VTK
vtkDistributedDataFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkDistributedDataFilter.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
15 /*----------------------------------------------------------------------------
16  Copyright (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
19 
89 #ifndef vtkDistributedDataFilter_h
90 #define vtkDistributedDataFilter_h
91 
92 #include "vtkFiltersParallelMPIModule.h" // For export macro
93 #include "vtkDataObjectAlgorithm.h"
94 
95 class vtkBSPCuts;
96 class vtkDataArray;
97 class vtkDistributedDataFilterSTLCloak;
98 class vtkFloatArray;
99 class vtkIdList;
100 class vtkIdTypeArray;
101 class vtkIntArray;
103 class vtkPKdTree;
104 class vtkUnstructuredGrid;
105 
106 class VTKFILTERSPARALLELMPI_EXPORT vtkDistributedDataFilter: public vtkDataObjectAlgorithm
107 {
108  vtkTypeMacro(vtkDistributedDataFilter,
110 
111 public:
112  void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
113 
114  static vtkDistributedDataFilter *New();
115 
117 
120  void SetController(vtkMultiProcessController *c);
121  vtkGetObjectMacro(Controller, vtkMultiProcessController);
123 
134  vtkPKdTree *GetKdtree();
135 
153  vtkBooleanMacro(RetainKdtree, int);
154  vtkGetMacro(RetainKdtree, int);
155  vtkSetMacro(RetainKdtree, int);
156 
168  vtkBooleanMacro(IncludeAllIntersectingCells, int);
169  vtkGetMacro(IncludeAllIntersectingCells, int);
170  vtkSetMacro(IncludeAllIntersectingCells, int);
171 
178  vtkBooleanMacro(ClipCells, int);
179  vtkGetMacro(ClipCells, int);
180  vtkSetMacro(ClipCells, int);
181 
183  ASSIGN_TO_ONE_REGION=0,
184  ASSIGN_TO_ALL_INTERSECTING_REGIONS=1,
185  SPLIT_BOUNDARY_CELLS=2
186  };
187 
189 
192  void SetBoundaryMode(int mode);
194  { this->SetBoundaryMode(vtkDistributedDataFilter::ASSIGN_TO_ONE_REGION); }
196  { this->SetBoundaryMode(
198  }
200  { this->SetBoundaryMode(vtkDistributedDataFilter::SPLIT_BOUNDARY_CELLS); }
201  int GetBoundaryMode();
203 
208 
219  vtkBooleanMacro(UseMinimalMemory, int);
220  vtkGetMacro(UseMinimalMemory, int);
221  vtkSetMacro(UseMinimalMemory, int);
222 
223 
228  vtkBooleanMacro(Timing, int);
229  vtkSetMacro(Timing, int);
230  vtkGetMacro(Timing, int);
231 
244  vtkBSPCuts* GetCuts() {return this->UserCuts;}
245  void SetCuts(vtkBSPCuts* cuts);
246 
256  void SetUserRegionAssignments(const int *map, int numRegions);
257 
258 protected:
261 
269  void AssignBoundaryCellsToOneRegionOn();
270  void AssignBoundaryCellsToOneRegionOff();
271  void SetAssignBoundaryCellsToOneRegion(int val);
272 
282  void AssignBoundaryCellsToAllIntersectingRegionsOn();
283  void AssignBoundaryCellsToAllIntersectingRegionsOff();
284  void SetAssignBoundaryCellsToAllIntersectingRegions(int val);
285 
294  void DivideBoundaryCellsOn();
295  void DivideBoundaryCellsOff();
296  void SetDivideBoundaryCells(int val);
297 
305  vtkInformationVector *) VTK_OVERRIDE;
306  void SingleProcessExecute(vtkDataSet *input, vtkUnstructuredGrid *output);
308  vtkInformationVector *) VTK_OVERRIDE;
309  virtual int FillInputPortInformation(int port, vtkInformation *info) VTK_OVERRIDE;
310 
316  virtual int RequestDataObject(vtkInformation*,
318  vtkInformationVector*) VTK_OVERRIDE;
319 
323  int RequestDataInternal(vtkDataSet* input, vtkUnstructuredGrid* output);
324 
325 private:
326 
327  enum{
328  DeleteNo = 0,
329  DeleteYes = 1
330  };
331 
332  enum{
333  DuplicateCellsNo = 0,
334  DuplicateCellsYes = 1
335  };
336 
337  enum{
338  GhostCellsNo = 0,
339  GhostCellsYes = 1
340  };
341 
342  enum{
343  UnsetGhostLevel = 99
344  };
345 
349  int PartitionDataAndAssignToProcesses(vtkDataSet *set);
350 
354  vtkUnstructuredGrid *RedistributeDataSet(vtkDataSet *set, vtkDataSet *input);
355 
359  int ClipGridCells(vtkUnstructuredGrid *grid);
360 
364  vtkUnstructuredGrid * AcquireGhostCells(vtkUnstructuredGrid *grid);
365 
369  void ComputeMyRegionBounds();
370 
374  int CheckFieldArrayTypes(vtkDataSet *set);
375 
381  vtkDataSet *TestFixTooFewInputFiles(vtkDataSet *input);
382 
386  vtkUnstructuredGrid *MPIRedistribute(vtkDataSet *in, vtkDataSet *input);
387 
391  vtkIdList **GetCellIdsForProcess(int proc, int *nlists);
392 
397  void SetUpPairWiseExchange();
398 
400 
403  void FreeIntArrays(vtkIdTypeArray **ar);
404  static void FreeIdLists(vtkIdList**lists, int nlists);
405  static vtkIdType GetIdListSize(vtkIdList**lists, int nlists);
407 
409 
412  vtkIdTypeArray *ExchangeCounts(vtkIdType myCount, int tag);
413  vtkIdTypeArray *ExchangeCountsLean(vtkIdType myCount, int tag);
414  vtkIdTypeArray *ExchangeCountsFast(vtkIdType myCount, int tag);
416 
418 
421  vtkIdTypeArray **ExchangeIdArrays(vtkIdTypeArray **arIn,
422  int deleteSendArrays, int tag);
423  vtkIdTypeArray **ExchangeIdArraysLean(vtkIdTypeArray **arIn,
424  int deleteSendArrays, int tag);
425  vtkIdTypeArray **ExchangeIdArraysFast(vtkIdTypeArray **arIn,
426  int deleteSendArrays, int tag);
428 
430 
433  vtkFloatArray **ExchangeFloatArrays(vtkFloatArray **myArray,
434  int deleteSendArrays, int tag);
435  vtkFloatArray **ExchangeFloatArraysLean(vtkFloatArray **myArray,
436  int deleteSendArrays, int tag);
437  vtkFloatArray **ExchangeFloatArraysFast(vtkFloatArray **myArray,
438  int deleteSendArrays, int tag);
440 
442 
445  vtkUnstructuredGrid *ExchangeMergeSubGrids(vtkIdList **cellIds, int deleteCellIds,
446  vtkDataSet *myGrid, int deleteMyGrid,
447  int filterOutDuplicateCells, int ghostCellFlag, int tag);
448  vtkUnstructuredGrid *ExchangeMergeSubGrids(vtkIdList ***cellIds, int *numLists,
449  int deleteCellIds,
450  vtkDataSet *myGrid, int deleteMyGrid,
451  int filterOutDuplicateCells, int ghostCellFlag, int tag);
452  vtkUnstructuredGrid *ExchangeMergeSubGridsLean(
453  vtkIdList ***cellIds, int *numLists,
454  int deleteCellIds,
455  vtkDataSet *myGrid, int deleteMyGrid,
456  int filterOutDuplicateCells, int ghostCellFlag, int tag);
457  vtkUnstructuredGrid *ExchangeMergeSubGridsFast(
458  vtkIdList ***cellIds, int *numLists,
459  int deleteCellIds,
460  vtkDataSet *myGrid, int deleteMyGrid,
461  int filterOutDuplicateCells, int ghostCellFlag, int tag);
463 
464 
466 
469  char *MarshallDataSet(vtkUnstructuredGrid *extractedGrid, int &size);
470  vtkUnstructuredGrid *UnMarshallDataSet(char *buf, int size);
472 
474 
477  void ClipCellsToSpatialRegion(vtkUnstructuredGrid *grid);
478 #if 0
479  void ClipWithVtkClipDataSet(vtkUnstructuredGrid *grid, double *bounds,
480  vtkUnstructuredGrid **outside, vtkUnstructuredGrid **inside);
481 #endif
482 
483 
484  void ClipWithBoxClipDataSet(vtkUnstructuredGrid *grid, double *bounds,
485  vtkUnstructuredGrid **outside, vtkUnstructuredGrid **inside);
486 
488 
494  vtkIdTypeArray *GetGlobalNodeIdArray(vtkDataSet *set);
495  vtkIdType *GetGlobalNodeIds(vtkDataSet *set);
496  vtkIdTypeArray *GetGlobalElementIdArray(vtkDataSet *set);
497  vtkIdType *GetGlobalElementIds(vtkDataSet *set);
498  int AssignGlobalNodeIds(vtkUnstructuredGrid *grid);
499  int AssignGlobalElementIds(vtkDataSet *in);
500  vtkIdTypeArray **FindGlobalPointIds(vtkFloatArray **ptarray,
501  vtkIdTypeArray *ids, vtkUnstructuredGrid *grid, vtkIdType &numUniqueMissingPoints);
503 
507  vtkIdTypeArray **MakeProcessLists(vtkIdTypeArray **pointIds,
508  vtkDistributedDataFilterSTLCloak *procs);
509 
513  vtkIdList **BuildRequestedGrids( vtkIdTypeArray **globalPtIds,
514  vtkUnstructuredGrid *grid,
515  vtkDistributedDataFilterSTLCloak *ptIdMap);
516 
518 
521  int InMySpatialRegion(float x, float y, float z);
522  int InMySpatialRegion(double x, double y, double z);
523  int StrictlyInsideMyBounds(float x, float y, float z);
524  int StrictlyInsideMyBounds(double x, double y, double z);
526 
528 
531  vtkIdTypeArray **GetGhostPointIds(int ghostLevel, vtkUnstructuredGrid *grid,
532  int AddCellsIAlreadyHave);
533  vtkUnstructuredGrid *AddGhostCellsUniqueCellAssignment(
534  vtkUnstructuredGrid *myGrid,
535  vtkDistributedDataFilterSTLCloak *globalToLocalMap);
536  vtkUnstructuredGrid *AddGhostCellsDuplicateCellAssignment(
537  vtkUnstructuredGrid *myGrid,
538  vtkDistributedDataFilterSTLCloak *globalToLocalMap);
539  vtkUnstructuredGrid *SetMergeGhostGrid(
540  vtkUnstructuredGrid *ghostCellGrid,
541  vtkUnstructuredGrid *incomingGhostCells,
542  int ghostLevel, vtkDistributedDataFilterSTLCloak *idMap);
544 
546 
549  vtkUnstructuredGrid *ExtractCells(vtkIdList *list,
550  int deleteCellLists, vtkDataSet *in);
551  vtkUnstructuredGrid *ExtractCells(vtkIdList **lists, int nlists,
552  int deleteCellLists, vtkDataSet *in);
553  vtkUnstructuredGrid *ExtractZeroCellGrid(vtkDataSet *in);
555 
557 
560  static int GlobalPointIdIsUsed(vtkUnstructuredGrid *grid,
561  int ptId, vtkDistributedDataFilterSTLCloak *globalToLocal);
562  static int LocalPointIdIsUsed(vtkUnstructuredGrid *grid, int ptId);
563  static vtkIdType FindId(vtkIdTypeArray *ids, vtkIdType gid, vtkIdType startLoc);
565 
569  static vtkIdTypeArray *AddPointAndCells(vtkIdType gid,
570  vtkIdType localId,
571  vtkUnstructuredGrid *grid,
572  vtkIdType *gidCells,
573  vtkIdTypeArray *ids);
574 
576 
579  static void AddConstantUnsignedCharPointArray(vtkUnstructuredGrid *grid,
580  const char *arrayName, unsigned char val);
581  static void AddConstantUnsignedCharCellArray(vtkUnstructuredGrid *grid,
582  const char *arrayName, unsigned char val);
584 
588  static void RemoveRemoteCellsFromList(vtkIdList *cellList,
589  vtkIdType *gidCells,
590  vtkIdType *remoteCells,
591  vtkIdType nRemoteCells);
592 
596  static vtkUnstructuredGrid *MergeGrids(vtkDataSet **sets, int nsets,
597  int deleteDataSets,
598  int useGlobalNodeIds, float pointMergeTolerance,
599  int useGlobalCellIds);
600 
601  vtkPKdTree *Kdtree;
602  vtkMultiProcessController *Controller;
603 
604  int NumProcesses;
605  int MyId;
606 
607  int *Target;
608  int *Source;
609 
610  int NumConvexSubRegions;
611  double *ConvexSubRegionBounds;
612 
613  int GhostLevel;
614 
615  int RetainKdtree;
616  int IncludeAllIntersectingCells;
617  int ClipCells;
618  int AssignBoundaryCellsToOneRegion;
619  int AssignBoundaryCellsToAllIntersectingRegions;
620  int DivideBoundaryCells;
621 
622  int Timing;
623 
624  int NextProgressStep;
625  double ProgressIncrement;
626 
627  int UseMinimalMemory;
628 
629  vtkBSPCuts* UserCuts;
630 
631  vtkDistributedDataFilter(const vtkDistributedDataFilter&) VTK_DELETE_FUNCTION;
632  void operator=(const vtkDistributedDataFilter&) VTK_DELETE_FUNCTION;
633 
634  class vtkInternals;
635  vtkInternals* Internals;
636 
637 };
638 #endif
virtual int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
Build a k-d tree decomposition of a list of points.
Definition: vtkPKdTree.h:55
Store vtkAlgorithm input/output information.
This class represents an axis-aligned Binary Spatial Partitioning of a 3D space.
Definition: vtkBSPCuts.h:44
abstract class to specify dataset behavior
Definition: vtkDataSet.h:56
static vtkDataObjectAlgorithm * New()
dynamic, self-adjusting array of float
Definition: vtkFloatArray.h:35
void SetBoundaryModeToSplitBoundaryCells()
Handling of ClipCells and IncludeAllIntersectingCells.
dynamic, self-adjusting array of vtkIdType
int vtkIdType
Definition: vtkType.h:345
void SetBoundaryModeToAssignToAllIntersectingRegions()
Handling of ClipCells and IncludeAllIntersectingCells.
vtkBSPCuts * GetCuts()
You can set the k-d tree decomposition, rather than have D3 compute it.
void SetBoundaryModeToAssignToOneRegion()
Handling of ClipCells and IncludeAllIntersectingCells.
Distribute data among processors.
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:39
a simple class to control print indentation
Definition: vtkIndent.h:33
list of point or cell ids
Definition: vtkIdList.h:30
dataset represents arbitrary combinations of all possible cell types
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
virtual int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
vtkSetMacro(IgnoreDriverBugs, bool)
When set known driver bugs are ignored during driver feature detection.
Superclass for algorithms that produce only data object as output.
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
virtual int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
Store zero or more vtkInformation instances.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkBooleanMacro(IgnoreDriverBugs, bool)
When set known driver bugs are ignored during driver feature detection.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
Multiprocessing communication superclass.