VTK  9.0.20210223
vtkBoxClipDataSet.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkBoxClipDataSet.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 =========================================================================*/
15 /*----------------------------------------------------------------------------
16  Copyright (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
19 
48 #ifndef vtkBoxClipDataSet_h
49 #define vtkBoxClipDataSet_h
50 
51 #include "vtkFiltersGeneralModule.h" // For export macro
53 
54 class vtkCell3D;
55 class vtkCellArray;
56 class vtkCellData;
57 class vtkDataArray;
59 class vtkIdList;
60 class vtkGenericCell;
61 class vtkPointData;
63 class vtkPoints;
64 
65 class VTKFILTERSGENERAL_EXPORT vtkBoxClipDataSet : public vtkUnstructuredGridAlgorithm
66 {
67 public:
69  void PrintSelf(ostream& os, vtkIndent indent) override;
70 
76  static vtkBoxClipDataSet* New();
77 
79 
84  void SetBoxClip(double xmin, double xmax, double ymin, double ymax, double zmin, double zmax);
85  void SetBoxClip(const double* n0, const double* o0, const double* n1, const double* o1,
86  const double* n2, const double* o2, const double* n3, const double* o3, const double* n4,
87  const double* o4, const double* n5, const double* o5);
89 
91 
95  vtkSetMacro(GenerateClipScalars, vtkTypeBool);
96  vtkGetMacro(GenerateClipScalars, vtkTypeBool);
97  vtkBooleanMacro(GenerateClipScalars, vtkTypeBool);
99 
101 
105  vtkSetMacro(GenerateClippedOutput, vtkTypeBool);
106  vtkGetMacro(GenerateClippedOutput, vtkTypeBool);
107  vtkBooleanMacro(GenerateClippedOutput, vtkTypeBool);
109 
123  vtkUnstructuredGrid* GetClippedOutput();
124  virtual int GetNumberOfOutputs();
126 
128 
132  void SetLocator(vtkIncrementalPointLocator* locator);
133  vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
135 
140  void CreateDefaultLocator();
141 
145  vtkMTimeType GetMTime() override;
146 
148 
152  vtkGetMacro(Orientation, unsigned int);
153  vtkSetMacro(Orientation, unsigned int);
155 
156  static void InterpolateEdge(vtkDataSetAttributes* attributes, vtkIdType toId, vtkIdType fromId1,
157  vtkIdType fromId2, double t);
158 
159  void MinEdgeF(const unsigned int* id_v, const vtkIdType* cellIds, unsigned int* edgF);
160  void PyramidToTetra(
161  const vtkIdType* pyramId, const vtkIdType* cellIds, vtkCellArray* newCellArray);
162  void WedgeToTetra(const vtkIdType* wedgeId, const vtkIdType* cellIds, vtkCellArray* newCellArray);
163  void CellGrid(
164  vtkIdType typeobj, vtkIdType npts, const vtkIdType* cellIds, vtkCellArray* newCellArray);
165  void CreateTetra(vtkIdType npts, const vtkIdType* cellIds, vtkCellArray* newCellArray);
166  void ClipBox(vtkPoints* newPoints, vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
167  vtkCellArray* tets, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
168  vtkIdType cellId, vtkCellData* outCD);
169  void ClipHexahedron(vtkPoints* newPoints, vtkGenericCell* cell,
171  vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
172  void ClipBoxInOut(vtkPoints* newPoints, vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
173  vtkCellArray** tets, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
174  vtkIdType cellId, vtkCellData** outCD);
175  void ClipHexahedronInOut(vtkPoints* newPoints, vtkGenericCell* cell,
176  vtkIncrementalPointLocator* locator, vtkCellArray** tets, vtkPointData* inPD,
177  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
178 
179  void ClipBox2D(vtkPoints* newPoints, vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
180  vtkCellArray* tets, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
181  vtkIdType cellId, vtkCellData* outCD);
182  void ClipBoxInOut2D(vtkPoints* newPoints, vtkGenericCell* cell,
183  vtkIncrementalPointLocator* locator, vtkCellArray** tets, vtkPointData* inPD,
184  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
185  void ClipHexahedron2D(vtkPoints* newPoints, vtkGenericCell* cell,
187  vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
188  void ClipHexahedronInOut2D(vtkPoints* newPoints, vtkGenericCell* cell,
189  vtkIncrementalPointLocator* locator, vtkCellArray** tets, vtkPointData* inPD,
190  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
191 
192  void ClipBox1D(vtkPoints* newPoints, vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
193  vtkCellArray* lines, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
194  vtkIdType cellId, vtkCellData* outCD);
195  void ClipBoxInOut1D(vtkPoints* newPoints, vtkGenericCell* cell,
196  vtkIncrementalPointLocator* locator, vtkCellArray** lines, vtkPointData* inPD,
197  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
198  void ClipHexahedron1D(vtkPoints* newPoints, vtkGenericCell* cell,
199  vtkIncrementalPointLocator* locator, vtkCellArray* lines, vtkPointData* inPD,
200  vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData* outCD);
201  void ClipHexahedronInOut1D(vtkPoints* newPoints, vtkGenericCell* cell,
202  vtkIncrementalPointLocator* locator, vtkCellArray** lines, vtkPointData* inPD,
203  vtkPointData** outPD, vtkCellData* inCD, vtkIdType cellId, vtkCellData** outCD);
204 
205  void ClipBox0D(vtkGenericCell* cell, vtkIncrementalPointLocator* locator, vtkCellArray* verts,
206  vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD, vtkIdType cellId,
207  vtkCellData* outCD);
208  void ClipBoxInOut0D(vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
209  vtkCellArray** verts, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
210  vtkIdType cellId, vtkCellData** outCD);
211  void ClipHexahedron0D(vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
212  vtkCellArray* verts, vtkPointData* inPD, vtkPointData* outPD, vtkCellData* inCD,
213  vtkIdType cellId, vtkCellData* outCD);
214  void ClipHexahedronInOut0D(vtkGenericCell* cell, vtkIncrementalPointLocator* locator,
215  vtkCellArray** verts, vtkPointData* inPD, vtkPointData** outPD, vtkCellData* inCD,
216  vtkIdType cellId, vtkCellData** outCD);
217 
218 protected:
220  ~vtkBoxClipDataSet() override;
221 
223  int FillInputPortInformation(int port, vtkInformation* info) override;
224 
227 
229 
230  // double MergeTolerance;
231 
232  double BoundBoxClip[3][2];
233  unsigned int Orientation;
234  double PlaneNormal[6][3]; // normal of each plane
235  double PlanePoint[6][3]; // point on the plane
236 
237 private:
238  vtkBoxClipDataSet(const vtkBoxClipDataSet&) = delete;
239  void operator=(const vtkBoxClipDataSet&) = delete;
240 };
241 
242 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:33
vtkBoxClipDataSet::GenerateClippedOutput
vtkTypeBool GenerateClippedOutput
Definition: vtkBoxClipDataSet.h:228
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:32
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkDataSetAttributes
represent and manipulate attribute data in a dataset
Definition: vtkDataSetAttributes.h:60
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:35
vtkBoxClipDataSet::Orientation
unsigned int Orientation
Definition: vtkBoxClipDataSet.h:233
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
vtkCell3D
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
vtkX3D::port
Definition: vtkX3D.h:453
vtkBoxClipDataSet
clip an unstructured grid
Definition: vtkBoxClipDataSet.h:65
vtkBoxClipDataSet::GenerateClipScalars
vtkTypeBool GenerateClipScalars
Definition: vtkBoxClipDataSet.h:226
vtkObject::GetMTime
virtual vtkMTimeType GetMTime()
Return this object's modified time.
vtkUnstructuredGridAlgorithm::RequestData
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:180
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkUnstructuredGridAlgorithm::FillInputPortInformation
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkBoxClipDataSet::Locator
vtkIncrementalPointLocator * Locator
Definition: vtkBoxClipDataSet.h:225
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:73
vtkX3D::info
Definition: vtkX3D.h:382
vtkUnstructuredGridAlgorithm::New
static vtkUnstructuredGridAlgorithm * New()
vtkUnstructuredGridAlgorithm::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkGenericCell
provides thread-safe access to cells
Definition: vtkGenericCell.h:36
vtkUnstructuredGridAlgorithm.h
vtkUnstructuredGridAlgorithm
Superclass for algorithms that produce only unstructured grid as output.
Definition: vtkUnstructuredGridAlgorithm.h:40
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:93
vtkTypeBool
int vtkTypeBool
Definition: vtkABI.h:69
vtkMTimeType
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:293