TTK
Loading...
Searching...
No Matches
ttkContourForests.h
Go to the documentation of this file.
1
36
37#pragma once
38
39// VTK includes
40#include <vtkNew.h>
41#include <vtkPolyData.h>
42#include <vtkSmartPointer.h>
43#include <vtkUnstructuredGrid.h>
44
45// VTK Module
46#include <ttkContourForestsModule.h>
47
48// vtk wrapper includes
49#include <ttkAlgorithm.h>
50
51// base code includes
52#include <Geometry.h>
53
54#include <ContourForests.h>
55#include <ContourForestsTree.h>
56#include <DeprecatedDataTypes.h>
57
58class TTKCONTOURFORESTS_EXPORT ttkContourForests
59 : public ttkAlgorithm,
60 protected ttk::cf::ContourForests {
61public:
64
65 vtkGetMacro(ForceInputOffsetScalarField, bool);
66 void SetForceInputOffsetScalarField(bool onOff);
67
68 void SetTreeType(int tree);
69
70 void ShowMin(bool state);
71 void ShowMax(bool state);
72 void ShowSaddle1(bool state);
73 void ShowSaddle2(bool state);
74
75 void ShowArc(bool state);
76 void SetArcResolution(int arcResolution);
77 void SetPartitionNumber(int partitionNum);
78 void SetLessPartition(bool l);
79
80 void SetSkeletonSmoothing(double skeletonSmooth);
81
82 void SetSimplificationType(int type);
83
84 void SetSimplificationThreshold(double simplificationThreshold);
85
86protected:
88
89 // VTK Interface //
90 int FillInputPortInformation(int port, vtkInformation *info) override;
91 int FillOutputPortInformation(int port, vtkInformation *info) override;
92 int RequestData(vtkInformation *request,
93 vtkInformationVector **inputVector,
94 vtkInformationVector *outputVector) override;
95 void Modified() override;
96
97 // Base //
98 bool isCoincident(float p1[], double p2[]);
99 bool isCoincident(double p1[], double p2[]);
100
101 // ContourForestsTree //
102 void getTree();
103 void updateTree();
104 ttk::CriticalType getNodeType(ttk::SimplexId id);
105 ttk::CriticalType getNodeType(ttk::SimplexId id,
107 ttk::cf::MergeTree *tree);
108 void getCriticalPoints();
109 void clearTree();
110
111 // Skeleton //
112 void getSkeleton();
113 void clearSkeleton();
114 void getSkeletonNodes();
115 void getSkeletonArcs();
116 int getSkeletonScalars(
117 const std::vector<double> &scalars,
118 std::vector<std::vector<double>> &skeletonScalars) const;
119
120 // Segmentation //
121 void getSegmentation(vtkDataSet *input);
122 void clearSegmentation();
123
124 int sample(unsigned int samplingLevel);
125
126 int computeBarycenters();
127 void computeSkeleton(unsigned int arcRes);
128 void smoothSkeleton(unsigned int skeletonSmoothing);
129 void smooth(const ttk::SimplexId idArc, bool order);
130
131private:
132 // Base //
133 bool isLoaded_{};
134 bool lessPartition_{true};
135 ttk::cf::MergeTree *tree_{};
136 vtkSmartPointer<vtkPolyData> skeletonNodes_{vtkPolyData::New()};
137 vtkSmartPointer<vtkPolyData> skeletonArcs_{vtkPolyData::New()};
138 vtkSmartPointer<vtkDataSet> segmentation_{};
139
140 // Void //
141 vtkNew<vtkUnstructuredGrid> voidUnstructuredGrid_{};
142 vtkNew<vtkPolyData> voidPolyData_{};
143
144 // Configuration //
145 bool ForceInputOffsetScalarField{false};
146 bool varyingMesh_{};
147 bool varyingDataValues_{};
149 bool showMin_{true};
150 bool showMax_{true};
151 bool showSaddle1_{true};
152 bool showSaddle2_{true};
153 bool showArc_{true};
154 unsigned int arcResolution_{20};
155 int partitionNum_{-1};
156 unsigned int skeletonSmoothing_{15};
157 int simplificationType_{};
158 double simplificationThreshold_{};
159 double simplificationThresholdBuffer_{};
160
161 // Computation handles
162 bool toUpdateVertexSoSoffsets_{true};
163 bool toComputeContourTree_{true};
164 bool toUpdateTree_{true};
165 bool toComputeSkeleton_{true};
166 bool toComputeSegmentation_{true};
167
168 // Convenient storage //
169 vtkDataArray *vtkInputScalars_{};
170 double deltaScalar_{};
171 ttk::SimplexId numberOfVertices_{};
172 ttk::Triangulation *triangulation_{};
173 ttk::SimplexId *vertexSoSoffsets_{};
174 std::vector<ttk::SimplexId> criticalPoints_{};
175 std::vector<double> vertexScalars_{};
176
177 // treeType, SuperArc, several vertices list.
178 std::vector<std::vector<std::vector<std::vector<ttk::SimplexId>>>> samples_{};
179 std::vector<std::vector<std::vector<std::vector<double>>>> barycenters_{};
180};
Baseclass of all VTK filters that wrap ttk modules.
virtual int RequestData(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
int FillInputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
int FillOutputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
TTK VTK-filter that efficiently computes the contour tree of scalar data and more (data segmentation,...
static ttkContourForests * New()
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
MergeTree * getTree(const TreeType &tt)
CriticalType
default value for critical index
Definition DataTypes.h:80
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22