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

Detailed Description

Definition at line 55 of file mergeCells.H.

Public Member Functions

 vtkTypeMacro (mergeCells, vtkObject) void PrintSelf(ostream &os
 
 vtkSetObjectMacro (UnstructuredGrid, vtkUnstructuredGrid)
 
 vtkGetObjectMacro (UnstructuredGrid, vtkUnstructuredGrid)
 
 vtkSetMacro (TotalNumberOfCells, vtkIdType)
 
 vtkGetMacro (TotalNumberOfCells, vtkIdType)
 
 vtkSetMacro (TotalNumberOfPoints, vtkIdType)
 
 vtkGetMacro (TotalNumberOfPoints, vtkIdType)
 
 vtkSetMacro (UseGlobalIds, int)
 
 vtkGetMacro (UseGlobalIds, int)
 
 vtkSetClampMacro (PointMergeTolerance, float, 0.0f, VTK_FLOAT_MAX)
 
 vtkGetMacro (PointMergeTolerance, float)
 
 vtkSetMacro (UseGlobalCellIds, int)
 
 vtkGetMacro (UseGlobalCellIds, int)
 
 vtkSetMacro (MergeDuplicatePoints, int)
 
 vtkGetMacro (MergeDuplicatePoints, int)
 
 vtkBooleanMacro (MergeDuplicatePoints, int)
 
 vtkSetMacro (TotalNumberOfDataSets, int)
 
 vtkGetMacro (TotalNumberOfDataSets, int)
 
vtkSmartPointer< vtkIdTypeArray > MergeDataSet (vtkDataSet *set)
 Provide a DataSet to be merged in to the final UnstructuredGrid. More...
 
void Finish ()
 

Static Public Member Functions

static mergeCellsNew ()
 

Public Attributes

vtkIndent indent override
 

Protected Member Functions

 mergeCells ()
 
 ~mergeCells () VTK_OVERRIDE
 
void FreeLists ()
 
void StartUGrid (vtkDataSet *set)
 
vtkIdTypeArray * MapPointsToIdsUsingGlobalIds (vtkDataSet *set)
 
vtkIdTypeArray * MapPointsToIdsUsingLocator (vtkDataSet *set)
 
vtkIdType AddNewCellsUnstructuredGrid (vtkDataSet *set, vtkIdTypeArray *idMap)
 
vtkIdType AddNewCellsDataSet (vtkDataSet *set, vtkIdTypeArray *idMap)
 
vtkIdType GlobalCellIdAccessGetId (vtkIdType idx)
 
int GlobalCellIdAccessStart (vtkDataSet *set)
 
vtkIdType GlobalNodeIdAccessGetId (vtkIdType idx)
 
int GlobalNodeIdAccessStart (vtkDataSet *set)
 
 mergeCells (const mergeCells &)=delete
 
void operator= (const mergeCells &)=delete
 

Protected Attributes

int TotalNumberOfDataSets
 
vtkIdType TotalNumberOfCells
 
vtkIdType TotalNumberOfPoints
 
vtkIdType NumberOfCells
 
vtkIdType NumberOfPoints
 
int UseGlobalIds
 
int GlobalIdArrayType
 
void * GlobalIdArray
 
int UseGlobalCellIds
 
int GlobalCellIdArrayType
 
void * GlobalCellIdArray
 
float PointMergeTolerance
 
int MergeDuplicatePoints
 
char InputIsUGrid
 
char InputIsPointSet
 
mergeCellsSTLCloak * GlobalIdMap
 
mergeCellsSTLCloak * GlobalCellIdMap
 
vtkDataSetAttributes::FieldList * ptList
 
vtkDataSetAttributes::FieldList * cellList
 
vtkUnstructuredGrid * UnstructuredGrid
 
int nextGrid
 

Inherits vtkObject.

Constructor & Destructor Documentation

◆ mergeCells() [1/2]

NEM::MSH::mergeCells::mergeCells ( )
protected

Definition at line 66 of file mergeCells.C.

66  {
67  this->TotalNumberOfDataSets = 0;
68  this->TotalNumberOfCells = 0;
69  this->TotalNumberOfPoints = 0;
70 
71  this->NumberOfCells = 0;
72  this->NumberOfPoints = 0;
73 
74  this->PointMergeTolerance = 10e-4f;
75  this->MergeDuplicatePoints = 1;
76 
77  this->InputIsUGrid = 0;
78  this->InputIsPointSet = 0;
79 
80  this->ptList = NULL;
81  this->cellList = NULL;
82 
83  this->UnstructuredGrid = NULL;
84 
85  this->GlobalIdMap = new mergeCellsSTLCloak;
86  this->GlobalCellIdMap = new mergeCellsSTLCloak;
87 
88  this->UseGlobalIds = 0;
89  this->UseGlobalCellIds = 0;
90 
91  this->nextGrid = 0;
92 }
vtkIdType TotalNumberOfCells
Definition: mergeCells.H:116
vtkDataSetAttributes::FieldList * cellList
Definition: mergeCells.H:140
vtkIdType NumberOfPoints
Definition: mergeCells.H:120
vtkUnstructuredGrid * UnstructuredGrid
Definition: mergeCells.H:142
vtkDataSetAttributes::FieldList * ptList
Definition: mergeCells.H:139
vtkIdType NumberOfCells
Definition: mergeCells.H:119
mergeCellsSTLCloak * GlobalIdMap
Definition: mergeCells.H:136
mergeCellsSTLCloak * GlobalCellIdMap
Definition: mergeCells.H:137
vtkIdType TotalNumberOfPoints
Definition: mergeCells.H:117

◆ ~mergeCells()

NEM::MSH::mergeCells::~mergeCells ( )
protected

Definition at line 94 of file mergeCells.C.

94  {
95  this->FreeLists();
96 
97  delete this->GlobalIdMap;
98  delete this->GlobalCellIdMap;
99 
100  this->SetUnstructuredGrid(0);
101 }
mergeCellsSTLCloak * GlobalIdMap
Definition: mergeCells.H:136
mergeCellsSTLCloak * GlobalCellIdMap
Definition: mergeCells.H:137

◆ mergeCells() [2/2]

NEM::MSH::mergeCells::mergeCells ( const mergeCells )
protecteddelete

Member Function Documentation

◆ AddNewCellsDataSet()

vtkIdType NEM::MSH::mergeCells::AddNewCellsDataSet ( vtkDataSet *  set,
vtkIdTypeArray *  idMap 
)
protected

Definition at line 221 of file mergeCells.C.

References id, and NEM::MSH::New().

222  {
223  vtkIdType oldCellId, id, newPtId, newCellId = 0, oldPtId;
224 
225  vtkUnstructuredGrid *ugrid = this->UnstructuredGrid;
226  vtkCellData *cellArrays = set->GetCellData();
227  vtkIdType numCells = set->GetNumberOfCells();
228 
229  vtkIdList *cellPoints = vtkIdList::New();
230  cellPoints->Allocate(VTK_CELL_SIZE);
231 
232  vtkIdType nextCellId = 0;
233 
234  int duplicateCellTest = 0;
235 
236  if (this->UseGlobalCellIds) {
237  int success = this->GlobalCellIdAccessStart(set);
238 
239  if (success) {
240  nextCellId = this->GlobalCellIdMap->IdTypeMap.size();
241  duplicateCellTest = 1;
242  }
243  }
244 
245  for (oldCellId = 0; oldCellId < numCells; oldCellId++) {
246  if (duplicateCellTest) {
247  vtkIdType globalId = this->GlobalCellIdAccessGetId(oldCellId);
248 
249  std::pair<std::map<vtkIdType, vtkIdType>::iterator, bool> inserted =
250 
251  this->GlobalCellIdMap->IdTypeMap.insert(
252  std::map<vtkIdType, vtkIdType>::value_type(globalId, nextCellId));
253 
254  if (inserted.second) {
255  nextCellId++;
256  } else {
257  continue; // skip it, we already have this cell
258  }
259  }
260 
261  set->GetCellPoints(oldCellId, cellPoints);
262 
263  for (id = 0; id < cellPoints->GetNumberOfIds(); id++) {
264  oldPtId = cellPoints->GetId(id);
265 
266  if (idMap) {
267  newPtId = idMap->GetValue(oldPtId);
268  } else {
269  newPtId = this->NumberOfPoints + oldPtId;
270  }
271  cellPoints->SetId(id, newPtId);
272  }
273 
274  newCellId = (vtkIdType)ugrid->InsertNextCell(set->GetCellType(oldCellId),
275  cellPoints);
276 
277  ugrid->GetCellData()->CopyData(*(this->cellList), cellArrays,
278  this->nextGrid, oldCellId, newCellId);
279  }
280 
281  cellPoints->Delete();
282 
283  return newCellId;
284 }
vtkIdType GlobalCellIdAccessGetId(vtkIdType idx)
Definition: mergeCells.C:773
vtkDataSetAttributes::FieldList * cellList
Definition: mergeCells.H:140
int GlobalCellIdAccessStart(vtkDataSet *set)
Definition: mergeCells.C:783
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkIdType NumberOfPoints
Definition: mergeCells.H:120
vtkUnstructuredGrid * UnstructuredGrid
Definition: mergeCells.H:142
vtkIdType id
id in .inp file
Definition: inpGeoMesh.C:128
mergeCellsSTLCloak * GlobalCellIdMap
Definition: mergeCells.H:137

◆ AddNewCellsUnstructuredGrid()

vtkIdType NEM::MSH::mergeCells::AddNewCellsUnstructuredGrid ( vtkDataSet *  set,
vtkIdTypeArray *  idMap 
)
protected

Definition at line 285 of file mergeCells.C.

References id, and NEM::MSH::New().

286  {
287  vtkIdType id;
288 
289  char firstSet = 0;
290 
291  if (this->nextGrid == 0) firstSet = 1;
292 
293  vtkUnstructuredGrid *newUgrid = vtkUnstructuredGrid::SafeDownCast(set);
294  vtkUnstructuredGrid *Ugrid = this->UnstructuredGrid;
295 
296  // connectivity information for the new data set
297 
298  vtkCellArray *newCellArray = newUgrid->GetCells();
299  vtkIdType *newCells = newCellArray->GetPointer();
300  vtkIdType *newLocs = newUgrid->GetCellLocationsArray()->GetPointer(0);
301  unsigned char *newTypes = newUgrid->GetCellTypesArray()->GetPointer(0);
302 
303  int newNumCells = newUgrid->GetNumberOfCells();
304  int newNumConnections = newCellArray->GetData()->GetNumberOfTuples();
305 
306  // If we are checking for duplicate cells, create a list now of
307  // any cells in the new data set that we already have.
308 
309  vtkIdList *duplicateCellIds = NULL;
310  int numDuplicateCells = 0;
311  int numDuplicateConnections = 0;
312 
313  if (this->UseGlobalCellIds) {
314  int success = this->GlobalCellIdAccessStart(set);
315 
316  if (success) {
317  vtkIdType nextLocalId = this->GlobalCellIdMap->IdTypeMap.size();
318 
319  duplicateCellIds = vtkIdList::New();
320 
321  for (id = 0; id < newNumCells; id++) {
322  vtkIdType globalId = this->GlobalCellIdAccessGetId(id);
323 
324  std::pair<std::map<vtkIdType, vtkIdType>::iterator, bool> inserted =
325 
326  this->GlobalCellIdMap->IdTypeMap.insert(
327  std::map<vtkIdType, vtkIdType>::value_type(globalId,
328  nextLocalId));
329 
330  if (inserted.second) {
331  nextLocalId++;
332  } else {
333  duplicateCellIds->InsertNextId(id);
334  numDuplicateCells++;
335 
336  int npoints = newCells[newLocs[id]];
337 
338  numDuplicateConnections += (npoints + 1);
339  }
340  }
341 
342  if (numDuplicateCells == 0) {
343  duplicateCellIds->Delete();
344  duplicateCellIds = NULL;
345  }
346  }
347  }
348 
349  // connectivity for the merged ugrid so far
350 
351  vtkCellArray *cellArray = NULL;
352  vtkIdType *cells = NULL;
353  vtkIdType *locs = NULL;
354  unsigned char *types = NULL;
355 
356  int numCells = 0;
357  int numConnections = 0;
358 
359  if (!firstSet) {
360  cellArray = Ugrid->GetCells();
361  cells = cellArray->GetPointer();
362  locs = Ugrid->GetCellLocationsArray()->GetPointer(0);
363  types = Ugrid->GetCellTypesArray()->GetPointer(0);
364  ;
365 
366  numCells = Ugrid->GetNumberOfCells();
367  numConnections = cellArray->GetData()->GetNumberOfTuples();
368  }
369 
370  // New output grid: merging of existing and incoming grids
371 
372  // CELL ARRAY
373 
374  int totalNumCells = numCells + newNumCells - numDuplicateCells;
375  int totalNumConnections =
376  numConnections + newNumConnections - numDuplicateConnections;
377 
378  vtkIdTypeArray *mergedcells = vtkIdTypeArray::New();
379  mergedcells->SetNumberOfValues(totalNumConnections);
380 
381  if (!firstSet) {
382  vtkIdType *idptr = mergedcells->GetPointer(0);
383  memcpy(idptr, cells, sizeof(vtkIdType) * numConnections);
384  }
385 
386  vtkCellArray *finalCellArray = vtkCellArray::New();
387  finalCellArray->SetCells(totalNumCells, mergedcells);
388 
389  // LOCATION ARRAY
390 
391  vtkIdTypeArray *locationArray = vtkIdTypeArray::New();
392  locationArray->SetNumberOfValues(totalNumCells);
393 
394  vtkIdType *iptr = locationArray->GetPointer(0); // new output dataset
395 
396  if (!firstSet) {
397  memcpy(iptr, locs, numCells * sizeof(vtkIdType)); // existing set
398  }
399 
400  // TYPE ARRAY
401 
402  vtkUnsignedCharArray *typeArray = vtkUnsignedCharArray::New();
403  typeArray->SetNumberOfValues(totalNumCells);
404 
405  unsigned char *cptr = typeArray->GetPointer(0);
406 
407  if (!firstSet) {
408  memcpy(cptr, types, numCells * sizeof(unsigned char));
409  }
410 
411  // set up new cell data
412 
413  vtkIdType finalCellId = numCells;
414  vtkIdType nextCellArrayIndex = static_cast<vtkIdType>(numConnections);
415  vtkCellData *cellArrays = set->GetCellData();
416 
417  vtkIdType oldPtId, finalPtId;
418 
419  int nextDuplicateCellId = 0;
420 
421  for (vtkIdType oldCellId = 0; oldCellId < newNumCells; oldCellId++) {
422  vtkIdType size = *newCells++;
423 
424  if (duplicateCellIds) {
425  vtkIdType skipId = duplicateCellIds->GetId(nextDuplicateCellId);
426 
427  if (skipId == oldCellId) {
428  newCells += size;
429  nextDuplicateCellId++;
430  continue;
431  }
432  }
433 
434  locationArray->SetValue(finalCellId, nextCellArrayIndex);
435 
436  typeArray->SetValue(finalCellId, newTypes[oldCellId]);
437 
438  mergedcells->SetValue(nextCellArrayIndex++, size);
439 
440  for (id = 0; id < size; id++) {
441  oldPtId = *newCells++;
442 
443  if (idMap) {
444  finalPtId = idMap->GetValue(oldPtId);
445  } else {
446  finalPtId = this->NumberOfPoints + oldPtId;
447  }
448 
449  mergedcells->SetValue(nextCellArrayIndex++, finalPtId);
450  }
451 
452  Ugrid->GetCellData()->CopyData(*(this->cellList), cellArrays,
453  this->nextGrid, oldCellId, finalCellId);
454 
455  finalCellId++;
456  }
457 
458  Ugrid->SetCells(typeArray, locationArray, finalCellArray);
459 
460  mergedcells->Delete();
461  typeArray->Delete();
462  locationArray->Delete();
463  finalCellArray->Delete();
464 
465  if (duplicateCellIds) {
466  duplicateCellIds->Delete();
467  }
468 
469  return finalCellId;
470 }
vtkIdType GlobalCellIdAccessGetId(vtkIdType idx)
Definition: mergeCells.C:773
vtkDataSetAttributes::FieldList * cellList
Definition: mergeCells.H:140
int GlobalCellIdAccessStart(vtkDataSet *set)
Definition: mergeCells.C:783
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkIdType NumberOfPoints
Definition: mergeCells.H:120
vtkUnstructuredGrid * UnstructuredGrid
Definition: mergeCells.H:142
vtkIdType id
id in .inp file
Definition: inpGeoMesh.C:128
mergeCellsSTLCloak * GlobalCellIdMap
Definition: mergeCells.H:137

◆ Finish()

void NEM::MSH::mergeCells::Finish ( )

Definition at line 524 of file mergeCells.C.

524  {
525  this->FreeLists();
526 
527  vtkUnstructuredGrid *ugrid = this->UnstructuredGrid;
528 
529  if (this->NumberOfPoints < this->TotalNumberOfPoints) {
530  // if we don't do this, ugrid->GetNumberOfPoints() gives
531  // the wrong value
532 
533  ugrid->GetPoints()->GetData()->Resize(this->NumberOfPoints);
534  }
535 
536  ugrid->Squeeze();
537 
538  return;
539 }
vtkIdType NumberOfPoints
Definition: mergeCells.H:120
vtkUnstructuredGrid * UnstructuredGrid
Definition: mergeCells.H:142
vtkIdType TotalNumberOfPoints
Definition: mergeCells.H:117

◆ FreeLists()

void NEM::MSH::mergeCells::FreeLists ( )
protected

Definition at line 103 of file mergeCells.C.

103  {
104  delete this->ptList;
105  this->ptList = NULL;
106 
107  delete this->cellList;
108  this->cellList = NULL;
109 }
vtkDataSetAttributes::FieldList * cellList
Definition: mergeCells.H:140
vtkDataSetAttributes::FieldList * ptList
Definition: mergeCells.H:139

◆ GlobalCellIdAccessGetId()

vtkIdType NEM::MSH::mergeCells::GlobalCellIdAccessGetId ( vtkIdType  idx)
protected

Definition at line 773 of file mergeCells.C.

773  {
774  if (this->GlobalCellIdArray) {
775  switch (this->GlobalCellIdArrayType) {
776  vtkTemplateMacro(VTK_TT *ids =
777  static_cast<VTK_TT *>(this->GlobalCellIdArray);
778  return static_cast<vtkIdType>(ids[idx]));
779  }
780  }
781  return 0;
782 }

◆ GlobalCellIdAccessStart()

int NEM::MSH::mergeCells::GlobalCellIdAccessStart ( vtkDataSet *  set)
protected

Definition at line 783 of file mergeCells.C.

783  {
784  if (this->UseGlobalCellIds) {
785  vtkDataArray *da = set->GetCellData()->GetGlobalIds();
786  if (da) {
787  this->GlobalCellIdArray = da->GetVoidPointer(0);
788  this->GlobalCellIdArrayType = da->GetDataType();
789  return 1;
790  }
791  }
792 
793  this->GlobalCellIdArray = 0;
794  return 0;
795 }

◆ GlobalNodeIdAccessGetId()

vtkIdType NEM::MSH::mergeCells::GlobalNodeIdAccessGetId ( vtkIdType  idx)
protected

Definition at line 797 of file mergeCells.C.

797  {
798  if (this->GlobalIdArray) {
799  switch (this->GlobalIdArrayType) {
800  vtkTemplateMacro(VTK_TT *ids = static_cast<VTK_TT *>(this->GlobalIdArray);
801  return static_cast<vtkIdType>(ids[idx]));
802  }
803  }
804  return 0;
805 }

◆ GlobalNodeIdAccessStart()

int NEM::MSH::mergeCells::GlobalNodeIdAccessStart ( vtkDataSet *  set)
protected

Definition at line 807 of file mergeCells.C.

807  {
808  if (this->UseGlobalIds) {
809  vtkDataArray *da = set->GetPointData()->GetGlobalIds();
810  if (da) {
811  this->GlobalIdArray = da->GetVoidPointer(0);
812  this->GlobalIdArrayType = da->GetDataType();
813  return 1;
814  }
815  }
816 
817  this->GlobalIdArray = 0;
818  return 0;
819 }

◆ MapPointsToIdsUsingGlobalIds()

vtkIdTypeArray * NEM::MSH::mergeCells::MapPointsToIdsUsingGlobalIds ( vtkDataSet *  set)
protected

Definition at line 544 of file mergeCells.C.

References NEM::MSH::New().

544  {
545  int success = this->GlobalNodeIdAccessStart(set);
546 
547  if (!success) {
548  vtkErrorMacro("global id array is not available");
549  return NULL;
550  }
551 
552  vtkIdType npoints = set->GetNumberOfPoints();
553 
554  auto idMap = vtkIdTypeArray::New();
555  idMap->SetNumberOfValues(npoints);
556 
557  vtkIdType nextNewLocalId = this->GlobalIdMap->IdTypeMap.size();
558 
559  // map global point Ids to Ids in the new data set
560 
561  for (vtkIdType oldId = 0; oldId < npoints; oldId++) {
562  vtkIdType globalId = this->GlobalNodeIdAccessGetId(oldId);
563 
564  std::pair<std::map<vtkIdType, vtkIdType>::iterator, bool> inserted =
565 
566  this->GlobalIdMap->IdTypeMap.insert(
567  std::map<vtkIdType, vtkIdType>::value_type(globalId,
568  nextNewLocalId));
569 
570  if (inserted.second) {
571  // this is a new global node Id
572 
573  idMap->SetValue(oldId, nextNewLocalId);
574 
575  nextNewLocalId++;
576  } else {
577  // a repeat, it was not inserted
578 
579  idMap->SetValue(oldId, inserted.first->second);
580  }
581  }
582 
583  return idMap;
584 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkIdType GlobalNodeIdAccessGetId(vtkIdType idx)
Definition: mergeCells.C:797
int GlobalNodeIdAccessStart(vtkDataSet *set)
Definition: mergeCells.C:807
mergeCellsSTLCloak * GlobalIdMap
Definition: mergeCells.H:136

◆ MapPointsToIdsUsingLocator()

vtkIdTypeArray * NEM::MSH::mergeCells::MapPointsToIdsUsingLocator ( vtkDataSet *  set)
protected

Definition at line 589 of file mergeCells.C.

References NEM::MSH::New().

589  {
590  vtkIdType ptId;
591 
592  vtkUnstructuredGrid *grid = this->UnstructuredGrid;
593  vtkPoints *points0 = grid->GetPoints();
594  vtkIdType npoints0 = (vtkIdType)this->NumberOfPoints;
595 
596  vtkPointSet *ps = vtkPointSet::SafeDownCast(set);
597  vtkPoints *points1;
598  vtkIdType npoints1 = set->GetNumberOfPoints();
599 
600  if (ps) {
601  points1 = ps->GetPoints();
602  } else {
603  points1 = vtkPoints::New();
604  points1->SetNumberOfPoints(npoints1);
605 
606  for (ptId = 0; ptId < npoints1; ptId++) {
607  points1->SetPoint(ptId, set->GetPoint(ptId));
608  }
609  }
610 
611  auto idMap = vtkIdTypeArray::New();
612  idMap->SetNumberOfValues(npoints1);
613 
614  vtkIdType nextNewLocalId = npoints0;
615 
616  if (this->PointMergeTolerance == 0.0) {
617  // testing shows vtkMergePoints is fastest when tolerance is 0
618 
619  vtkMergePoints *locator = vtkMergePoints::New();
620 
621  vtkPoints *ptarray = vtkPoints::New();
622 
623  double bounds[6];
624 
625  set->GetBounds(bounds);
626 
627  if (npoints0 > 0) {
628  double tmpbounds[6];
629 
630  // Prior to MapPointsToIdsUsingLocator(), points0->SetNumberOfPoints()
631  // has been called to set the number of points to the upper bound on the
632  // points TO BE merged and now points0->GetNumberOfPoints() does not
633  // refer to the number of the points merged so far. Thus we need to
634  // temporarily set the number to the latter such that grid->GetBounds()
635  // is able to return the correct bounding information. This is a fix to
636  // bug #0009626.
637 
638  points0->GetData()->SetNumberOfTuples(npoints0);
639  grid->GetBounds(tmpbounds); // safe to call GetBounds() for real info
640  points0->GetData()->SetNumberOfTuples(this->TotalNumberOfPoints);
641 
642  bounds[0] = ((tmpbounds[0] < bounds[0]) ? tmpbounds[0] : bounds[0]);
643  bounds[2] = ((tmpbounds[2] < bounds[2]) ? tmpbounds[2] : bounds[2]);
644  bounds[4] = ((tmpbounds[4] < bounds[4]) ? tmpbounds[4] : bounds[4]);
645 
646  bounds[1] = ((tmpbounds[1] > bounds[1]) ? tmpbounds[1] : bounds[1]);
647  bounds[3] = ((tmpbounds[3] > bounds[3]) ? tmpbounds[3] : bounds[3]);
648  bounds[5] = ((tmpbounds[5] > bounds[5]) ? tmpbounds[5] : bounds[5]);
649  }
650 
651  locator->InitPointInsertion(ptarray, bounds);
652 
653  vtkIdType newId;
654  double x[3];
655 
656  for (ptId = 0; ptId < npoints0; ptId++) {
657  // We already know there are no duplicates in this array.
658  // Just add them to the locator's point array.
659 
660  points0->GetPoint(ptId, x);
661  locator->InsertUniquePoint(x, newId);
662  }
663  for (ptId = 0; ptId < npoints1; ptId++) {
664  points1->GetPoint(ptId, x);
665  locator->InsertUniquePoint(x, newId);
666 
667  idMap->SetValue(ptId, newId);
668  }
669 
670  locator->Delete();
671  ptarray->Delete();
672  } else {
673  // testing shows vtkKdTree is fastest when tolerance is > 0
674 
675  vtkKdTree *kd = vtkKdTree::New();
676 
677  vtkPoints *ptArrays[2];
678  int numArrays;
679 
680  if (npoints0 > 0) {
681  // points0->GetNumberOfPoints() is equal to the upper bound
682  // on the points in the final merged grid. We need to temporarily
683  // set it to the number of points added to the merged grid so far.
684 
685  points0->GetData()->SetNumberOfTuples(npoints0);
686 
687  ptArrays[0] = points0;
688  ptArrays[1] = points1;
689  numArrays = 2;
690  } else {
691  ptArrays[0] = points1;
692  numArrays = 1;
693  }
694 
695  kd->BuildLocatorFromPoints(ptArrays, numArrays);
696 
697  vtkIdTypeArray *pointToEquivClassMap =
698  kd->BuildMapForDuplicatePoints(this->PointMergeTolerance);
699 
700  kd->Delete();
701 
702  if (npoints0 > 0) {
703  points0->GetData()->SetNumberOfTuples(this->TotalNumberOfPoints);
704  }
705 
706  // The map we get back isn't quite what we need. The range of
707  // the map is a subset of original point IDs which each
708  // represent an equivalence class of duplicate points. But the
709  // point chosen to represent the class could be any one of the
710  // equivalent points. We need to create a map that uses IDs
711  // of points in the points0 array as the representative, and
712  // then new logical contiguous point IDs
713  // (npoints0, npoints0+1, ..., numUniquePoints-1) for the
714  // points in the new set that are not duplicates of points
715  // in the points0 array.
716 
717  std::map<vtkIdType, vtkIdType> newIdMap;
718 
719  if (npoints0 > 0) // these were already a unique set
720  {
721  for (ptId = 0; ptId < npoints0; ptId++) {
722  vtkIdType EqClassRep = pointToEquivClassMap->GetValue(ptId);
723 
724  if (EqClassRep != ptId) {
725  newIdMap.insert(
726  std::map<vtkIdType, vtkIdType>::value_type(EqClassRep, ptId));
727  }
728  }
729  }
730  for (ptId = 0; ptId < npoints1; ptId++) {
731  vtkIdType EqClassRep = pointToEquivClassMap->GetValue(ptId + npoints0);
732 
733  if (EqClassRep < npoints0) {
734  idMap->SetValue(ptId,
735  EqClassRep); // a duplicate of a point in the first set
736  continue;
737  }
738 
739  std::pair<std::map<vtkIdType, vtkIdType>::iterator, bool> inserted =
740 
741  newIdMap.insert(std::map<vtkIdType, vtkIdType>::value_type(
742  EqClassRep, nextNewLocalId));
743 
744  bool newEqClassRep = inserted.second;
745  vtkIdType existingMappedId = inserted.first->second;
746 
747  if (newEqClassRep) {
748  idMap->SetValue(ptId, nextNewLocalId); // here's a new unique point
749 
750  nextNewLocalId++;
751  } else {
752  idMap->SetValue(
753  ptId, existingMappedId); // a duplicate of a point in the new set
754  }
755  }
756 
757  pointToEquivClassMap->Delete();
758  newIdMap.clear();
759  }
760 
761  if (!ps) {
762  points1->Delete();
763  }
764 
765  return idMap;
766 }
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkIdType NumberOfPoints
Definition: mergeCells.H:120
vtkUnstructuredGrid * UnstructuredGrid
Definition: mergeCells.H:142
vtkIdType TotalNumberOfPoints
Definition: mergeCells.H:117

◆ MergeDataSet()

vtkSmartPointer< vtkIdTypeArray > NEM::MSH::mergeCells::MergeDataSet ( vtkDataSet *  set)

This call returns after the merge has completed. Be sure to call SetTotalNumberOfCells, SetTotalNumberOfPoints, and SetTotalNumberOfDataSets before making this call. Return the idMap from old and new points.

Definition at line 111 of file mergeCells.C.

111  {
112  vtkIdType newPtId, oldPtId, newCellId;
113  vtkSmartPointer<vtkIdTypeArray> idMap;
114 
115  vtkUnstructuredGrid *ugrid = this->UnstructuredGrid;
116 
117  if (!ugrid) {
118  vtkErrorMacro(<< "SetUnstructuredGrid first");
119  return nullptr;
120  }
121 
122  if (this->TotalNumberOfDataSets <= 0) {
123  // TotalNumberOfCells and TotalNumberOfPoints may both be zero
124  // if all data sets to be merged are empty
125 
126  vtkErrorMacro(<< "Must SetTotalNumberOfCells, SetTotalNumberOfPoints and "
127  "SetTotalNumberOfDataSets (upper bounds at least)"
128  " before starting to MergeDataSets");
129 
130  return nullptr;
131  }
132 
133  vtkPointData *pointArrays = set->GetPointData();
134  vtkCellData *cellArrays = set->GetCellData();
135 
136  // Since mergeCells is to be used only on distributed vtkDataSets,
137  // each DataSet should have the same field arrays. However I've been
138  // told that the field arrays may get rearranged in the process of
139  // Marshalling/UnMarshalling. So we use a
140  // vtkDataSetAttributes::FieldList to ensure the field arrays are
141  // merged in the right order.
142 
143  if (ugrid->GetNumberOfCells() == 0) {
144  vtkPointSet *check1 = vtkPointSet::SafeDownCast(set);
145 
146  if (check1) {
147  this->InputIsPointSet = 1;
148  vtkUnstructuredGrid *check2 = vtkUnstructuredGrid::SafeDownCast(set);
149  this->InputIsUGrid = (check2 != nullptr);
150  }
151 
152  this->StartUGrid(set);
153  } else {
154  this->ptList->IntersectFieldList(pointArrays);
155  this->cellList->IntersectFieldList(cellArrays);
156  }
157 
158  vtkIdType numPoints = set->GetNumberOfPoints();
159  vtkIdType numCells = set->GetNumberOfCells();
160 
161  if (numCells == 0) {
162  return nullptr;
163  }
164 
165  if (this->MergeDuplicatePoints) {
166  if (this->UseGlobalIds) // faster by far
167  {
168  // Note: It has been observed that an input dataset may
169  // have an invalid global ID array. Using the array to
170  // merge points results in bad geometry. It may be
171  // worthwhile to do a quick sanity check when merging
172  // points. Downside is that will slow down this filter.
173 
174  idMap = vtkSmartPointer<vtkIdTypeArray>::Take(
175  this->MapPointsToIdsUsingGlobalIds(set));
176  } else {
177  idMap = vtkSmartPointer<vtkIdTypeArray>::Take(
178  this->MapPointsToIdsUsingLocator(set));
179  }
180  } else {
181  idMap = NULL;
182  }
183 
184  vtkIdType nextPt = (vtkIdType)this->NumberOfPoints;
185 
186  vtkPoints *pts = ugrid->GetPoints();
187 
188  for (oldPtId = 0; oldPtId < numPoints; oldPtId++) {
189  if (idMap) {
190  newPtId = idMap->GetValue(oldPtId);
191  } else {
192  newPtId = nextPt;
193  }
194 
195  if (newPtId == nextPt) {
196  pts->SetPoint(nextPt, set->GetPoint(oldPtId));
197 
198  ugrid->GetPointData()->CopyData(*this->ptList, pointArrays,
199  this->nextGrid, oldPtId, nextPt);
200 
201  nextPt++;
202  }
203  }
204 
205  pts->Modified(); // so that subsequent GetBounds will be correct
206 
207  if (this->InputIsUGrid) {
208  newCellId = this->AddNewCellsUnstructuredGrid(set, idMap);
209  } else {
210  newCellId = this->AddNewCellsDataSet(set, idMap);
211  }
212 
213  this->NumberOfPoints = nextPt;
214  this->NumberOfCells = newCellId;
215 
216  this->nextGrid++;
217 
218  return idMap;
219 }
vtkIdType AddNewCellsDataSet(vtkDataSet *set, vtkIdTypeArray *idMap)
Definition: mergeCells.C:221
vtkIdTypeArray * MapPointsToIdsUsingGlobalIds(vtkDataSet *set)
Definition: mergeCells.C:544
vtkDataSetAttributes::FieldList * cellList
Definition: mergeCells.H:140
vtkIdTypeArray * MapPointsToIdsUsingLocator(vtkDataSet *set)
Definition: mergeCells.C:589
vtkIdType NumberOfPoints
Definition: mergeCells.H:120
vtkUnstructuredGrid * UnstructuredGrid
Definition: mergeCells.H:142
vtkDataSetAttributes::FieldList * ptList
Definition: mergeCells.H:139
vtkIdType AddNewCellsUnstructuredGrid(vtkDataSet *set, vtkIdTypeArray *idMap)
Definition: mergeCells.C:285
vtkIdType NumberOfCells
Definition: mergeCells.H:119
void StartUGrid(vtkDataSet *set)
Definition: mergeCells.C:472

◆ New()

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

◆ operator=()

void NEM::MSH::mergeCells::operator= ( const mergeCells )
protecteddelete

◆ StartUGrid()

void NEM::MSH::mergeCells::StartUGrid ( vtkDataSet *  set)
protected

Definition at line 472 of file mergeCells.C.

References NEM::MSH::New().

472  {
473  vtkPointData *PD = set->GetPointData();
474  vtkCellData *CD = set->GetCellData();
475 
476  vtkUnstructuredGrid *ugrid = this->UnstructuredGrid;
477 
478  if (!this->InputIsUGrid) {
479  ugrid->Allocate(this->TotalNumberOfCells);
480  }
481 
482  vtkPoints *pts = vtkPoints::New();
483 
484  // If the input has a vtkPoints object, we'll make the merged output
485  // grid have a vtkPoints object of the same data type. Otherwise,
486  // the merged output grid will have the default of points of type float.
487 
488  if (this->InputIsPointSet) {
489  vtkPointSet *ps = vtkPointSet::SafeDownCast(set);
490  pts->SetDataType(ps->GetPoints()->GetDataType());
491  }
492 
493  pts->SetNumberOfPoints(
494  this->TotalNumberOfPoints); // allocate for upper bound
495 
496  ugrid->SetPoints(pts);
497 
498  pts->Delete();
499 
500  // Order of field arrays may get changed when data sets are
501  // marshalled/sent/unmarshalled. So we need to re-index the
502  // field arrays before copying them using a FieldList
503 
504  this->ptList =
505  new vtkDataSetAttributes::FieldList(this->TotalNumberOfDataSets);
506  this->cellList =
507  new vtkDataSetAttributes::FieldList(this->TotalNumberOfDataSets);
508 
509  this->ptList->InitializeFieldList(PD);
510  this->cellList->InitializeFieldList(CD);
511 
512  if (this->UseGlobalIds) {
513  ugrid->GetPointData()->CopyGlobalIdsOn();
514  }
515  ugrid->GetPointData()->CopyAllocate(*ptList, this->TotalNumberOfPoints);
516  if (this->UseGlobalCellIds) {
517  ugrid->GetCellData()->CopyGlobalIdsOn();
518  }
519  ugrid->GetCellData()->CopyAllocate(*cellList, this->TotalNumberOfCells);
520 
521  return;
522 }
vtkIdType TotalNumberOfCells
Definition: mergeCells.H:116
vtkDataSetAttributes::FieldList * cellList
Definition: mergeCells.H:140
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
vtkUnstructuredGrid * UnstructuredGrid
Definition: mergeCells.H:142
vtkDataSetAttributes::FieldList * ptList
Definition: mergeCells.H:139
vtkIdType TotalNumberOfPoints
Definition: mergeCells.H:117

◆ vtkBooleanMacro()

NEM::MSH::mergeCells::vtkBooleanMacro ( MergeDuplicatePoints  ,
int   
)

◆ vtkGetMacro() [1/7]

NEM::MSH::mergeCells::vtkGetMacro ( TotalNumberOfCells  ,
vtkIdType   
)

◆ vtkGetMacro() [2/7]

NEM::MSH::mergeCells::vtkGetMacro ( TotalNumberOfPoints  ,
vtkIdType   
)

◆ vtkGetMacro() [3/7]

NEM::MSH::mergeCells::vtkGetMacro ( UseGlobalIds  ,
int   
)

◆ vtkGetMacro() [4/7]

NEM::MSH::mergeCells::vtkGetMacro ( PointMergeTolerance  ,
float   
)

◆ vtkGetMacro() [5/7]

NEM::MSH::mergeCells::vtkGetMacro ( UseGlobalCellIds  ,
int   
)

◆ vtkGetMacro() [6/7]

NEM::MSH::mergeCells::vtkGetMacro ( MergeDuplicatePoints  ,
int   
)

◆ vtkGetMacro() [7/7]

NEM::MSH::mergeCells::vtkGetMacro ( TotalNumberOfDataSets  ,
int   
)

◆ vtkGetObjectMacro()

NEM::MSH::mergeCells::vtkGetObjectMacro ( UnstructuredGrid  ,
vtkUnstructuredGrid   
)

◆ vtkSetClampMacro()

NEM::MSH::mergeCells::vtkSetClampMacro ( PointMergeTolerance  ,
float  ,
0.  0f,
VTK_FLOAT_MAX   
)

◆ vtkSetMacro() [1/6]

NEM::MSH::mergeCells::vtkSetMacro ( TotalNumberOfCells  ,
vtkIdType   
)

◆ vtkSetMacro() [2/6]

NEM::MSH::mergeCells::vtkSetMacro ( TotalNumberOfPoints  ,
vtkIdType   
)

◆ vtkSetMacro() [3/6]

NEM::MSH::mergeCells::vtkSetMacro ( UseGlobalIds  ,
int   
)

◆ vtkSetMacro() [4/6]

NEM::MSH::mergeCells::vtkSetMacro ( UseGlobalCellIds  ,
int   
)

◆ vtkSetMacro() [5/6]

NEM::MSH::mergeCells::vtkSetMacro ( MergeDuplicatePoints  ,
int   
)

◆ vtkSetMacro() [6/6]

NEM::MSH::mergeCells::vtkSetMacro ( TotalNumberOfDataSets  ,
int   
)

◆ vtkSetObjectMacro()

NEM::MSH::mergeCells::vtkSetObjectMacro ( UnstructuredGrid  ,
vtkUnstructuredGrid   
)

◆ vtkTypeMacro()

NEM::MSH::mergeCells::vtkTypeMacro ( mergeCells  ,
vtkObject   
) &

Member Data Documentation

◆ cellList

vtkDataSetAttributes::FieldList* NEM::MSH::mergeCells::cellList
protected

Definition at line 140 of file mergeCells.H.

◆ GlobalCellIdArray

void* NEM::MSH::mergeCells::GlobalCellIdArray
protected

Definition at line 128 of file mergeCells.H.

◆ GlobalCellIdArrayType

int NEM::MSH::mergeCells::GlobalCellIdArrayType
protected

Definition at line 127 of file mergeCells.H.

◆ GlobalCellIdMap

mergeCellsSTLCloak* NEM::MSH::mergeCells::GlobalCellIdMap
protected

Definition at line 137 of file mergeCells.H.

◆ GlobalIdArray

void* NEM::MSH::mergeCells::GlobalIdArray
protected

Definition at line 124 of file mergeCells.H.

◆ GlobalIdArrayType

int NEM::MSH::mergeCells::GlobalIdArrayType
protected

Definition at line 123 of file mergeCells.H.

◆ GlobalIdMap

mergeCellsSTLCloak* NEM::MSH::mergeCells::GlobalIdMap
protected

Definition at line 136 of file mergeCells.H.

◆ InputIsPointSet

char NEM::MSH::mergeCells::InputIsPointSet
protected

Definition at line 134 of file mergeCells.H.

◆ InputIsUGrid

char NEM::MSH::mergeCells::InputIsUGrid
protected

Definition at line 133 of file mergeCells.H.

◆ MergeDuplicatePoints

int NEM::MSH::mergeCells::MergeDuplicatePoints
protected

Definition at line 131 of file mergeCells.H.

◆ nextGrid

int NEM::MSH::mergeCells::nextGrid
protected

Definition at line 144 of file mergeCells.H.

◆ NumberOfCells

vtkIdType NEM::MSH::mergeCells::NumberOfCells
protected

Definition at line 119 of file mergeCells.H.

◆ NumberOfPoints

vtkIdType NEM::MSH::mergeCells::NumberOfPoints
protected

Definition at line 120 of file mergeCells.H.

◆ override

vtkIndent indent NEM::MSH::mergeCells::override

Definition at line 58 of file mergeCells.H.

◆ PointMergeTolerance

float NEM::MSH::mergeCells::PointMergeTolerance
protected

Definition at line 130 of file mergeCells.H.

◆ ptList

vtkDataSetAttributes::FieldList* NEM::MSH::mergeCells::ptList
protected

Definition at line 139 of file mergeCells.H.

◆ TotalNumberOfCells

vtkIdType NEM::MSH::mergeCells::TotalNumberOfCells
protected

Definition at line 116 of file mergeCells.H.

◆ TotalNumberOfDataSets

int NEM::MSH::mergeCells::TotalNumberOfDataSets
protected

Definition at line 114 of file mergeCells.H.

◆ TotalNumberOfPoints

vtkIdType NEM::MSH::mergeCells::TotalNumberOfPoints
protected

Definition at line 117 of file mergeCells.H.

◆ UnstructuredGrid

vtkUnstructuredGrid* NEM::MSH::mergeCells::UnstructuredGrid
protected

Definition at line 142 of file mergeCells.H.

◆ UseGlobalCellIds

int NEM::MSH::mergeCells::UseGlobalCellIds
protected

Definition at line 126 of file mergeCells.H.

◆ UseGlobalIds

int NEM::MSH::mergeCells::UseGlobalIds
protected

Definition at line 122 of file mergeCells.H.


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