VTK  9.6.20260204
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
22
23#ifndef vtkAdaptiveDataSetSurfaceFilter_h
24#define vtkAdaptiveDataSetSurfaceFilter_h
25
26#include "vtkFiltersHybridModule.h" // For export macro
27#include "vtkGeometryFilter.h"
28
29VTK_ABI_NAMESPACE_BEGIN
30class vtkBitArray;
31class vtkCamera;
33class vtkMatrix4x4;
34class vtkRenderer;
35
38
39class VTKFILTERSHYBRID_EXPORT vtkAdaptiveDataSetSurfaceFilter : public vtkGeometryFilter
40{
41public:
44 void PrintSelf(ostream& os, vtkIndent indent) override;
45
47
51 vtkGetObjectMacro(Renderer, vtkRenderer);
53
58
60
65 vtkSetMacro(ViewPointDepend, bool);
66 vtkGetMacro(ViewPointDepend, bool);
68
70
75 vtkSetMacro(FixedLevelMax, int);
76 vtkGetMacro(FixedLevelMax, int);
78
79protected:
82
83 int RequestData(vtkInformation* vtkNotUsed(request), vtkInformationVector** inputVector,
84 vtkInformationVector* outputVector) override;
86 int FillInputPortInformation(int port, vtkInformation* info) override;
87
88private:
90 void operator=(const vtkAdaptiveDataSetSurfaceFilter&) = delete;
91
92 enum class ShapeState : uint8_t;
93
100 template <int N>
101 ShapeState IsShapeVisible(const std::array<std::array<double, 3>, N>& points, int level);
102
106 void ProcessTrees(vtkHyperTreeGrid* input, vtkPolyData* output);
107
111 void RecursivelyProcessTree1D(vtkHyperTreeGridNonOrientedGeometryCursor*, int);
112 void RecursivelyProcessTree2D(vtkHyperTreeGridNonOrientedGeometryCursor*, int);
113 void RecursivelyProcessTree3D(vtkHyperTreeGridNonOrientedVonNeumannSuperCursorLight*, int);
114
119
124
129
133 void AddFace(vtkIdType, const double*, const double*, int, unsigned int);
134
135 vtkDataSetAttributes* InData = nullptr;
136 vtkDataSetAttributes* OutData = nullptr;
137
141 unsigned int Dimension = 0;
142
146 unsigned int Orientation = 0;
147
151 vtkBitArray* Mask;
152
156 vtkPoints* Points = nullptr;
157
161 vtkCellArray* Cells = nullptr;
162
166 vtkRenderer* Renderer = nullptr;
167
171 unsigned int Axis1;
172
176 unsigned int Axis2;
177
181 int LastRendererSize[2] = { 0, 0 };
182
186 bool ViewPointDepend = true;
187
191 int FixedLevelMax = -1;
192
196 bool IsParallel = false;
197
201 int MaxLevel = VTK_INT_MAX;
202
203 vtkSmartPointer<vtkMatrix4x4> ModelViewMatrix;
204 vtkSmartPointer<vtkMatrix4x4> ProjectionMatrix;
205};
206
207VTK_ABI_NAMESPACE_END
208#endif // vtkAdaptiveDataSetSurfaceFilter_h
vtkMTimeType GetMTime() override
Get the mtime of this object.
void SetRenderer(vtkRenderer *ren)
Set/Get the renderer attached to this adaptive surface extractor.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
int DataObjectExecute(vtkDataObject *input, vtkPolyData *output)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
static vtkAdaptiveDataSetSurfaceFilter * New()
dynamic, self-adjusting array of bits
Definition vtkBitArray.h:31
a virtual camera for 3D rendering
Definition vtkCamera.h:151
object to represent cell connectivity
general representation of visualization data
represent and manipulate attribute data in a dataset
A dataset containing a grid of vtkHyperTree instances arranged as a rectilinear grid.
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 4x4 transformation matrices
represent and manipulate 3D points
Definition vtkPoints.h:139
concrete dataset represents vertices, lines, polygons, and triangle strips
abstract specification for renderers
Hold a reference to a vtkObjectBase instance.
int vtkIdType
Definition vtkType.h:354
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:309
#define VTK_INT_MAX
Definition vtkType.h:183