VTK  9.2.6
vtkReebGraph.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkReebGraph.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
15/*----------------------------------------------------------------------------
16 Copyright (c) Sandia Corporation
17 See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18----------------------------------------------------------------------------*/
19
121#ifndef vtkReebGraph_h
122#define vtkReebGraph_h
123
124#include "vtkCommonDataModelModule.h" // For export macro
126
127class vtkDataArray;
128class vtkDataSet;
129class vtkIdList;
130class vtkPolyData;
133
134class VTKCOMMONDATAMODEL_EXPORT vtkReebGraph : public vtkMutableDirectedGraph
135{
136
137public:
138 static vtkReebGraph* New();
139
141 void PrintSelf(ostream& os, vtkIndent indent) override;
142 void PrintNodeData(ostream& os, vtkIndent indent);
143
150 int GetDataObjectType() override { return VTK_REEB_GRAPH; }
151
152 enum
153 {
154 ERR_INCORRECT_FIELD = -1,
155 ERR_NO_SUCH_FIELD = -2,
156 ERR_NOT_A_SIMPLICIAL_MESH = -3
157 };
158
172 int Build(vtkPolyData* mesh, vtkDataArray* scalarField);
173
186 int Build(vtkUnstructuredGrid* mesh, vtkDataArray* scalarField);
187
204 int Build(vtkPolyData* mesh, vtkIdType scalarFieldId);
205
221 int Build(vtkUnstructuredGrid* mesh, vtkIdType scalarFieldId);
222
239 int Build(vtkPolyData* mesh, const char* scalarFieldName);
240
256 int Build(vtkUnstructuredGrid* mesh, const char* scalarFieldName);
257
271 int StreamTriangle(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1,
272 vtkIdType vertex2Id, double scalar2);
273
288 int StreamTetrahedron(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1,
289 vtkIdType vertex2Id, double scalar2, vtkIdType vertex3Id, double scalar3);
290
302
303 // Description:
304 // Implements deep copy
305 void DeepCopy(vtkDataObject* src) override;
306
349 double simplificationThreshold, vtkReebGraphSimplificationMetric* simplificationMetric);
350
356
357protected:
359 ~vtkReebGraph() override;
360
361 class Implementation;
362 Implementation* Storage;
363
364private:
365 vtkReebGraph(const vtkReebGraph&) = delete;
366 void operator=(const vtkReebGraph&) = delete;
367};
368
369#endif
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:56
general representation of visualization data
Definition: vtkDataObject.h:66
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
list of point or cell ids
Definition: vtkIdList.h:34
a simple class to control print indentation
Definition: vtkIndent.h:40
An editable directed graph.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:91
abstract class for custom Reeb graph simplification metric design.
Reeb graph computation for PL scalar fields.
Definition: vtkReebGraph.h:135
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.
Definition: vtkReebGraph.h:150
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
Definition: vtkReebGraph.h:362
~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
#define VTK_REEB_GRAPH
Definition: vtkType.h:105