VTK  9.1.0
vtkRedistributeDataSetFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkRedistributeDataSetFilter.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 =========================================================================*/
70 #ifndef vtkRedistributeDataSetFilter_h
71 #define vtkRedistributeDataSetFilter_h
72 
73 #include "vtkDataObjectAlgorithm.h"
74 #include "vtkFiltersParallelDIY2Module.h" // for export macros
75 #include "vtkSmartPointer.h" // for vtkSmartPointer
76 
77 #include <memory> // for std::shared_ptr
78 #include <vector> // for std::vector
79 
80 // clang-format off
81 #include "vtk_diy2.h" // for DIY2 APIs
82 #include VTK_DIY2(diy/assigner.hpp)
83 // clang-format on
84 
86 class vtkBoundingBox;
90 class vtkDataObjectTree;
91 
92 class VTKFILTERSPARALLELDIY2_EXPORT vtkRedistributeDataSetFilter : public vtkDataObjectAlgorithm
93 {
94 public:
97  void PrintSelf(ostream& os, vtkIndent indent) override;
98 
100 
105  vtkGetObjectMacro(Controller, vtkMultiProcessController);
107 
109  {
110  ASSIGN_TO_ONE_REGION = 0,
111  ASSIGN_TO_ALL_INTERSECTING_REGIONS = 1,
112  SPLIT_BOUNDARY_CELLS = 2
113  };
114 
116 
128  vtkSetClampMacro(BoundaryMode, int, ASSIGN_TO_ONE_REGION, SPLIT_BOUNDARY_CELLS);
129  vtkGetMacro(BoundaryMode, int);
130  void SetBoundaryModeToAssignToOneRegion() { this->SetBoundaryMode(ASSIGN_TO_ONE_REGION); }
132  {
133  this->SetBoundaryMode(ASSIGN_TO_ALL_INTERSECTING_REGIONS);
134  }
135  void SetBoundaryModeToSplitBoundaryCells() { this->SetBoundaryMode(SPLIT_BOUNDARY_CELLS); }
137 
139 
144  vtkSetMacro(UseExplicitCuts, bool);
145  vtkGetMacro(UseExplicitCuts, bool);
146  vtkBooleanMacro(UseExplicitCuts, bool);
148 
150 
153  void SetExplicitCuts(const std::vector<vtkBoundingBox>& boxes);
154  const std::vector<vtkBoundingBox>& GetExplicitCuts() const { return this->ExplicitCuts; }
156  void AddExplicitCut(const vtkBoundingBox& bbox);
157  void AddExplicitCut(const double bbox[6]);
161 
163 
168  void SetAssigner(std::shared_ptr<diy::Assigner> assigner);
169  std::shared_ptr<diy::Assigner> GetAssigner();
170  std::shared_ptr<const diy::Assigner> GetAssigner() const;
171 
173 
183  vtkSetMacro(ExpandExplicitCuts, bool);
184  vtkGetMacro(ExpandExplicitCuts, bool);
185  vtkBooleanMacro(ExpandExplicitCuts, bool);
187 
189 
193  const std::vector<vtkBoundingBox>& GetCuts() const { return this->Cuts; }
194 
196 
214  vtkSetClampMacro(NumberOfPartitions, int, 0, VTK_INT_MAX);
215  vtkGetMacro(NumberOfPartitions, int);
217 
219 
232  vtkSetMacro(PreservePartitionsInOutput, bool);
233  vtkGetMacro(PreservePartitionsInOutput, bool);
234  vtkBooleanMacro(PreservePartitionsInOutput, bool);
236 
238 
242  vtkSetMacro(GenerateGlobalCellIds, bool);
243  vtkGetMacro(GenerateGlobalCellIds, bool);
244  vtkBooleanMacro(GenerateGlobalCellIds, bool);
246 
253  std::vector<vtkBoundingBox> ExpandCuts(
254  const std::vector<vtkBoundingBox>& cuts, const vtkBoundingBox& bounds);
255 
257 
264  vtkSetMacro(EnableDebugging, bool);
265  vtkGetMacro(EnableDebugging, bool);
266  vtkBooleanMacro(EnableDebugging, bool);
268 
270 
278  vtkSetMacro(LoadBalanceAcrossAllBlocks, bool);
279  vtkGetMacro(LoadBalanceAcrossAllBlocks, bool);
280  vtkBooleanMacro(LoadBalanceAcrossAllBlocks, bool);
282 
283 protected:
286 
290 
300  virtual std::vector<vtkBoundingBox> GenerateCuts(vtkDataObject* data);
301 
313  vtkDataSet* dataset, const std::vector<vtkBoundingBox>& cuts);
314 
315 private:
317  void operator=(const vtkRedistributeDataSetFilter&) = delete;
318 
319  bool InitializeCuts(vtkDataObjectTree* input);
320  bool Redistribute(vtkPartitionedDataSet* inputDO, vtkPartitionedDataSet* outputPDS,
321  const std::vector<vtkBoundingBox>& cuts, vtkIdType* mb_offset = nullptr);
322  bool RedistributeDataSet(
323  vtkDataSet* inputDS, vtkPartitionedDataSet* outputPDS, const std::vector<vtkBoundingBox>& cuts);
324  vtkSmartPointer<vtkDataSet> ClipDataSet(vtkDataSet* dataset, const vtkBoundingBox& bbox);
325 
326  void MarkGhostCells(vtkPartitionedDataSet* pieces);
327 
328  vtkSmartPointer<vtkPartitionedDataSet> AssignGlobalCellIds(
329  vtkPartitionedDataSet* input, vtkIdType* mb_offset = nullptr);
330  vtkSmartPointer<vtkDataSet> AssignGlobalCellIds(
331  vtkDataSet* input, vtkIdType* mb_offset = nullptr);
332 
333  void MarkValidDimensions(const vtkBoundingBox& gbounds);
334 
335  std::vector<vtkBoundingBox> ExplicitCuts;
336  std::vector<vtkBoundingBox> Cuts;
337  std::shared_ptr<diy::Assigner> Assigner;
338 
339  vtkMultiProcessController* Controller;
340  int BoundaryMode;
341  int NumberOfPartitions;
342  bool PreservePartitionsInOutput;
343  bool GenerateGlobalCellIds;
344  bool UseExplicitCuts;
345  bool ExpandExplicitCuts;
346  bool EnableDebugging;
347  bool ValidDim[3];
348  bool LoadBalanceAcrossAllBlocks;
349 };
350 
351 #endif
Fast, simple class for representing and operating on 3D bounds.
Superclass for algorithms that produce only data object as output.
provides implementation for most abstract methods in the superclass vtkCompositeDataSet
general representation of visualization data
Definition: vtkDataObject.h:60
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Composite dataset that organizes datasets into blocks.
composite dataset to encapsulates pieces of dataset.
Multiprocessing communication superclass.
composite dataset to encapsulates a dataset consisting of partitions.
redistributes input dataset into requested number of partitions
static vtkRedistributeDataSetFilter * New()
const std::vector< vtkBoundingBox > & GetExplicitCuts() const
Specify the cuts to use when UseExplicitCuts is true.
virtual vtkSmartPointer< vtkPartitionedDataSet > SplitDataSet(vtkDataSet *dataset, const std::vector< vtkBoundingBox > &cuts)
This method is called to split a vtkDataSet into multiple datasets by the vector of vtkBoundingBox pa...
virtual std::vector< vtkBoundingBox > GenerateCuts(vtkDataObject *data)
This method is called to generate the partitions for the input dataset.
const vtkBoundingBox & GetExplicitCut(int index) const
Specify the cuts to use when UseExplicitCuts is true.
void SetBoundaryModeToAssignToAllIntersectingRegions()
Specify how cells on the boundaries are handled.
int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void AddExplicitCut(const double bbox[6])
Specify the cuts to use when UseExplicitCuts is true.
std::vector< vtkBoundingBox > ExpandCuts(const std::vector< vtkBoundingBox > &cuts, const vtkBoundingBox &bounds)
Helper function to expand a collection of bounding boxes to include the bounds specified.
int GetNumberOfExplicitCuts() const
Specify the cuts to use when UseExplicitCuts is true.
void SetController(vtkMultiProcessController *)
Get/Set the controller to use.
void SetExplicitCuts(const std::vector< vtkBoundingBox > &boxes)
Specify the cuts to use when UseExplicitCuts is true.
~vtkRedistributeDataSetFilter() override
std::shared_ptr< diy::Assigner > GetAssigner()
Specify the DIY assigner used for distributing cuts.
void RemoveAllExplicitCuts()
Specify the cuts to use when UseExplicitCuts is true.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
void SetBoundaryModeToAssignToOneRegion()
Specify how cells on the boundaries are handled.
std::shared_ptr< const diy::Assigner > GetAssigner() const
Specify the DIY assigner used for distributing cuts.
void AddExplicitCut(const vtkBoundingBox &bbox)
Specify the cuts to use when UseExplicitCuts is true.
void SetAssigner(std::shared_ptr< diy::Assigner > assigner)
Specify the DIY assigner used for distributing cuts.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetBoundaryModeToSplitBoundaryCells()
Specify how cells on the boundaries are handled.
const std::vector< vtkBoundingBox > & GetCuts() const
Returns the cuts used by the most recent RequestData call.
Hold a reference to a vtkObjectBase instance.
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
@ index
Definition: vtkX3D.h:252
@ data
Definition: vtkX3D.h:321
int vtkIdType
Definition: vtkType.h:332
#define VTK_INT_MAX
Definition: vtkType.h:155