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 
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 
136  static vtkReebGraph *New();
137 
139  void PrintSelf(ostream& os, vtkIndent indent) override;
140  void PrintNodeData(ostream& os, vtkIndent indent);
141 
148  int GetDataObjectType() override {return VTK_REEB_GRAPH;}
149 
150 
151  enum
152  {
153  ERR_INCORRECT_FIELD = -1,
154  ERR_NO_SUCH_FIELD = -2,
155  ERR_NOT_A_SIMPLICIAL_MESH = -3
156  };
157 
171  int Build(vtkPolyData *mesh, vtkDataArray *scalarField);
172 
185  int Build(vtkUnstructuredGrid *mesh, vtkDataArray *scalarField);
186 
187 
204  int Build(vtkPolyData *mesh, vtkIdType scalarFieldId);
205 
221  int Build(vtkUnstructuredGrid *mesh, vtkIdType scalarFieldId);
222 
223 
240  int Build(vtkPolyData *mesh, const char* scalarFieldName);
241 
257  int Build(vtkUnstructuredGrid *mesh, const char* scalarFieldName);
258 
272  int StreamTriangle( vtkIdType vertex0Id, double scalar0,
273  vtkIdType vertex1Id, double scalar1,
274  vtkIdType vertex2Id, double scalar2);
275 
290  int StreamTetrahedron( vtkIdType vertex0Id, double scalar0,
291  vtkIdType vertex1Id, double scalar1,
292  vtkIdType vertex2Id, double scalar2,
293  vtkIdType vertex3Id, double scalar3);
294 
305  void CloseStream();
306 
307  // Description:
308  // Implements deep copy
309  void DeepCopy(vtkDataObject *src) override;
310 
352  int Simplify(double simplificationThreshold,
353  vtkReebGraphSimplificationMetric *simplificationMetric);
354 
359  void Set(vtkMutableDirectedGraph *g);
360 
361 protected:
362 
363  vtkReebGraph();
364  ~vtkReebGraph() override;
365 
366  class Implementation;
367  Implementation* Storage;
368 
369 private:
370  vtkReebGraph(const vtkReebGraph&) = delete;
371  void operator=(const vtkReebGraph&) = delete;
372 
373 };
374 
375 #endif
vtkMutableDirectedGraph
An editable directed graph.
Definition: vtkMutableDirectedGraph.h:45
vtkIdType
int vtkIdType
Definition: vtkType.h:347
vtkGraph::DeepCopy
void DeepCopy(vtkDataObject *obj) override
Deep copies the data object into this graph.
VTK_REEB_GRAPH
#define VTK_REEB_GRAPH
Definition: vtkType.h:119
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
vtkMutableDirectedGraph.h
vtkReebGraphSimplificationMetric
abstract class for custom Reeb graph simplification metric design.
Definition: vtkReebGraphSimplificationMetric.h:56
vtkReebGraph
Reeb graph computation for PL scalar fields.
Definition: vtkReebGraph.h:131
vtkReebGraph::Storage
Implementation * Storage
Definition: vtkReebGraph.h:366
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkDataSet
abstract class to specify dataset behavior
Definition: vtkDataSet.h:56
vtkPolyData
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:79
vtkMutableDirectedGraph::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:81
vtkReebGraph::GetDataObjectType
int GetDataObjectType() override
Return class name of data type.
Definition: vtkReebGraph.h:148
vtkDataObject
general representation of visualization data
Definition: vtkDataObject.h:58
vtkMutableDirectedGraph::New
static vtkMutableDirectedGraph * New()