VTK  9.3.0
vtkTessellatorFilter.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright 2003 Sandia Corporation
3 // SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-NVIDIA-USGov
4 #ifndef vtkTessellatorFilter_h
5 #define vtkTessellatorFilter_h
6 
48 #include "vtkFiltersGeneralModule.h" // For export macro
50 
51 VTK_ABI_NAMESPACE_BEGIN
52 class vtkDataArray;
53 class vtkDataSet;
55 class vtkPointLocator;
56 class vtkPoints;
60 
61 class VTKFILTERSGENERAL_EXPORT vtkTessellatorFilter : public vtkUnstructuredGridAlgorithm
62 {
63 public:
65  void PrintSelf(ostream& os, vtkIndent indent) override;
66 
68 
70  vtkGetObjectMacro(Tessellator, vtkStreamingTessellator);
71 
73  vtkGetObjectMacro(Subdivider, vtkDataSetEdgeSubdivisionCriterion);
74 
75  vtkMTimeType GetMTime() override;
76 
78 
86  vtkSetClampMacro(OutputDimension, int, 1, 3);
87  vtkGetMacro(OutputDimension, int);
89 
90 // With VTK_USE_FUTURE_CONST, vtkGetMacro already makes the member const.
91 #if !VTK_USE_FUTURE_CONST
92  int GetOutputDimension() const;
93 #endif
94 
96 
101  virtual void SetMaximumNumberOfSubdivisions(int num_subdiv_in);
103  virtual void SetChordError(double ce);
104  double GetChordError();
106 
108 
111  virtual void ResetFieldCriteria();
112  virtual void SetFieldCriterion(int field, double err);
114 
116 
122  vtkGetMacro(MergePoints, vtkTypeBool);
123  vtkSetMacro(MergePoints, vtkTypeBool);
124  vtkBooleanMacro(MergePoints, vtkTypeBool);
126 
127 protected:
130 
132 
139 
144 
148  void Teardown();
149 
153  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
154  vtkInformationVector* outputVector) override;
155 
161 
163 
172 
173  static void AddAPoint(const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
174  static void AddALine(
175  const double*, const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
176  static void AddATriangle(
177  const double*, const double*, const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
178  static void AddATetrahedron(const double*, const double*, const double*, const double*,
179  vtkEdgeSubdivisionCriterion*, void*, const void*);
180  void OutputPoint(const double*);
181  void OutputLine(const double*, const double*);
182  void OutputTriangle(const double*, const double*, const double*);
183  void OutputTetrahedron(const double*, const double*, const double*, const double*);
184 
185 private:
187  void operator=(const vtkTessellatorFilter&) = delete;
188 };
189 
190 // With VTK_USE_FUTURE_CONST, vtkGetMacro already makes the member const.
191 #if !VTK_USE_FUTURE_CONST
193 {
194  return this->OutputDimension;
195 }
196 #endif
197 
198 VTK_ABI_NAMESPACE_END
199 #endif // vtkTessellatorFilter_h
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
a subclass of vtkEdgeSubdivisionCriterion for vtkDataSet objects.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:53
how to decide whether a linear approximation to nonlinear geometry or field should be subdivided
a simple class to control print indentation
Definition: vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
quickly locate points in 3-space
represent and manipulate 3D points
Definition: vtkPoints.h:29
An algorithm that refines an initial simplicial tessellation using edge subdivision.
approximate nonlinear FEM elements with simplices
void OutputTetrahedron(const double *, const double *, const double *, const double *)
void MergeOutputPoints(vtkUnstructuredGrid *input, vtkUnstructuredGrid *output)
Called by RequestData to merge output points.
vtkPointLocator * Locator
virtual void SetSubdivider(vtkDataSetEdgeSubdivisionCriterion *)
virtual void SetChordError(double ce)
These are convenience routines for setting properties maintained by the tessellator and subdivider.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual void SetTessellator(vtkStreamingTessellator *)
~vtkTessellatorFilter() override
void SetupOutput(vtkDataSet *input, vtkUnstructuredGrid *output)
Called by RequestData to set up a multitude of member variables used by the per-primitive output func...
vtkPoints * OutputPoints
These member variables are set by SetupOutput for use inside the callback members OutputLine and Outp...
vtkDataArray ** OutputAttributes
These member variables are set by SetupOutput for use inside the callback members OutputLine and Outp...
static void AddAPoint(const double *, vtkEdgeSubdivisionCriterion *, void *, const void *)
virtual void SetFieldCriterion(int field, double err)
These methods are for the ParaView client.
static void AddATriangle(const double *, const double *, const double *, vtkEdgeSubdivisionCriterion *, void *, const void *)
static void AddATetrahedron(const double *, const double *, const double *, const double *, vtkEdgeSubdivisionCriterion *, void *, const void *)
vtkMTimeType GetMTime() override
Return this object's modified time.
static vtkTessellatorFilter * New()
virtual void ResetFieldCriteria()
These methods are for the ParaView client.
int GetMaximumNumberOfSubdivisions()
These are convenience routines for setting properties maintained by the tessellator and subdivider.
vtkDataSetEdgeSubdivisionCriterion * Subdivider
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void Teardown()
Reset the temporary variables used during the filter's RequestData() method.
int * OutputAttributeIndices
These member variables are set by SetupOutput for use inside the callback members OutputLine and Outp...
virtual int GetOutputDimension()
Set the dimension of the output tessellation.
void OutputTriangle(const double *, const double *, const double *)
vtkUnstructuredGrid * OutputMesh
These member variables are set by SetupOutput for use inside the callback members OutputLine and Outp...
void OutputLine(const double *, const double *)
void OutputPoint(const double *)
static void AddALine(const double *, const double *, vtkEdgeSubdivisionCriterion *, void *, const void *)
double GetChordError()
These are convenience routines for setting properties maintained by the tessellator and subdivider.
virtual void SetMaximumNumberOfSubdivisions(int num_subdiv_in)
These are convenience routines for setting properties maintained by the tessellator and subdivider.
vtkStreamingTessellator * Tessellator
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Run the filter; produce a polygonal approximation to the grid.
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
@ field
Definition: vtkX3D.h:177
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
int vtkTypeBool
Definition: vtkABI.h:64
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270