VTK  9.3.0
vtkAdaptiveDataSetSurfaceFilter.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
24 #ifndef vtkAdaptiveDataSetSurfaceFilter_h
25 #define vtkAdaptiveDataSetSurfaceFilter_h
26 
27 #include "vtkFiltersHybridModule.h" // For export macro
28 #include "vtkGeometryFilter.h"
29 
30 VTK_ABI_NAMESPACE_BEGIN
31 class vtkBitArray;
32 class vtkCamera;
33 class vtkHyperTreeGrid;
34 class vtkRenderer;
35 
38 
39 class VTKFILTERSHYBRID_EXPORT vtkAdaptiveDataSetSurfaceFilter : public vtkGeometryFilter
40 {
41 public:
44  void PrintSelf(ostream& os, vtkIndent indent) override;
45 
47 
51  vtkGetObjectMacro(Renderer, vtkRenderer);
53 
57  vtkMTimeType GetMTime() override;
58 
60 
63  vtkSetMacro(CircleSelection, bool);
64  vtkGetMacro(CircleSelection, bool);
66 
68 
73  vtkSetMacro(BBSelection, bool);
74  vtkGetMacro(BBSelection, bool);
76 
78 
81  vtkSetMacro(ViewPointDepend, bool);
82  vtkGetMacro(ViewPointDepend, bool);
84 
86 
89  vtkSetMacro(FixedLevelMax, int);
90  vtkGetMacro(FixedLevelMax, int);
92 
94 
99  vtkSetMacro(Scale, double);
100  vtkGetMacro(Scale, double);
102 
104 
109  vtkSetMacro(DynamicDecimateLevelMax, int);
110  vtkGetMacro(DynamicDecimateLevelMax, int);
112 
113 protected:
116 
117  int RequestData(vtkInformation* vtkNotUsed(request), vtkInformationVector** inputVector,
118  vtkInformationVector* outputVector) override;
121 
126 
132 
137 
142 
147 
151  void AddFace(vtkIdType, const double*, const double*, int, unsigned int);
152 
155 
159  unsigned int Dimension;
160 
164  unsigned int Orientation;
165 
170 
175 
180 
185 
189  unsigned int Axis1;
190 
194  unsigned int Axis2;
195 
199  int LevelMax;
200 
205 
209  int LastRendererSize[2];
210 
215 
219  double LastCameraFocalPoint[3];
220 
225 
229  double WindowBounds[4];
230 
235 
239  double Radius;
240 
245 
250 
254  double Scale;
255 
260 
261 private:
263  void operator=(const vtkAdaptiveDataSetSurfaceFilter&) = delete;
264 };
265 
266 VTK_ABI_NAMESPACE_END
267 #endif // vtkAdaptiveDataSetSurfaceFilter_h
Adaptively extract dataset surface.
void RecursivelyProcessTreeNot3D(vtkHyperTreeGridNonOrientedGeometryCursor *, int)
Recursively descend into tree down to leaves.
void ProcessTrees(vtkHyperTreeGrid *input, vtkPolyData *output)
Main routine to generate external boundary.
vtkPoints * Points
Storage for points of output unstructured mesh.
double Scale
Scale factor for adaptive view.
static vtkAdaptiveDataSetSurfaceFilter * New()
vtkMTimeType GetMTime() override
Get the mtime of this object.
unsigned int Orientation
Orientation of input grid when dimension < 3.
bool BBSelection
Product cell when in nounding box selection.
double Radius
Radius parameter for adaptive view.
bool CircleSelection
Product cell when in circle selection.
vtkCellArray * Cells
Storage for cells of output unstructured mesh.
void SetRenderer(vtkRenderer *ren)
Set/Get the renderer attached to this adaptive surface extractor.
bool ViewPointDepend
JB Activation de la dependance au point de vue.
vtkRenderer * Renderer
Pointer to the renderer in use.
void ProcessLeaf2D(vtkHyperTreeGridNonOrientedGeometryCursor *)
Process 2D leaves and issue corresponding faces (quads)
unsigned int Dimension
Dimension of input grid.
double LastCameraParallelScale
Last camera parallel scale for adaptive view.
unsigned int Axis2
Second axis parameter for adaptive view.
int LevelMax
Maximum depth parameter for adaptive view.
int FixedLevelMax
JB Forced, fixed the level depth, ignored automatic determination.
void AddFace(vtkIdType, const double *, const double *, int, unsigned int)
Helper method to generate a face based on its normal and offset from cursor origin.
bool ParallelProjection
Parallel projection parameter for adaptive view.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
void ProcessLeaf3D(vtkHyperTreeGridNonOrientedVonNeumannSuperCursorLight *)
Process 3D leaves and issue corresponding cells (voxels)
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void RecursivelyProcessTree3D(vtkHyperTreeGridNonOrientedVonNeumannSuperCursorLight *, int)
void ProcessLeaf1D(vtkHyperTreeGridNonOrientedGeometryCursor *)
Process 1D leaves and issue corresponding edges (lines)
unsigned int Axis1
First axis parameter for adaptive view.
int DataObjectExecute(vtkDataObject *input, vtkPolyData *output)
int RequestData(vtkInformation *vtkNotUsed(request), vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int DynamicDecimateLevelMax
JB Decimate level max after automatic determination.
dynamic, self-adjusting array of bits
Definition: vtkBitArray.h:29
a virtual camera for 3D rendering
Definition: vtkCamera.h:50
object to represent cell connectivity
Definition: vtkCellArray.h:185
general representation of visualization data
Definition: vtkDataObject.h:64
represent and manipulate attribute data in a dataset
extract boundary geometry from dataset (or convert data to polygonal type)
A dataset containing a grid of vtkHyperTree instances arranged as a rectilinear grid.
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 3D points
Definition: vtkPoints.h:38
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:89
abstract specification for renderers
Definition: vtkRenderer.h:71
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
int vtkIdType
Definition: vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270