VTK  9.3.1
vtkAMRCutPlane.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
15 #ifndef vtkAMRCutPlane_h
16 #define vtkAMRCutPlane_h
17 
18 #include "vtkFiltersAMRModule.h" // For export macro
20 
21 #include <map> // For STL map
22 #include <vector> // For STL vector
23 
24 VTK_ABI_NAMESPACE_BEGIN
26 class vtkOverlappingAMR;
28 class vtkInformation;
30 class vtkIndent;
31 class vtkPlane;
32 class vtkUniformGrid;
33 class vtkCell;
34 class vtkPoints;
35 class vtkCellArray;
36 class vtkPointData;
37 class vtkCellData;
38 
39 class VTKFILTERSAMR_EXPORT vtkAMRCutPlane : public vtkMultiBlockDataSetAlgorithm
40 {
41 public:
42  static vtkAMRCutPlane* New();
44  void PrintSelf(ostream& oss, vtkIndent indent) override;
45 
47 
50  vtkSetVector3Macro(Center, double);
52 
54 
57  vtkSetVector3Macro(Normal, double);
59 
61 
64  vtkSetMacro(LevelOfResolution, int);
65  vtkGetMacro(LevelOfResolution, int);
67 
69 
74  vtkSetMacro(UseNativeCutter, bool);
75  vtkGetMacro(UseNativeCutter, bool);
76  vtkBooleanMacro(UseNativeCutter, bool);
78 
80 
84  virtual void SetController(vtkMultiProcessController*);
85  vtkGetObjectMacro(Controller, vtkMultiProcessController);
87 
88  // Standard pipeline routines
89 
93 
99  vtkInformationVector* outputVector) override;
100 
105 
109  vtkSetMacro(InitialRequest, bool);
110 
111 protected:
112  vtkAMRCutPlane();
113  ~vtkAMRCutPlane() override;
114 
119  vtkPlane* GetCutPlane(vtkOverlappingAMR* metadata);
120 
124  void ExtractCellFromGrid(vtkUniformGrid* grid, vtkCell* cell,
125  std::map<vtkIdType, vtkIdType>& gridPntMapping, vtkPoints* nodes, vtkCellArray* cells);
126 
131  void ExtractPointDataFromGrid(vtkUniformGrid* grid,
132  std::map<vtkIdType, vtkIdType>& gridPntMapping, vtkIdType NumNodes, vtkPointData* PD);
133 
138  void ExtractCellDataFromGrid(
139  vtkUniformGrid* grid, std::vector<vtkIdType>& cellIdxList, vtkCellData* CD);
140 
147  void ComputeAMRBlocksToLoad(vtkPlane* p, vtkOverlappingAMR* m);
148 
149  // Descriription:
150  // Initializes the cut-plane center given the min/max bounds.
151  void InitializeCenter(double min[3], double max[3]);
152 
154 
157  bool PlaneIntersectsAMRBox(vtkPlane* pl, double bounds[6]);
158  bool PlaneIntersectsAMRBox(double plane[4], double bounds[6]);
160 
164  bool PlaneIntersectsCell(vtkPlane* pl, vtkCell* cell);
165 
169  bool IsAMRData2D(vtkOverlappingAMR* input);
170 
174  void CutAMRBlock(
175  vtkPlane* cutPlane, unsigned int blockIdx, vtkUniformGrid* grid, vtkMultiBlockDataSet* dataSet);
176 
178  double Center[3];
179  double Normal[3];
183 
184  std::vector<int> BlocksToLoad;
185 
186 private:
187  vtkAMRCutPlane(const vtkAMRCutPlane&) = delete;
188  void operator=(const vtkAMRCutPlane&) = delete;
189 };
190 
191 VTK_ABI_NAMESPACE_END
192 #endif /* vtkAMRCutPlane_h */
vtkMultiProcessController * Controller
represent and manipulate point attribute data
Definition: vtkPointData.h:29
Store vtkAlgorithm input/output information.
represent and manipulate cell attribute data
Definition: vtkCellData.h:30
A concrete instance of vtkMultiBlockDataSet that provides functionality for cutting an AMR dataset (a...
int vtkIdType
Definition: vtkType.h:315
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
static vtkMultiBlockDataSetAlgorithm * New()
virtual int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
abstract class to specify cell behavior
Definition: vtkCell.h:49
a simple class to control print indentation
Definition: vtkIndent.h:28
perform various plane computations
Definition: vtkPlane.h:25
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
image data with blanking
object to represent cell connectivity
Definition: vtkCellArray.h:175
Composite dataset that organizes datasets into blocks.
hierarchical dataset of vtkUniformGrids
Store zero or more vtkInformation instances.
#define max(a, b)
represent and manipulate 3D points
Definition: vtkPoints.h:28
std::vector< int > BlocksToLoad
virtual int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
Multiprocessing communication superclass.