VTK
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
139 static vtkReebGraph *New();
140
142 void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
143 void PrintNodeData(ostream& os, vtkIndent indent);
144
151 int GetDataObjectType() VTK_OVERRIDE {return VTK_REEB_GRAPH;}
152
153
154 enum
155 {
156 ERR_INCORRECT_FIELD = -1,
157 ERR_NO_SUCH_FIELD = -2,
158 ERR_NOT_A_SIMPLICIAL_MESH = -3
159 };
160
174 int Build(vtkPolyData *mesh, vtkDataArray *scalarField);
175
188 int Build(vtkUnstructuredGrid *mesh, vtkDataArray *scalarField);
189
190
207 int Build(vtkPolyData *mesh, vtkIdType scalarFieldId);
208
224 int Build(vtkUnstructuredGrid *mesh, vtkIdType scalarFieldId);
225
226
243 int Build(vtkPolyData *mesh, const char* scalarFieldName);
244
260 int Build(vtkUnstructuredGrid *mesh, const char* scalarFieldName);
261
275 int StreamTriangle( vtkIdType vertex0Id, double scalar0,
276 vtkIdType vertex1Id, double scalar1,
277 vtkIdType vertex2Id, double scalar2);
278
293 int StreamTetrahedron( vtkIdType vertex0Id, double scalar0,
294 vtkIdType vertex1Id, double scalar1,
295 vtkIdType vertex2Id, double scalar2,
296 vtkIdType vertex3Id, double scalar3);
297
309
310 // Descrition:
311 // Implements deep copy
312 void DeepCopy(vtkDataObject *src) VTK_OVERRIDE;
313
355 int Simplify(double simplificationThreshold,
356 vtkReebGraphSimplificationMetric *simplificationMetric);
357
363
364protected:
365
367 ~vtkReebGraph() VTK_OVERRIDE;
368
369 class Implementation;
370 Implementation* Storage;
371
372private:
373 vtkReebGraph(const vtkReebGraph&) VTK_DELETE_FUNCTION;
374 void operator=(const vtkReebGraph&) VTK_DELETE_FUNCTION;
375
376};
377
378#endif
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
general representation of visualization data
Definition: vtkDataObject.h:65
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
list of point or cell ids
Definition: vtkIdList.h:37
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:86
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:151
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.
~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:287
#define VTK_REEB_GRAPH
Definition: vtkType.h:115