VTK  9.4.20250303
vtkReebGraph.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3// SPDX-License-Identifier: BSD-3-Clause
4
106#ifndef vtkReebGraph_h
107#define vtkReebGraph_h
108
109#include "vtkCommonDataModelModule.h" // For export macro
111
112VTK_ABI_NAMESPACE_BEGIN
113class vtkDataArray;
114class vtkDataSet;
115class vtkIdList;
116class vtkPolyData;
119
120class VTKCOMMONDATAMODEL_EXPORT vtkReebGraph : public vtkMutableDirectedGraph
121{
122
123public:
124 static vtkReebGraph* New();
125
127 void PrintSelf(ostream& os, vtkIndent indent) override;
128 void PrintNodeData(ostream& os, vtkIndent indent);
129
136 int GetDataObjectType() override { return VTK_REEB_GRAPH; }
137
138 enum
139 {
140 ERR_INCORRECT_FIELD = -1,
141 ERR_NO_SUCH_FIELD = -2,
142 ERR_NOT_A_SIMPLICIAL_MESH = -3
143 };
144
158 int Build(vtkPolyData* mesh, vtkDataArray* scalarField);
159
172 int Build(vtkUnstructuredGrid* mesh, vtkDataArray* scalarField);
173
190 int Build(vtkPolyData* mesh, vtkIdType scalarFieldId);
191
207 int Build(vtkUnstructuredGrid* mesh, vtkIdType scalarFieldId);
208
225 int Build(vtkPolyData* mesh, const char* scalarFieldName);
226
242 int Build(vtkUnstructuredGrid* mesh, const char* scalarFieldName);
243
257 int StreamTriangle(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1,
258 vtkIdType vertex2Id, double scalar2);
259
274 int StreamTetrahedron(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1,
275 vtkIdType vertex2Id, double scalar2, vtkIdType vertex3Id, double scalar3);
276
288
289 // Description:
290 // Implements deep copy
291 void DeepCopy(vtkDataObject* src) override;
292
335 double simplificationThreshold, vtkReebGraphSimplificationMetric* simplificationMetric);
336
342
343protected:
345 ~vtkReebGraph() override;
346
347 class Implementation;
348 Implementation* Storage;
349
350private:
351 vtkReebGraph(const vtkReebGraph&) = delete;
352 void operator=(const vtkReebGraph&) = delete;
353};
354
355VTK_ABI_NAMESPACE_END
356#endif
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
list of point or cell ids
Definition vtkIdList.h:133
a simple class to control print indentation
Definition vtkIndent.h:108
An editable directed graph.
concrete dataset represents vertices, lines, polygons, and triangle strips
abstract class for custom Reeb graph simplification metric design.
Reeb graph computation for PL scalar fields.
int Simplify(double simplificationThreshold, vtkReebGraphSimplificationMetric *simplificationMetric)
Simplify the Reeb graph given a threshold 'simplificationThreshold' (between 0 and 1).
void PrintNodeData(ostream &os, vtkIndent indent)
void CloseStream()
Finalize internal data structures, in the case of streaming computations (with StreamTriangle or Stre...
int Build(vtkPolyData *mesh, vtkIdType scalarFieldId)
Build the Reeb graph of the field given by the Id 'scalarFieldId', defined on the surface mesh 'mesh'...
int Build(vtkPolyData *mesh, const char *scalarFieldName)
Build the Reeb graph of the field given by the name 'scalarFieldName', defined on the surface mesh 'm...
int StreamTriangle(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1, vtkIdType vertex2Id, double scalar2)
Streaming Reeb graph computation.
int Build(vtkUnstructuredGrid *mesh, vtkIdType scalarFieldId)
Build the Reeb graph of the field given by the Id 'scalarFieldId', defined on the volume mesh 'mesh'.
static vtkReebGraph * New()
int Build(vtkUnstructuredGrid *mesh, vtkDataArray *scalarField)
Build the Reeb graph of the field 'scalarField' defined on the volume mesh 'mesh'.
int GetDataObjectType() override
Return class name of data type.
int Build(vtkPolyData *mesh, vtkDataArray *scalarField)
Build the Reeb graph of the field 'scalarField' defined on the surface mesh 'mesh'.
void Set(vtkMutableDirectedGraph *g)
Use a pre-defined Reeb graph (post-processing).
void DeepCopy(vtkDataObject *src) override
Deep copies the data object into this graph.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int StreamTetrahedron(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1, vtkIdType vertex2Id, double scalar2, vtkIdType vertex3Id, double scalar3)
Streaming Reeb graph computation.
Implementation * Storage
~vtkReebGraph() override
int Build(vtkUnstructuredGrid *mesh, const char *scalarFieldName)
Build the Reeb graph of the field given by the name 'scalarFieldName', defined on the volume mesh 'me...
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition vtkType.h:332
@ VTK_REEB_GRAPH
Definition vtkType.h:108