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::dataSetRegionBoundaryFilter Class Reference

Like vtkDataSetRegionSurfaceFilter, this filter extracts boundaries of materials, including both interfaces and external boundaries. More...

Detailed Description

Unlike vtkDataSetRegionSurfaceFilter, can extract boundary edges for 2d data sets.

Output edges/faces are represented once for interfaces, with the orientation chosen by the input cell with lower material number. If input is a vtkPointSet, output points are shallow copy; otherwise, input points are deep copy - in either case, all points are copied, even those unused in output. Output point data is shallow copy of input point data; input cell data is not copied. Unlike vtkDataSetRegionSurfaceFilter, point ordering should match the edge/face described in OrigCellArrayName and CellSideArrayName arrays.

Definition at line 58 of file dataSetRegionBoundaryFilter.H.

Public Member Functions

 vtkTypeMacro (dataSetRegionBoundaryFilter, vtkPolyDataAlgorithm) vtkSetClampMacro(Dimension
 If 2-D, finds edges; if 3-D, finds faces. More...
 
 vtkGetMacro (Dimension, int)
 
 vtkSetMacro (MaterialArrayName, const std::string &)
 Input material/region array name. More...
 
 vtkGetMacro (MaterialArrayName, const std::string &)
 
 vtkSetMacro (BoundaryRegionArrayName, const std::string &)
 Output region array name. More...
 
 vtkGetMacro (BoundaryRegionArrayName, const std::string &)
 
 vtkSetMacro (RegionToMaterialArrayName, const std::string &)
 Output region to material array name. More...
 
 vtkGetMacro (RegionToMaterialArrayName, const std::string &)
 
 vtkSetMacro (OrigCellArrayName, const std::string &)
 Output original cell id array name. More...
 
 vtkGetMacro (OrigCellArrayName, const std::string &)
 
 vtkSetMacro (CellSideArrayName, const std::string &)
 Output original cell edge/face id array name. More...
 
 vtkGetMacro (CellSideArrayName, const std::string &)
 

Static Public Member Functions

static dataSetRegionBoundaryFilterNew ()
 

Public Attributes

 int
 

Protected Member Functions

 dataSetRegionBoundaryFilter ()=default
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 
int FillInputPortInformation (int port, vtkInformation *info) override
 

Protected Attributes

int Dimension {3}
 
std::string MaterialArrayName {"Material"}
 
std::string RegionToMaterialArrayName {"RegionToMaterial"}
 
std::string OrigCellArrayName {"OrigCellIds"}
 
std::string CellSideArrayName {"CellFaceIds"}
 
std::string BoundaryRegionArrayName {"Region"}
 

Inherits vtkPolyDataAlgorithm.

Constructor & Destructor Documentation

◆ dataSetRegionBoundaryFilter()

NEM::MSH::dataSetRegionBoundaryFilter::dataSetRegionBoundaryFilter ( )
protecteddefault

Member Function Documentation

◆ FillInputPortInformation()

int NEM::MSH::dataSetRegionBoundaryFilter::FillInputPortInformation ( int  port,
vtkInformation *  info 
)
overrideprotected

◆ New()

static dataSetRegionBoundaryFilter* NEM::MSH::dataSetRegionBoundaryFilter::New ( )
static

◆ RequestData()

int NEM::MSH::dataSetRegionBoundaryFilter::RequestData ( vtkInformation *  request,
vtkInformationVector **  inputVector,
vtkInformationVector *  outputVector 
)
overrideprotected

Definition at line 254 of file dataSetRegionBoundaryFilter.C.

References BoundaryCellSet, NEM::MSH::New(), and points.

256  {
257  // get the info objects
258  vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
259  vtkInformation *outInfo = outputVector->GetInformationObject(0);
260 
261  // get the input and output
262  auto input =
263  vtkDataSet::SafeDownCast(inInfo->Get(vtkDataObject::DATA_OBJECT()));
264  auto output =
265  vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
266 
267  if (input->CheckAttributes()) {
268  return 0;
269  }
270 
271  const vtkIdType numCells = input->GetNumberOfCells();
272  if (numCells == 0) {
273  vtkDebugMacro(<< "Number of cells is zero, no data to process.");
274  return 1;
275  }
276  vtkNew<vtkCellArray> newCells{};
277  // Approximate
278  newCells->Allocate(numCells);
279 
280  BoundaryCellSet boundaryCellSet(input->GetNumberOfPoints());
281  auto regionArr = vtkArrayDownCast<vtkIntArray>(
282  input->GetCellData()->GetAbstractArray(this->MaterialArrayName.c_str()));
283  for (vtkIdType i = 0; i < numCells; ++i) {
284  auto cell = input->GetCell(i);
285  auto region = regionArr->GetTypedComponent(i, 0);
286 
287  if (cell->GetCellDimension() != this->GetDimension()) {
288  continue;
289  }
290 
291  if (this->GetDimension() == 2) {
292  for (int j = 0; j < cell->GetNumberOfEdges(); ++j) {
293  auto edge = cell->GetEdge(j);
294  boundaryCellSet.insert(i, j, edge->GetPointIds(), 2, region);
295  }
296  } else { // this->GetDimension() == 3
297  switch (static_cast<VTKCellType>(cell->GetCellType())) {
298  case VTK_TETRA:
299  case VTK_QUADRATIC_TETRA:
300  case VTK_HIGHER_ORDER_TETRAHEDRON:
301  case VTK_LAGRANGE_TETRAHEDRON: {
302  auto points = cell->GetPointIds()->GetPointer(0);
303  boundaryCellSet.insertTetFaces(i, points, region);
304  break;
305  }
306  case VTK_HEXAHEDRON:
307  case VTK_QUADRATIC_HEXAHEDRON:
308  case VTK_TRIQUADRATIC_HEXAHEDRON:
309  case VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON:
310  case VTK_HIGHER_ORDER_HEXAHEDRON:
311  case VTK_LAGRANGE_HEXAHEDRON: {
312  auto points = cell->GetPointIds()->GetPointer(0);
313  boundaryCellSet.insertHexFaces(i, points, region);
314  break;
315  }
316  case VTK_WEDGE:
317  case VTK_QUADRATIC_WEDGE:
318  case VTK_QUADRATIC_LINEAR_WEDGE:
319  case VTK_BIQUADRATIC_QUADRATIC_WEDGE:
320  case VTK_HIGHER_ORDER_WEDGE:
321  case VTK_LAGRANGE_WEDGE: {
322  auto points = cell->GetPointIds()->GetPointer(0);
323  boundaryCellSet.insertWedgeFaces(i, points, region);
324  break;
325  }
326  default: {
327  for (int j = 0; j < cell->GetNumberOfFaces(); ++j) {
328  auto face = cell->GetFace(j);
329  // Treat non-linear faces as linear ones by only looking at
330  // beginning of point id array. Number of edges should be equal to
331  // number of vertices even for non-linear faces.
332  boundaryCellSet.insert(i, j, face->GetPointIds(),
333  face->GetNumberOfEdges(), region);
334  }
335  break;
336  }
337  }
338  }
339  }
340 
341  vtkSmartPointer<vtkIdTypeArray> origCellId{};
342  vtkSmartPointer<vtkIntArray> cellSideId{}, surfRegionId{};
343  if (!this->GetOrigCellArrayName().empty()) {
345  origCellId->SetName(this->GetOrigCellArrayName().c_str());
346  origCellId->SetNumberOfComponents(2);
347  origCellId->SetNumberOfTuples(boundaryCellSet.size());
348  }
349  if (!this->GetCellSideArrayName().empty()) {
350  cellSideId = vtkSmartPointer<vtkIntArray>::New();
351  cellSideId->SetName(this->GetCellSideArrayName().c_str());
352  cellSideId->SetNumberOfComponents(2);
353  cellSideId->SetNumberOfTuples(boundaryCellSet.size());
354  }
355  if (!this->GetBoundaryRegionArrayName().empty()) {
356  surfRegionId = vtkSmartPointer<vtkIntArray>::New();
357  surfRegionId->SetName(this->GetBoundaryRegionArrayName().c_str());
358  surfRegionId->SetNumberOfTuples(boundaryCellSet.size());
359  }
360 
361  // Map from regionArr to surfRegionId.
362  std::map<std::array<int, 2>, int> oldRegions2SurfRegions{};
363  for (const auto &boundaryCell : boundaryCellSet) {
364  auto origCellToInsert = input->GetCell(boundaryCell.sources[0]);
365  auto sideToInsert = this->GetDimension() == 2
366  ? origCellToInsert->GetEdge(boundaryCell.sides[0])
367  : origCellToInsert->GetFace(boundaryCell.sides[0]);
368  auto sideIdx = newCells->InsertNextCell(
369  this->GetDimension() == 2 ? 2 : sideToInsert->GetNumberOfEdges(),
370  sideToInsert->GetPointIds()->GetPointer(0));
371  if (origCellId) {
372  origCellId->SetTypedTuple(sideIdx, boundaryCell.sources.data());
373  }
374  if (cellSideId) {
375  cellSideId->SetTypedTuple(sideIdx, boundaryCell.sides.data());
376  }
377  auto surfRegionIter = oldRegions2SurfRegions.emplace(
378  boundaryCell.regions, oldRegions2SurfRegions.size());
379  if (surfRegionId) {
380  surfRegionId->SetTypedComponent(sideIdx, 0, surfRegionIter.first->second);
381  }
382  }
383  newCells->Squeeze();
384 
385  if (!this->GetRegionToMaterialArrayName().empty()) {
386  vtkNew<vtkIntArray> region2MaterialId;
387  region2MaterialId->SetName(this->GetRegionToMaterialArrayName().c_str());
388  region2MaterialId->SetNumberOfComponents(2);
389  region2MaterialId->SetNumberOfTuples(oldRegions2SurfRegions.size());
390  for (const auto &old2newRegion : oldRegions2SurfRegions) {
391  auto newRegion = old2newRegion.second;
392  region2MaterialId->SetComponent(newRegion, 0, old2newRegion.first[0]);
393  region2MaterialId->SetComponent(newRegion, 1, old2newRegion.first[1]);
394  }
395  output->GetFieldData()->AddArray(region2MaterialId);
396  }
397 
398  if (auto pointSet = vtkPointSet::SafeDownCast(input)) {
399  output->SetPoints(pointSet->GetPoints());
400  } else if (auto rectGrid = vtkRectilinearGrid::SafeDownCast(input)) {
401  rectGrid->GetPoints(output->GetPoints());
402  } else {
403  auto outPoints = output->GetPoints();
404  outPoints->Initialize();
405  outPoints->SetNumberOfPoints(input->GetNumberOfPoints());
406  for (vtkIdType i = 0; i < input->GetNumberOfPoints(); ++i) {
407  outPoints->SetPoint(i, input->GetPoint(i));
408  }
409  }
410  output->GetPointData()->ShallowCopy(input->GetPointData());
411 
412  if (this->GetDimension() == 2) {
413  output->SetLines(newCells);
414  } else { // this->GetDimension() == 3
415  output->SetPolys(newCells);
416  }
417 
418  vtkCellData *outputCD = output->GetCellData();
419  outputCD->AddArray(origCellId);
420  outputCD->AddArray(cellSideId);
421  outputCD->AddArray(surfRegionId);
422 
423  if (output->CheckAttributes()) {
424  return 0;
425  } else {
426  return 1;
427  }
428 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
friend BoundaryCellSet
std::vector< vtkIdType > points
points given by id in .inp file
Definition: inpGeoMesh.C:133

◆ vtkGetMacro() [1/6]

NEM::MSH::dataSetRegionBoundaryFilter::vtkGetMacro ( Dimension  ,
int   
)

◆ vtkGetMacro() [2/6]

NEM::MSH::dataSetRegionBoundaryFilter::vtkGetMacro ( MaterialArrayName  ,
const std::string &   
)

◆ vtkGetMacro() [3/6]

NEM::MSH::dataSetRegionBoundaryFilter::vtkGetMacro ( BoundaryRegionArrayName  ,
const std::string &   
)

◆ vtkGetMacro() [4/6]

NEM::MSH::dataSetRegionBoundaryFilter::vtkGetMacro ( RegionToMaterialArrayName  ,
const std::string &   
)

◆ vtkGetMacro() [5/6]

NEM::MSH::dataSetRegionBoundaryFilter::vtkGetMacro ( OrigCellArrayName  ,
const std::string &   
)

◆ vtkGetMacro() [6/6]

NEM::MSH::dataSetRegionBoundaryFilter::vtkGetMacro ( CellSideArrayName  ,
const std::string &   
)

◆ vtkSetMacro() [1/5]

NEM::MSH::dataSetRegionBoundaryFilter::vtkSetMacro ( MaterialArrayName  ,
const std::string &   
)

Default is "Material".

Assumed to be vtkIntArray in input cell data. Interfaces between materials are included in output. Should not use the value -1.

◆ vtkSetMacro() [2/5]

NEM::MSH::dataSetRegionBoundaryFilter::vtkSetMacro ( BoundaryRegionArrayName  ,
const std::string &   
)

Default is "Region".

If not empty, array will be added as vtkIntArray in output cell data. To recover relation to material, see RegionToMaterialArrayName

◆ vtkSetMacro() [3/5]

NEM::MSH::dataSetRegionBoundaryFilter::vtkSetMacro ( RegionToMaterialArrayName  ,
const std::string &   
)

Default is "RegionToMaterial".

If not empty, array will be added as vtkIntArray with two components in output field data. Indexed by region (value in BoundaryRegionArrayName array); value corresponds to material(s), with -1 corresponding to no material (for regions bounding only one material).

◆ vtkSetMacro() [4/5]

NEM::MSH::dataSetRegionBoundaryFilter::vtkSetMacro ( OrigCellArrayName  ,
const std::string &   
)

Default is "OrigCellIds".

If not empty, array will be added as vtkIdTypeArray with two components in output cell data. For each output edge/face, contains the index of the input cell. Second component is -1 if external edge/face.

◆ vtkSetMacro() [5/5]

NEM::MSH::dataSetRegionBoundaryFilter::vtkSetMacro ( CellSideArrayName  ,
const std::string &   
)

Default is "CellFaceIds".

If not empty, array will be added as vtkIntArray with two components in output cell data. For each output edge/face, contains the index of the corresponding edge/face of the cell given by the OrigCellArrayName array. Second component is -1 if external edge/face.

◆ vtkTypeMacro()

NEM::MSH::dataSetRegionBoundaryFilter::vtkTypeMacro ( dataSetRegionBoundaryFilter  ,
vtkPolyDataAlgorithm   
)

Default is 3.

Member Data Documentation

◆ BoundaryRegionArrayName

std::string NEM::MSH::dataSetRegionBoundaryFilter::BoundaryRegionArrayName {"Region"}
protected

Definition at line 139 of file dataSetRegionBoundaryFilter.H.

◆ CellSideArrayName

std::string NEM::MSH::dataSetRegionBoundaryFilter::CellSideArrayName {"CellFaceIds"}
protected

Definition at line 138 of file dataSetRegionBoundaryFilter.H.

◆ Dimension

int NEM::MSH::dataSetRegionBoundaryFilter::Dimension {3}
protected

Definition at line 134 of file dataSetRegionBoundaryFilter.H.

◆ int

NEM::MSH::dataSetRegionBoundaryFilter::int

Definition at line 67 of file dataSetRegionBoundaryFilter.H.

◆ MaterialArrayName

std::string NEM::MSH::dataSetRegionBoundaryFilter::MaterialArrayName {"Material"}
protected

Definition at line 135 of file dataSetRegionBoundaryFilter.H.

◆ OrigCellArrayName

std::string NEM::MSH::dataSetRegionBoundaryFilter::OrigCellArrayName {"OrigCellIds"}
protected

Definition at line 137 of file dataSetRegionBoundaryFilter.H.

◆ RegionToMaterialArrayName

std::string NEM::MSH::dataSetRegionBoundaryFilter::RegionToMaterialArrayName {"RegionToMaterial"}
protected

Definition at line 136 of file dataSetRegionBoundaryFilter.H.


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