VTK  9.0.1
vtkAdaptiveDataSetSurfaceFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkAdaptiveDataSetSurfaceFilter.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 =========================================================================*/
33 #ifndef vtkAdaptiveDataSetSurfaceFilter_h
34 #define vtkAdaptiveDataSetSurfaceFilter_h
35 
36 #include "vtkFiltersHybridModule.h" // For export macro
37 #include "vtkGeometryFilter.h"
38 
39 class vtkBitArray;
40 class vtkCamera;
41 class vtkHyperTreeGrid;
42 class vtkRenderer;
43 
46 
47 class VTKFILTERSHYBRID_EXPORT vtkAdaptiveDataSetSurfaceFilter : public vtkGeometryFilter
48 {
49 public:
52  void PrintSelf(ostream&, vtkIndent) override;
53 
55 
59  vtkGetObjectMacro(Renderer, vtkRenderer);
61 
65  vtkMTimeType GetMTime() override;
66 
68 
71  vtkSetMacro(CircleSelection, bool);
72  vtkGetMacro(CircleSelection, bool);
74 
76 
81  vtkSetMacro(BBSelection, bool);
82  vtkGetMacro(BBSelection, bool);
84 
86 
89  vtkSetMacro(ViewPointDepend, bool);
90  vtkGetMacro(ViewPointDepend, bool);
92 
94 
97  vtkSetMacro(FixedLevelMax, int);
98  vtkGetMacro(FixedLevelMax, int);
100 
102 
107  vtkSetMacro(Scale, double);
108  vtkGetMacro(Scale, double);
110 
112 
117  vtkSetMacro(DynamicDecimateLevelMax, int);
118  vtkGetMacro(DynamicDecimateLevelMax, int);
120 
121 protected:
124 
125  int RequestData(vtkInformation* vtkNotUsed(request), vtkInformationVector** inputVector,
126  vtkInformationVector* outputVector) override;
127  int DataSetExecute(vtkDataObject* input, vtkPolyData* output) /*override*/;
129 
134 
140 
145 
150 
155 
159  void AddFace(vtkIdType, const double*, const double*, int, unsigned int);
160 
163 
167  unsigned int Dimension;
168 
172  unsigned int Orientation;
173 
178 
183 
188 
193 
197  unsigned int Axis1;
198 
202  unsigned int Axis2;
203 
207  int LevelMax;
208 
213 
217  int LastRendererSize[2];
218 
223 
227  double LastCameraFocalPoint[3];
228 
233 
237  double WindowBounds[4];
238 
243 
247  double Radius;
248 
253 
254 #ifndef NDEBUG
259  long int NbRejectByBB;
260 #endif
261 
266 
270  double Scale;
271 
276 
277 private:
279  void operator=(const vtkAdaptiveDataSetSurfaceFilter&) = delete;
280 };
281 
282 #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.
int DataSetExecute(vtkDataObject *input, vtkPolyData *output)
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)
void PrintSelf(ostream &, vtkIndent) override
Methods invoked by print to print information about the object including superclasses.
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.
long int NbRejectByCircle
Effect of options selection.
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 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:34
a virtual camera for 3D rendering
Definition: vtkCamera.h:46
object to represent cell connectivity
Definition: vtkCellArray.h:180
general representation of visualization data
Definition: vtkDataObject.h:60
represent and manipulate attribute data in a dataset
extract geometry from data (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:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition: vtkPoints.h:34
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
abstract specification for renderers
Definition: vtkRenderer.h:59
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
int vtkIdType
Definition: vtkType.h:338
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:293