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
57#include "vtkFiltersGeneralModule.h" // For export macro
59
60VTK_ABI_NAMESPACE_BEGIN
61class vtkDataArray;
62class vtkDataSet;
64class vtkPointLocator;
65class vtkPoints;
69
70class VTKFILTERSGENERAL_EXPORT vtkTessellatorFilter : public vtkUnstructuredGridAlgorithm
71{
72public:
74 void PrintSelf(ostream& os, vtkIndent indent) override;
75
77
79 vtkGetObjectMacro(Tessellator, vtkStreamingTessellator);
80
82 vtkGetObjectMacro(Subdivider, vtkDataSetEdgeSubdivisionCriterion);
83
85
87
95 vtkSetClampMacro(OutputDimension, int, 1, 3);
96 vtkGetMacro(OutputDimension, int);
98
99// With VTK_USE_FUTURE_CONST, vtkGetMacro already makes the member const.
100#if !VTK_USE_FUTURE_CONST
101 int GetOutputDimension() const;
102#endif
103
105
110 virtual void SetMaximumNumberOfSubdivisions(int num_subdiv_in);
112 virtual void SetChordError(double ce);
115
117
120 virtual void ResetFieldCriteria();
121 virtual void SetFieldCriterion(int field, double err);
123
125
131 vtkGetMacro(MergePoints, vtkTypeBool);
132 vtkSetMacro(MergePoints, vtkTypeBool);
133 vtkBooleanMacro(MergePoints, vtkTypeBool);
135
136protected:
139
140 int FillInputPortInformation(int port, vtkInformation* info) override;
141
148
153
157 void Teardown();
158
163 vtkInformationVector* outputVector) override;
164
170
172
181
182 static void AddAPoint(const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
183 static void AddALine(
184 const double*, const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
185 static void AddATriangle(
186 const double*, const double*, const double*, vtkEdgeSubdivisionCriterion*, void*, const void*);
187 static void AddATetrahedron(const double*, const double*, const double*, const double*,
188 vtkEdgeSubdivisionCriterion*, void*, const void*);
189 void OutputPoint(const double*);
190 void OutputLine(const double*, const double*);
191 void OutputTriangle(const double*, const double*, const double*);
192 void OutputTetrahedron(const double*, const double*, const double*, const double*);
193
194private:
196 void operator=(const vtkTessellatorFilter&) = delete;
197};
198
199// With VTK_USE_FUTURE_CONST, vtkGetMacro already makes the member const.
200#if !VTK_USE_FUTURE_CONST
202{
203 return this->OutputDimension;
204}
205#endif
206
207VTK_ABI_NAMESPACE_END
208#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:62
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:38
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:38
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:270