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
45#ifndef vtkBoxClipDataSet_h
46#define vtkBoxClipDataSet_h
47
48#include "vtkFiltersGeneralModule.h" // For export macro
50
51VTK_ABI_NAMESPACE_BEGIN
52class vtkCell3D;
53class vtkCellArray;
54class vtkCellData;
55class vtkDataArray;
57class vtkIdList;
58class vtkGenericCell;
59class vtkPointData;
61class vtkPoints;
62
63class VTKFILTERSGENERAL_EXPORT vtkBoxClipDataSet : public vtkUnstructuredGridAlgorithm
64{
65public:
67 void PrintSelf(ostream& os, vtkIndent indent) override;
68
75
77
82 void SetBoxClip(double xmin, double xmax, double ymin, double ymax, double zmin, double zmax);
83 void SetBoxClip(const double* n0, const double* o0, const double* n1, const double* o1,
84 const double* n2, const double* o2, const double* n3, const double* o3, const double* n4,
85 const double* o4, const double* n5, const double* o5);
87
89
93 vtkSetMacro(GenerateClipScalars, vtkTypeBool);
94 vtkGetMacro(GenerateClipScalars, vtkTypeBool);
95 vtkBooleanMacro(GenerateClipScalars, vtkTypeBool);
97
99
103 vtkSetMacro(GenerateClippedOutput, vtkTypeBool);
104 vtkGetMacro(GenerateClippedOutput, vtkTypeBool);
105 vtkBooleanMacro(GenerateClippedOutput, vtkTypeBool);
107
118
122 virtual int GetNumberOfOutputs();
124
126
131 vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
133
139
144
146
150 vtkGetMacro(Orientation, unsigned int);
151 vtkSetMacro(Orientation, unsigned int);
153
154 static void InterpolateEdge(vtkDataSetAttributes* attributes, vtkIdType toId, vtkIdType fromId1,
155 vtkIdType fromId2, double t);
156
157 void MinEdgeF(const unsigned int* id_v, const vtkIdType* cellIds, unsigned int* edgF);
159 const vtkIdType* pyramId, const vtkIdType* cellIds, vtkCellArray* newCellArray);
160 void WedgeToTetra(const vtkIdType* wedgeId, const vtkIdType* cellIds, vtkCellArray* newCellArray);
162 vtkIdType typeobj, vtkIdType npts, const vtkIdType* cellIds, vtkCellArray* newCellArray);
163 void CreateTetra(vtkIdType npts, const vtkIdType* cellIds, vtkCellArray* newCellArray);
165 vtkCellArray* tets, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
166 vtkIdType cellId, vtkCellData* outCD);
169 vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
171 vtkCellArray** tets, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
172 vtkIdType cellId, vtkCellData** outCD);
175 vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
176
178 vtkCellArray* tets, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
179 vtkIdType cellId, vtkCellData* outCD);
182 vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
185 vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
188 vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
189
191 vtkCellArray* lines, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
192 vtkIdType cellId, vtkCellData* outCD);
195 vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
198 vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
201 vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
202
204 vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId,
205 vtkCellData* outCD);
207 vtkCellArray** verts, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
208 vtkIdType cellId, vtkCellData** outCD);
210 vtkCellArray* verts, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
211 vtkIdType cellId, vtkCellData* outCD);
213 vtkCellArray** verts, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
214 vtkIdType cellId, vtkCellData** outCD);
215
216protected:
219
221 int FillInputPortInformation(int port, vtkInformation* info) override;
222
225
227
228 // double MergeTolerance;
229
230 double BoundBoxClip[3][2];
231 unsigned int Orientation;
232 double PlaneNormal[6][3]; // normal of each plane
233 double PlanePoint[6][3]; // point on the plane
234
235private:
236 vtkBoxClipDataSet(const vtkBoxClipDataSet&) = delete;
237 void operator=(const vtkBoxClipDataSet&) = delete;
238};
239
240VTK_ABI_NAMESPACE_END
241#endif
clip an unstructured grid
vtkIncrementalPointLocator * Locator
static vtkBoxClipDataSet * New()
Constructor of the clipping box.
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)
vtkUnstructuredGrid * GetClippedOutput()
Set the tolerance for merging clip intersection points that are near the vertices of cells.
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)
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.
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
represent and manipulate cell attribute data
Definition vtkCellData.h:40
abstract superclass for arrays of numeric data
represent and manipulate attribute data in a dataset
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:38
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270