NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
NEM::MSH::geoMeshBase::GeoMesh Struct Reference

Detailed Description

Definition at line 418 of file geoMeshBase.H.

Public Member Functions

void findSide2OrigCell ()
 For each edge/polygon in sideSet, find the corresponding mesh element in mesh. More...
 
int getDimension () const
 Get dimension of mesh. More...
 

Public Attributes

vtkSmartPointer< vtkUnstructuredGrid > mesh
 
std::string geo
 
std::string link
 
SideSet sideSet
 

Member Function Documentation

◆ findSide2OrigCell()

void NEM::MSH::geoMeshBase::GeoMesh::findSide2OrigCell ( )

Does nothing if "OrigCellIds"/"CellFaceIds" arrays already present. Sets the sideSet's "OrigCellIds"/"CellFaceIds" arrays. Orientation of side only checked if multiple candidate original cells. Both "OrigCellIds" and "CellFaceIds" are set to -1 for any cells in the sideSet that can't be matched.

Definition at line 314 of file geoMeshBase.C.

References mesh.

314  {
316  vtkIdType numSides = sideSet.sides->GetNumberOfCells();
317  auto entArr = sideSet.getGeoEntArr();
318  vtkNew<vtkIdTypeArray> origCellArr;
319  origCellArr->SetNumberOfComponents(2);
320  origCellArr->SetNumberOfTuples(numSides);
321  origCellArr->FillTypedComponent(0, -1);
322  origCellArr->FillTypedComponent(1, -1);
323  vtkNew<vtkIntArray> cellFaceArr;
324  cellFaceArr->SetNumberOfComponents(2);
325  cellFaceArr->SetNumberOfTuples(numSides);
326  cellFaceArr->FillTypedComponent(0, -1);
327  cellFaceArr->FillTypedComponent(1, -1);
328  // Using points->cell, intersect the set of cells that have the same points
329  // as the side cell
330  mesh->BuildLinks();
331  auto pt2Cell = mesh->GetCellLinks();
332  for (vtkIdType i = 0; i < sideSet.sides->GetNumberOfCells(); ++i) {
333  auto sidePoints = sideSet.sides->GetCell(i)->GetPointIds();
334  auto cellList = pt2Cell->GetLink(sidePoints->GetId(0));
335  std::vector<vtkIdType> candidates{cellList.cells,
336  cellList.cells + cellList.ncells};
337  std::sort(candidates.begin(), candidates.end());
338  for (auto j = 1; j < sidePoints->GetNumberOfIds(); ++j) {
339  std::vector<vtkIdType> newCandidates{};
340  cellList = pt2Cell->GetLink(sidePoints->GetId(j));
341  for (auto iterCell = cellList.cells;
342  iterCell < cellList.cells + cellList.ncells; ++iterCell) {
343  auto findIter =
344  std::lower_bound(candidates.begin(), candidates.end(), *iterCell);
345  if (findIter != candidates.end() && *iterCell == *findIter) {
346  newCandidates.emplace_back(*iterCell);
347  }
348  }
349  candidates = std::move(newCandidates);
350  std::sort(candidates.begin(), candidates.end());
351  }
352  // For each candidate cell, check each edge/face for a match
353  bool foundMatch = false;
354  if (candidates.size() == 1) {
355  // If there is only one candidate (side is on boundary, not material
356  // interface), check by sorting points - note no check for orientation
357  std::vector<vtkIdType> sortedSidePoints(sidePoints->GetNumberOfIds());
358  auto pointsBegin = sidePoints->GetPointer(0);
359  std::partial_sort_copy(
360  pointsBegin, pointsBegin + sidePoints->GetNumberOfIds(),
361  sortedSidePoints.begin(), sortedSidePoints.end());
362  auto candOrigCell = mesh->GetCell(candidates.at(0));
363  auto dim = candOrigCell->GetCellDimension();
364  for (int j = 0; j < (dim == 3 ? candOrigCell->GetNumberOfFaces()
365  : candOrigCell->GetNumberOfEdges());
366  ++j) {
367  auto candSide =
368  dim == 3 ? candOrigCell->GetFace(j) : candOrigCell->GetEdge(j);
369  std::vector<vtkIdType> sortedCandPoints(sortedSidePoints.size());
370  auto candPointsBegin = candSide->GetPointIds()->GetPointer(0);
371  // Note use of sidePoints->GetNumberOfIds() in case, eg, candSide is
372  // a vtkQuadraticTriangle, but side is a vtkTriangle
373  auto candPointsEnd = candPointsBegin + sidePoints->GetNumberOfIds();
374  std::partial_sort_copy(candPointsBegin, candPointsEnd,
375  sortedCandPoints.begin(),
376  sortedCandPoints.end());
377  if (sortedSidePoints == sortedCandPoints) {
378  foundMatch = true;
379  origCellArr->SetTypedComponent(i, 0, candidates.at(0));
380  cellFaceArr->SetTypedComponent(i, 0, j);
381  break;
382  }
383  }
384  } else {
385  // If there are two candidates, match the orientation of the side
386  assert(candidates.size() == 2);
387  for (const auto &candIdx : candidates) {
388  auto candOrigCell = mesh->GetCell(candIdx);
389  auto dim = candOrigCell->GetCellDimension();
390  bool foundSide = false;
391  for (int j = 0; j < (dim == 3 ? candOrigCell->GetNumberOfFaces()
392  : candOrigCell->GetNumberOfEdges()) &&
393  !foundSide;
394  ++j) {
395  auto candSide =
396  dim == 3 ? candOrigCell->GetFace(j) : candOrigCell->GetEdge(j);
397  auto candPointsBegin = candSide->GetPointIds()->GetPointer(0);
398  // Note use of sidePoints->GetNumberOfIds() in case, eg, candSide is
399  // a vtkQuadraticTriangle, but side is a vtkTriangle
400  auto candPointsEnd = candPointsBegin + sidePoints->GetNumberOfIds();
401  auto findFirstPoint =
402  std::find(candPointsBegin, candPointsEnd, sidePoints->GetId(0));
403  if (findFirstPoint != candPointsEnd) {
404  std::vector<vtkIdType> rotatedCandPoints(
405  sidePoints->GetNumberOfIds());
406  std::rotate_copy(candPointsBegin, findFirstPoint, candPointsEnd,
407  rotatedCandPoints.begin());
408  if (std::equal(rotatedCandPoints.begin(), rotatedCandPoints.end(),
409  sidePoints->GetPointer(0))) {
410  foundMatch = true;
411  foundSide = true;
412  origCellArr->SetTypedComponent(i, 0, candIdx);
413  cellFaceArr->SetTypedComponent(i, 0, j);
414  } else if (std::equal(rotatedCandPoints.rbegin(),
415  rotatedCandPoints.rend(),
416  sidePoints->GetPointer(0))) {
417  foundSide = true;
418  origCellArr->SetTypedComponent(i, 1, candIdx);
419  cellFaceArr->SetTypedComponent(i, 1, j);
420  }
421  }
422  }
423  if (foundMatch) { break; }
424  }
425  }
426  if (!foundMatch) {
427  std::cerr << "No cell in mesh found to correspond to side set cell "
428  << i << ".\n";
429  }
430  }
431  sideSet.setOrigCellArr(origCellArr);
432  sideSet.setCellFaceArr(cellFaceArr);
433  }
434 }
vtkSmartPointer< vtkUnstructuredGrid > mesh
Definition: geoMeshBase.H:419
vtkSmartPointer< vtkIntArray > getCellFaceArr() const
Definition: geoMeshBase.C:281
vtkSmartPointer< vtkIdTypeArray > getOrigCellArr() const
Definition: geoMeshBase.C:264
void setCellFaceArr(vtkIntArray *arr)
Definition: geoMeshBase.C:287
void setOrigCellArr(vtkIdTypeArray *arr)
Definition: geoMeshBase.C:270
vtkSmartPointer< vtkPolyData > sides
Cells represent edges/faces of some GeoMesh.
Definition: geoMeshBase.H:399
vtkSmartPointer< vtkIntArray > getGeoEntArr() const
Definition: geoMeshBase.C:252

◆ getDimension()

int NEM::MSH::geoMeshBase::GeoMesh::getDimension ( ) const

All cells in mesh assumed to have same dimension.

Returns
dimension of mesh

Definition at line 436 of file geoMeshBase.C.

References mesh.

Referenced by NEM::MSH::exoGeoMesh::stitch().

436  {
437  if (this->mesh->GetNumberOfCells() > 0) {
438  return this->mesh->GetCell(0)->GetCellDimension();
439  } else {
440  return -1;
441  }
442 }
vtkSmartPointer< vtkUnstructuredGrid > mesh
Definition: geoMeshBase.H:419

Member Data Documentation

◆ geo

◆ link

◆ mesh

◆ sideSet


The documentation for this struct was generated from the following files: