VTK  9.5.20250531
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
107#include "vtkFiltersGeneralModule.h" // For export macro
109
110VTK_ABI_NAMESPACE_BEGIN
111class vtkDataArray;
112class vtkDataSet;
114class vtkPointLocator;
115class vtkPoints;
119
120class VTKFILTERSGENERAL_EXPORT vtkTessellatorFilter : public vtkUnstructuredGridAlgorithm
121{
122public:
124 void PrintSelf(ostream& os, vtkIndent indent) override;
125
127
129 vtkGetObjectMacro(Tessellator, vtkStreamingTessellator);
130
132 vtkGetObjectMacro(Subdivider, vtkDataSetEdgeSubdivisionCriterion);
133
135
137
145 vtkSetClampMacro(OutputDimension, int, 1, 3);
146 vtkGetMacro(OutputDimension, int);
148
149// With VTK_USE_FUTURE_CONST, vtkGetMacro already makes the member const.
150#if !VTK_USE_FUTURE_CONST
151 int GetOutputDimension() const;
152#endif
153
155
160 virtual void SetMaximumNumberOfSubdivisions(int num_subdiv_in);
162 virtual void SetChordError(double ce);
165
167
170 virtual void ResetFieldCriteria();
171 virtual void SetFieldCriterion(int field, double err);
173
175
181 vtkGetMacro(MergePoints, vtkTypeBool);
182 vtkSetMacro(MergePoints, vtkTypeBool);
183 vtkBooleanMacro(MergePoints, vtkTypeBool);
185
186protected:
189
190 int FillInputPortInformation(int port, vtkInformation* info) override;
191
198
203
207 void Teardown();
208
213 vtkInformationVector* outputVector) override;
214
220
222
231
232 static void AddAPoint(const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
233 static void AddALine(
234 const double*, const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
235 static void AddATriangle(
236 const double*, const double*, const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
237 static void AddATetrahedron(const double*, const double*, const double*, const double*,
238 vtkEdgeSubdivisionCriterion*, void*, const void*);
239 void OutputPoint(const double*);
240 void OutputLine(const double*, const double*);
241 void OutputTriangle(const double*, const double*, const double*);
242 void OutputTetrahedron(const double*, const double*, const double*, const double*);
243
244private:
246 void operator=(const vtkTessellatorFilter&) = delete;
247};
248
249// With VTK_USE_FUTURE_CONST, vtkGetMacro already makes the member const.
250#if !VTK_USE_FUTURE_CONST
252{
253 return this->OutputDimension;
254}
255#endif
256
257VTK_ABI_NAMESPACE_END
258#endif // vtkTessellatorFilter_h
abstract superclass for arrays of numeric data
a subclass of vtkEdgeSubdivisionCriterion for vtkDataSet objects.
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
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:108
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:139
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.
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.
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.
static vtkTessellatorFilter * New()
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
int vtkTypeBool
Definition vtkABI.h:64
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287