VTK  9.1.0
vtkAMRCutPlane.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkAMRCutPlane.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 =========================================================================*/
24#ifndef vtkAMRCutPlane_h
25#define vtkAMRCutPlane_h
26
27#include "vtkFiltersAMRModule.h" // For export macro
29
30#include <map> // For STL map
31#include <vector> // For STL vector
32
36class vtkInformation;
38class vtkIndent;
39class vtkPlane;
40class vtkUniformGrid;
41class vtkCell;
42class vtkPoints;
43class vtkCellArray;
44class vtkPointData;
45class vtkCellData;
46
47class VTKFILTERSAMR_EXPORT vtkAMRCutPlane : public vtkMultiBlockDataSetAlgorithm
48{
49public:
52 void PrintSelf(ostream& oss, vtkIndent indent) override;
53
55
58 vtkSetVector3Macro(Center, double);
60
62
65 vtkSetVector3Macro(Normal, double);
67
69
72 vtkSetMacro(LevelOfResolution, int);
73 vtkGetMacro(LevelOfResolution, int);
75
77
80 vtkSetMacro(UseNativeCutter, bool);
81 vtkGetMacro(UseNativeCutter, bool);
82 vtkBooleanMacro(UseNativeCutter, bool);
84
86
90 vtkSetMacro(Controller, vtkMultiProcessController*);
91 vtkGetMacro(Controller, vtkMultiProcessController*);
93
94 // Standard pipeline routines
95
99
105 vtkInformationVector* outputVector) override;
106
111
112protected:
114 ~vtkAMRCutPlane() override;
115
121
126 std::map<vtkIdType, vtkIdType>& gridPntMapping, vtkPoints* nodes, vtkCellArray* cells);
127
133 std::map<vtkIdType, vtkIdType>& gridPntMapping, vtkIdType NumNodes, vtkPointData* PD);
134
140 vtkUniformGrid* grid, std::vector<vtkIdType>& cellIdxList, vtkCellData* CD);
141
149
150 // Descriription:
151 // Initializes the cut-plane center given the min/max bounds.
152 void InitializeCenter(double min[3], double max[3]);
153
155
158 bool PlaneIntersectsAMRBox(vtkPlane* pl, double bounds[6]);
159 bool PlaneIntersectsAMRBox(double plane[4], double bounds[6]);
161
166
171
176 vtkPlane* cutPlane, unsigned int blockIdx, vtkUniformGrid* grid, vtkMultiBlockDataSet* dataSet);
177
179 double Center[3];
180 double Normal[3];
184
185 std::vector<int> BlocksToLoad;
186
187private:
188 vtkAMRCutPlane(const vtkAMRCutPlane&) = delete;
189 void operator=(const vtkAMRCutPlane&) = delete;
190};
191
192#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.
std::vector< int > BlocksToLoad
void InitializeCenter(double min[3], double max[3])
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.
void CutAMRBlock(vtkPlane *cutPlane, unsigned int blockIdx, vtkUniformGrid *grid, vtkMultiBlockDataSet *dataSet)
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.
vtkMultiProcessController * Controller
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
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
a simple class to control print indentation
Definition: vtkIndent.h:113
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:143
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
image data with blanking
@ Normal
Definition: vtkX3D.h:51
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
int vtkIdType
Definition: vtkType.h:332
#define max(a, b)