VTK  9.0.1
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 
118 #ifndef vtkReebGraph_h
119 #define vtkReebGraph_h
120 
121 #include "vtkCommonDataModelModule.h" // For export macro
122 #include "vtkMutableDirectedGraph.h"
123 
124 class vtkDataArray;
125 class vtkDataSet;
126 class vtkIdList;
127 class vtkPolyData;
129 class vtkUnstructuredGrid;
130 
131 class VTKCOMMONDATAMODEL_EXPORT vtkReebGraph : public vtkMutableDirectedGraph
132 {
133 
134 public:
135  static vtkReebGraph* New();
136 
138  void PrintSelf(ostream& os, vtkIndent indent) override;
139  void PrintNodeData(ostream& os, vtkIndent indent);
140 
147  int GetDataObjectType() override { return VTK_REEB_GRAPH; }
148 
149  enum
150  {
151  ERR_INCORRECT_FIELD = -1,
152  ERR_NO_SUCH_FIELD = -2,
153  ERR_NOT_A_SIMPLICIAL_MESH = -3
154  };
155 
169  int Build(vtkPolyData* mesh, vtkDataArray* scalarField);
170 
183  int Build(vtkUnstructuredGrid* mesh, vtkDataArray* scalarField);
184 
201  int Build(vtkPolyData* mesh, vtkIdType scalarFieldId);
202 
218  int Build(vtkUnstructuredGrid* mesh, vtkIdType scalarFieldId);
219 
236  int Build(vtkPolyData* mesh, const char* scalarFieldName);
237 
253  int Build(vtkUnstructuredGrid* mesh, const char* scalarFieldName);
254 
268  int StreamTriangle(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1,
269  vtkIdType vertex2Id, double scalar2);
270 
285  int StreamTetrahedron(vtkIdType vertex0Id, double scalar0, vtkIdType vertex1Id, double scalar1,
286  vtkIdType vertex2Id, double scalar2, vtkIdType vertex3Id, double scalar3);
287 
298  void CloseStream();
299 
300  // Description:
301  // Implements deep copy
302  void DeepCopy(vtkDataObject* src) override;
303 
345  int Simplify(
346  double simplificationThreshold, vtkReebGraphSimplificationMetric* simplificationMetric);
347 
353 
354 protected:
356  ~vtkReebGraph() override;
357 
358  class Implementation;
359  Implementation* Storage;
360 
361 private:
362  vtkReebGraph(const vtkReebGraph&) = delete;
363  void operator=(const vtkReebGraph&) = delete;
364 };
365 
366 #endif
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
general representation of visualization data
Definition: vtkDataObject.h:60
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
list of point or cell ids
Definition: vtkIdList.h:31
a simple class to control print indentation
Definition: vtkIndent.h:34
An editable directed graph.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
abstract class for custom Reeb graph simplification metric design.
Reeb graph computation for PL scalar fields.
Definition: vtkReebGraph.h:132
int Simplify(double simplificationThreshold, vtkReebGraphSimplificationMetric *simplificationMetric)
Simplify the Reeb graph given a threshold 'simplificationThreshold' (between 0 and 1).
static vtkReebGraph * New()
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'.
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:147
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:358
~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:338
#define VTK_REEB_GRAPH
Definition: vtkType.h:113