VTK  9.3.0
vtkBoxClipDataSet.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3 // SPDX-License-Identifier: BSD-3-Clause
4 
36 #ifndef vtkBoxClipDataSet_h
37 #define vtkBoxClipDataSet_h
38 
39 #include "vtkFiltersGeneralModule.h" // For export macro
41 
42 VTK_ABI_NAMESPACE_BEGIN
43 class vtkCell3D;
44 class vtkCellArray;
45 class vtkCellData;
46 class vtkDataArray;
48 class vtkIdList;
49 class vtkGenericCell;
50 class vtkPointData;
52 class vtkPoints;
53 
54 class VTKFILTERSGENERAL_EXPORT vtkBoxClipDataSet : public vtkUnstructuredGridAlgorithm
55 {
56 public:
58  void PrintSelf(ostream& os, vtkIndent indent) override;
59 
66 
68 
73  void SetBoxClip(double xmin, double xmax, double ymin, double ymax, double zmin, double zmax);
74  void SetBoxClip(const double* n0, const double* o0, const double* n1, const double* o1,
75  const double* n2, const double* o2, const double* n3, const double* o3, const double* n4,
76  const double* o4, const double* n5, const double* o5);
78 
80 
84  vtkSetMacro(GenerateClipScalars, vtkTypeBool);
85  vtkGetMacro(GenerateClipScalars, vtkTypeBool);
86  vtkBooleanMacro(GenerateClipScalars, vtkTypeBool);
88 
90 
94  vtkSetMacro(GenerateClippedOutput, vtkTypeBool);
95  vtkGetMacro(GenerateClippedOutput, vtkTypeBool);
96  vtkBooleanMacro(GenerateClippedOutput, vtkTypeBool);
98 
109 
113  virtual int GetNumberOfOutputs();
115 
117 
122  vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
124 
130 
134  vtkMTimeType GetMTime() override;
135 
137 
141  vtkGetMacro(Orientation, unsigned int);
142  vtkSetMacro(Orientation, unsigned int);
144 
145  static void InterpolateEdge(vtkDataSetAttributes* attributes, vtkIdType toId, vtkIdType fromId1,
146  vtkIdType fromId2, double t);
147 
148  void MinEdgeF(const unsigned int* id_v, const vtkIdType* cellIds, unsigned int* edgF);
150  const vtkIdType* pyramId, const vtkIdType* cellIds, vtkCellArray* newCellArray);
151  void WedgeToTetra(const vtkIdType* wedgeId, const vtkIdType* cellIds, vtkCellArray* newCellArray);
152  void CellGrid(
153  vtkIdType typeobj, vtkIdType npts, const vtkIdType* cellIds, vtkCellArray* newCellArray);
154  void CreateTetra(vtkIdType npts, const vtkIdType* cellIds, vtkCellArray* newCellArray);
155  void ClipBox(vtkPoints* newPoints, vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
156  vtkCellArray* tets, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
157  vtkIdType cellId, vtkCellData* outCD);
158  void ClipHexahedron(vtkPoints* newPoints, vtkGenericCell* cell,
160  vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
162  vtkCellArray** tets, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
163  vtkIdType cellId, vtkCellData** outCD);
165  vtkIncrementalPointLocator* locator, vtkCellArray** tets, vtkPointData* inPD,
166  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
167 
169  vtkCellArray* tets, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
170  vtkIdType cellId, vtkCellData* outCD);
171  void ClipBoxInOut2D(vtkPoints* newPoints, vtkGenericCell* cell,
172  vtkIncrementalPointLocator* locator, vtkCellArray** tets, vtkPointData* inPD,
173  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
174  void ClipHexahedron2D(vtkPoints* newPoints, vtkGenericCell* cell,
176  vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
178  vtkIncrementalPointLocator* locator, vtkCellArray** tets, vtkPointData* inPD,
179  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
180 
182  vtkCellArray* lines, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
183  vtkIdType cellId, vtkCellData* outCD);
184  void ClipBoxInOut1D(vtkPoints* newPoints, vtkGenericCell* cell,
185  vtkIncrementalPointLocator* locator, vtkCellArray** lines, vtkPointData* inPD,
186  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
187  void ClipHexahedron1D(vtkPoints* newPoints, vtkGenericCell* cell,
188  vtkIncrementalPointLocator* locator, vtkCellArray* lines, vtkPointData* inPD,
189  vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
191  vtkIncrementalPointLocator* locator, vtkCellArray** lines, vtkPointData* inPD,
192  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
193 
195  vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId,
196  vtkCellData* outCD);
198  vtkCellArray** verts, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
199  vtkIdType cellId, vtkCellData** outCD);
201  vtkCellArray* verts, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
202  vtkIdType cellId, vtkCellData* outCD);
204  vtkCellArray** verts, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
205  vtkIdType cellId, vtkCellData** outCD);
206 
207 protected:
209  ~vtkBoxClipDataSet() override;
210 
213 
216 
218 
219  // double MergeTolerance;
220 
221  double BoundBoxClip[3][2];
222  unsigned int Orientation;
223  double PlaneNormal[6][3]; // normal of each plane
224  double PlanePoint[6][3]; // point on the plane
225 
226 private:
227  vtkBoxClipDataSet(const vtkBoxClipDataSet&) = delete;
228  void operator=(const vtkBoxClipDataSet&) = delete;
229 };
230 
231 VTK_ABI_NAMESPACE_END
232 #endif
clip an unstructured grid
vtkIncrementalPointLocator * Locator
void ClipBox2D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *tets, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
void ClipBoxInOut(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **tets, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
void CreateDefaultLocator()
Create default locator.
void ClipHexahedron(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *tets, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
void ClipHexahedron1D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *lines, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
void CellGrid(vtkIdType typeobj, vtkIdType npts, const vtkIdType *cellIds, vtkCellArray *newCellArray)
void CreateTetra(vtkIdType npts, const vtkIdType *cellIds, vtkCellArray *newCellArray)
void ClipHexahedronInOut(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **tets, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
void ClipBoxInOut1D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **lines, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
void ClipHexahedronInOut2D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **tets, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
static vtkBoxClipDataSet * New()
Constructor of the clipping box.
vtkTypeBool GenerateClipScalars
void ClipBoxInOut2D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **tets, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
static void InterpolateEdge(vtkDataSetAttributes *attributes, vtkIdType toId, vtkIdType fromId1, vtkIdType fromId2, double t)
void ClipHexahedronInOut0D(vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **verts, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
void ClipBox1D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *lines, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
void SetLocator(vtkIncrementalPointLocator *locator)
Specify a spatial locator for merging points.
~vtkBoxClipDataSet() override
vtkTypeBool GenerateClippedOutput
void WedgeToTetra(const vtkIdType *wedgeId, const vtkIdType *cellIds, vtkCellArray *newCellArray)
void ClipBox0D(vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void MinEdgeF(const unsigned int *id_v, const vtkIdType *cellIds, unsigned int *edgF)
void ClipBoxInOut0D(vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **verts, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
void ClipHexahedronInOut1D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray **lines, vtkPointData *inPD, vtkPointData **outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData **outCD)
void ClipHexahedron0D(vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
void SetBoxClip(double xmin, double xmax, double ymin, double ymax, double zmin, double zmax)
Specify the Box with which to perform the clipping.
vtkMTimeType GetMTime() override
Return the mtime also considering the locator.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void ClipHexahedron2D(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *tets, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
void PyramidToTetra(const vtkIdType *pyramId, const vtkIdType *cellIds, vtkCellArray *newCellArray)
virtual int GetNumberOfOutputs()
Set the tolerance for merging clip intersection points that are near the vertices of cells.
vtkUnstructuredGrid * GetClippedOutput()
Set the tolerance for merging clip intersection points that are near the vertices of cells.
unsigned int Orientation
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void SetBoxClip(const double *n0, const double *o0, const double *n1, const double *o1, const double *n2, const double *o2, const double *n3, const double *o3, const double *n4, const double *o4, const double *n5, const double *o5)
Specify the Box with which to perform the clipping.
void ClipBox(vtkPoints *newPoints, vtkGenericCell *cell, vtkIncrementalPointLocator *locator, vtkCellArray *tets, vtkPointData *inPD, vtkPointData *outPD, vtkCellData *inCD, vtkIdType cellId, vtkCellData *outCD)
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:28
object to represent cell connectivity
Definition: vtkCellArray.h:176
represent and manipulate cell attribute data
Definition: vtkCellData.h:31
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
represent and manipulate attribute data in a dataset
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:23
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:30
represent and manipulate 3D points
Definition: vtkPoints.h:29
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270