VTK  9.4.20241221
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 "vtkDeprecation.h" // For VTK_DEPRECATED
19#include "vtkFiltersAMRModule.h" // For export macro
21
22#include <map> // For STL map
23#include <vector> // For STL vector
24
25VTK_ABI_NAMESPACE_BEGIN
29class vtkInformation;
31class vtkIndent;
32class vtkPlane;
33class vtkUniformGrid;
35class vtkCell;
36class vtkPoints;
37class vtkCellArray;
38class vtkPointData;
39class vtkCellData;
40
41class VTKFILTERSAMR_EXPORT vtkAMRCutPlane : public vtkMultiBlockDataSetAlgorithm
42{
43public:
46 void PrintSelf(ostream& oss, vtkIndent indent) override;
47
49
52 vtkSetVector3Macro(Center, double);
54
56
59 vtkSetVector3Macro(Normal, double);
61
63
66 vtkSetMacro(LevelOfResolution, int);
67 vtkGetMacro(LevelOfResolution, int);
69
71
76 vtkSetMacro(UseNativeCutter, bool);
77 vtkGetMacro(UseNativeCutter, bool);
78 vtkBooleanMacro(UseNativeCutter, bool);
80
82
87 vtkGetObjectMacro(Controller, vtkMultiProcessController);
89
90 // Standard pipeline routines
91
93 int FillInputPortInformation(int port, vtkInformation* info) override;
94 int FillOutputPortInformation(int port, vtkInformation* info) override;
95
101 vtkInformationVector* outputVector) override;
102
107
111 vtkSetMacro(InitialRequest, bool);
112
113protected:
115 ~vtkAMRCutPlane() override;
116
122
127 std::map<vtkIdType, vtkIdType>& gridPntMapping, vtkPoints* nodes, vtkCellArray* cells);
128
134 std::map<vtkIdType, vtkIdType>& gridPntMapping, vtkIdType NumNodes, vtkPointData* PD);
135
141 vtkUniformGrid* grid, std::vector<vtkIdType>& cellIdxList, vtkCellData* CD);
142
150
151 // Descriription:
152 // Initializes the cut-plane center given the min/max bounds.
153 void InitializeCenter(double min[3], double max[3]);
154
156
159 bool PlaneIntersectsAMRBox(vtkPlane* pl, double bounds[6]);
160 bool PlaneIntersectsAMRBox(double plane[4], double bounds[6]);
162
167
172
177
181 VTK_DEPRECATED_IN_9_4_0("Use CutAMRBlock(vtkPlane*, vtkUniformGrid*) instead.")
182 void CutAMRBlock(
183 vtkPlane* cutPlane, unsigned int blockIdx, vtkUniformGrid* grid, vtkMultiBlockDataSet* dataSet);
184
185 int LevelOfResolution;
186 double Center[3];
187 double Normal[3];
188 bool InitialRequest;
189 bool UseNativeCutter;
191
192 std::vector<int> BlocksToLoad;
193
194private:
195 vtkAMRCutPlane(const vtkAMRCutPlane&) = delete;
196 void operator=(const vtkAMRCutPlane&) = delete;
197};
198
199VTK_ABI_NAMESPACE_END
200#endif /* vtkAMRCutPlane_h */
A concrete instance of vtkMultiBlockDataSet that provides functionality for cutting an AMR dataset (a...
bool IsAMRData2D(vtkOverlappingAMR *input)
A utility function that checks if the input AMR data is 2-D.
vtkPlane * GetCutPlane(vtkOverlappingAMR *metadata)
Returns the cut-plane defined by a vtkCutPlane instance based on the user-supplied center and normal.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void ExtractCellFromGrid(vtkUniformGrid *grid, vtkCell *cell, std::map< vtkIdType, vtkIdType > &gridPntMapping, vtkPoints *nodes, vtkCellArray *cells)
Extracts cell.
void InitializeCenter(double min[3], double max[3])
virtual void SetController(vtkMultiProcessController *)
Set/Get a multiprocess controller for parallel processing.
void ComputeAMRBlocksToLoad(vtkPlane *p, vtkOverlappingAMR *m)
Given a cut-plane, p, and the metadata, m, this method computes which blocks need to be loaded.
static vtkAMRCutPlane * New()
~vtkAMRCutPlane() override
bool PlaneIntersectsAMRBox(double plane[4], double bounds[6])
Determines if a plane intersects with an AMR box.
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
Performs upstream requests to the reader.
vtkSmartPointer< vtkUnstructuredGrid > CutAMRBlock(vtkPlane *cutPlane, vtkUniformGrid *grid)
Applies cutting to an AMR block.
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
int RequestInformation(vtkInformation *rqst, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Gets the metadata from upstream module and determines which blocks should be loaded by this instance.
void ExtractPointDataFromGrid(vtkUniformGrid *grid, std::map< vtkIdType, vtkIdType > &gridPntMapping, vtkIdType NumNodes, vtkPointData *PD)
Given the grid and a subset ID pair, grid IDs mapping to the extracted grid IDs, extract the point da...
void PrintSelf(ostream &oss, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
bool PlaneIntersectsAMRBox(vtkPlane *pl, double bounds[6])
Determines if a plane intersects with an AMR box.
bool PlaneIntersectsCell(vtkPlane *pl, vtkCell *cell)
Determines if a plane intersects with a grid cell.
void ExtractCellDataFromGrid(vtkUniformGrid *grid, std::vector< vtkIdType > &cellIdxList, vtkCellData *CD)
Given the grid and the list of cells that are extracted, extract the corresponding cell data.
object to represent cell connectivity
represent and manipulate cell attribute data
abstract class to specify cell behavior
Definition vtkCell.h:130
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
Composite dataset that organizes datasets into blocks.
Multiprocessing communication superclass.
hierarchical dataset of vtkUniformGrids
perform various plane computations
Definition vtkPlane.h:138
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:139
Hold a reference to a vtkObjectBase instance.
image data with blanking
dataset represents arbitrary combinations of all possible cell types
#define VTK_DEPRECATED_IN_9_4_0(reason)
int vtkIdType
Definition vtkType.h:315
#define max(a, b)