TTK
Loading...
Searching...
No Matches
ttkMergeTreeClustering.h
Go to the documentation of this file.
1
42
43#pragma once
44
45// VTK Module
46#include <ttkMergeTreeClusteringModule.h>
47
48// VTK Includes
49#include <ttkAlgorithm.h>
50#include <vtkMultiBlockDataSet.h>
51#include <vtkSmartPointer.h>
52#include <vtkUnstructuredGrid.h>
53
54// TTK Base Includes
55#include <MergeTreeBarycenter.h>
56#include <MergeTreeClustering.h>
57#include <MergeTreeDistance.h>
58
59class TTKMERGETREECLUSTERING_EXPORT ttkMergeTreeClustering
60 : public ttkAlgorithm // we inherit from the generic ttkAlgorithm class
61//, protected ttk::MergeTreeDistance // and we inherit from the base class
62{
63private:
68 // Input Options
69 bool Epsilon1UseFarthestSaddle = false;
70 double EpsilonTree1 = 5.;
71 double EpsilonTree2 = 5.;
72 double Epsilon2Tree1 = 95.;
73 double Epsilon2Tree2 = 95.;
74 double Epsilon3Tree1 = 90.;
75 double Epsilon3Tree2 = 90.;
76 double PersistenceThreshold = 0.;
77 bool DeleteMultiPersPairs = false;
78 bool UseMinMaxPair = true;
79 bool IsPersistenceDiagram = false;
80
81 // Execution Options
82 int Backend = 0;
83 double Alpha = 0.5;
84 int AssignmentSolver = 0;
85 bool BranchDecomposition = true;
86 bool NormalizedWasserstein = true;
87 bool KeepSubtree = false;
88 bool oldBD = BranchDecomposition;
89 bool oldNW = NormalizedWasserstein;
90 bool oldKS = KeepSubtree;
91 double JoinSplitMixtureCoefficient = 0.5;
92 bool ComputeBarycenter = false;
93 unsigned int NumberOfBarycenters = 1;
94 double BarycenterSizeLimitPercent = 0.0;
95 bool Deterministic = false;
96 int pathMetric = 0;
97 int branchMetric = 0;
98 int baseModule = 0;
99 double NonMatchingWeight = 1.0;
100
101 // Output Options
102 bool OutputTrees = true;
103 bool OutputSegmentation = false;
104 bool PlanarLayout = false;
105 bool BranchDecompositionPlanarLayout = false;
106 double BranchSpacing = 1.;
107 bool RescaleTreesIndividually = false;
108 double DimensionSpacing = 1.;
109 int DimensionToShift = 0;
110 double ImportantPairs = 50.;
111 int MaximumImportantPairs = 0;
112 int MinimumImportantPairs = 0;
113 double ImportantPairsSpacing = 1.;
114 double NonImportantPairsSpacing = 1.;
115 double NonImportantPairsProximity = 0.05;
116 bool BarycenterPositionAlpha = false;
117 std::string ExcludeImportantPairsLower = "";
118 std::string ExcludeImportantPairsHigher = "";
119
120 // ----------------------
121 // Data for visualization
122 // ----------------------
123 // Trees
124 std::vector<ttk::ftm::MergeTree<double>> intermediateSTrees,
125 intermediateSTrees2;
126 std::vector<vtkUnstructuredGrid *> treesNodes, treesNodes2;
127 std::vector<vtkUnstructuredGrid *> treesArcs, treesArcs2;
128 std::vector<vtkDataSet *> treesSegmentation, treesSegmentation2;
129
130 // Matching
131 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>
132 outputMatching;
133 std::vector<std::vector<
134 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>>>
135 outputMatchingBarycenter, outputMatchingBarycenter2;
136
137 // Barycenter
138 std::vector<ttk::ftm::MergeTree<double>> barycentersS, barycentersS2;
139 std::vector<int> clusteringAssignment;
140
141 // Node correspondence
142 std::vector<std::vector<int>> trees1NodeCorrMesh, trees2NodeCorrMesh;
143
144 // Distances
145 std::vector<double> finalDistances;
146
147 void setDataVisualization(int numInputs, int numInputs2) {
148 // Trees
149 intermediateSTrees = std::vector<ttk::ftm::MergeTree<double>>(numInputs);
150 intermediateSTrees2 = std::vector<ttk::ftm::MergeTree<double>>(numInputs2);
151 treesNodes = std::vector<vtkUnstructuredGrid *>(numInputs);
152 treesNodes2 = std::vector<vtkUnstructuredGrid *>(numInputs2);
153 treesArcs = std::vector<vtkUnstructuredGrid *>(numInputs);
154 treesArcs2 = std::vector<vtkUnstructuredGrid *>(numInputs2);
155 treesSegmentation = std::vector<vtkDataSet *>(numInputs);
156 treesSegmentation2 = std::vector<vtkDataSet *>(numInputs2);
157
158 // Matching
159 outputMatchingBarycenter = std::vector<std::vector<
160 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>>>(
161 NumberOfBarycenters,
162 std::vector<
163 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>>(
164 numInputs));
165
166 outputMatchingBarycenter2 = std::vector<std::vector<
167 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>>>(
168 NumberOfBarycenters,
169 std::vector<
170 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>>(
171 numInputs2));
172
173 // Barycenter
174 barycentersS
175 = std::vector<ttk::ftm::MergeTree<double>>(NumberOfBarycenters);
176 clusteringAssignment = std::vector<int>(numInputs, 0);
177 }
178
179 void resetDataVisualization() {
180 setDataVisualization(0, 0);
181 trees1NodeCorrMesh = std::vector<std::vector<int>>();
182 trees2NodeCorrMesh = std::vector<std::vector<int>>();
183 finalDistances = std::vector<double>();
184 }
185
186 bool isDataVisualizationFilled() {
187 return trees1NodeCorrMesh.size() != 0 and finalDistances.size() != 0;
188 }
189
190public:
195 // Input Options
196 void SetEpsilon1UseFarthestSaddle(bool epsilon1UseFarthestSaddle) {
197 Epsilon1UseFarthestSaddle = epsilon1UseFarthestSaddle;
198 Modified();
199 resetDataVisualization();
200 }
201 vtkGetMacro(Epsilon1UseFarthestSaddle, bool);
202
203 void SetEpsilonTree1(double epsilonTree1) {
204 EpsilonTree1 = epsilonTree1;
205 Modified();
206 resetDataVisualization();
207 }
208 vtkGetMacro(EpsilonTree1, double);
209
210 void SetEpsilon2Tree1(double epsilon2Tree1) {
211 Epsilon2Tree1 = epsilon2Tree1;
212 Modified();
213 resetDataVisualization();
214 }
215 vtkGetMacro(Epsilon2Tree1, double);
216
217 void SetEpsilon3Tree1(double epsilon3Tree1) {
218 Epsilon3Tree1 = epsilon3Tree1;
219 Modified();
220 resetDataVisualization();
221 }
222 vtkGetMacro(Epsilon3Tree1, double);
223
224 void SetPersistenceThreshold(double persistenceThreshold) {
225 PersistenceThreshold = persistenceThreshold;
226 Modified();
227 resetDataVisualization();
228 }
229 vtkGetMacro(PersistenceThreshold, double);
230
231 void SetUseMinMaxPair(bool useMinMaxPair) {
232 UseMinMaxPair = useMinMaxPair;
233 Modified();
234 resetDataVisualization();
235 }
236 vtkGetMacro(UseMinMaxPair, bool);
237
238 void SetDeleteMultiPersPairs(bool deleteMultiPersPairs) {
239 DeleteMultiPersPairs = deleteMultiPersPairs;
240 Modified();
241 resetDataVisualization();
242 }
243 vtkGetMacro(DeleteMultiPersPairs, bool);
244
245 // Execution Options
246 void SetBackend(int newBackend) {
247 if(Backend == 2) { // Custom
248 oldBD = BranchDecomposition;
249 oldNW = NormalizedWasserstein;
250 oldKS = KeepSubtree;
251 }
252 if(newBackend == 2) { // Custom
253 BranchDecomposition = oldBD;
254 NormalizedWasserstein = oldNW;
255 KeepSubtree = oldKS;
256 }
257 Backend = newBackend;
258 Modified();
259 resetDataVisualization();
260 }
261 vtkGetMacro(Backend, int);
262
263 void SetAlpha(double alpha) {
264 Alpha = 1 - alpha;
265 Alpha = std::min(1 - 1e-6, Alpha);
266 Alpha = std::max(1e-6, Alpha);
267 Modified();
268 resetDataVisualization();
269 }
270 vtkGetMacro(Alpha, double);
271
272 void SetAssignmentSolver(int assignmentSolver) {
273 AssignmentSolver = assignmentSolver;
274 Modified();
275 resetDataVisualization();
276 }
277 vtkGetMacro(AssignmentSolver, int);
278
279 void SetBranchDecomposition(bool branchDecomposition) {
280 BranchDecomposition = branchDecomposition;
281 Modified();
282 resetDataVisualization();
283 }
284 vtkGetMacro(BranchDecomposition, bool);
285
286 void SetDeterministic(bool deterministic) {
287 Deterministic = deterministic;
288 Modified();
289 resetDataVisualization();
290 }
291 vtkGetMacro(Deterministic, bool);
292
293 void SetNormalizedWasserstein(bool normalizedWasserstein) {
294 NormalizedWasserstein = normalizedWasserstein;
295 Modified();
296 resetDataVisualization();
297 }
298 vtkGetMacro(NormalizedWasserstein, bool);
299
300 void SetKeepSubtree(bool keepSubtree) {
301 KeepSubtree = keepSubtree;
302 Modified();
303 resetDataVisualization();
304 }
305 vtkGetMacro(KeepSubtree, bool);
306
307 void SetJoinSplitMixtureCoefficient(double joinSplitMixtureCoefficient) {
308 JoinSplitMixtureCoefficient = joinSplitMixtureCoefficient;
309 Modified();
310 resetDataVisualization();
311 }
312 vtkGetMacro(JoinSplitMixtureCoefficient, double);
313
314 void SetComputeBarycenter(bool computeBarycenter) {
315 ComputeBarycenter = computeBarycenter;
316 Modified();
317 resetDataVisualization();
318 }
319 vtkGetMacro(ComputeBarycenter, bool);
320
321 void SetNumberOfBarycenters(unsigned int numberOfBarycenters) {
322 NumberOfBarycenters = numberOfBarycenters;
323 Modified();
324 resetDataVisualization();
325 }
326 vtkGetMacro(NumberOfBarycenters, unsigned int);
327
328 void SetBarycenterSizeLimitPercent(double percent) {
329 BarycenterSizeLimitPercent = percent;
330 Modified();
331 resetDataVisualization();
332 }
333 vtkGetMacro(BarycenterSizeLimitPercent, double);
334
335 void SetBranchMetric(int m) {
336 branchMetric = m;
337 Modified();
338 }
339
340 void SetPathMetric(int m) {
341 pathMetric = m;
342 Modified();
343 }
344
345 void SetNonMatchingWeight(double weight) {
346 NonMatchingWeight = weight;
347 Modified();
348 resetDataVisualization();
349 }
350 vtkGetMacro(NonMatchingWeight, double);
351
352 // Output Options
353 vtkSetMacro(BarycenterPositionAlpha, bool);
354 vtkGetMacro(BarycenterPositionAlpha, bool);
355
356 vtkSetMacro(OutputTrees, bool);
357 vtkGetMacro(OutputTrees, bool);
358
359 vtkSetMacro(OutputSegmentation, bool);
360 vtkGetMacro(OutputSegmentation, bool);
361
362 vtkSetMacro(PlanarLayout, bool);
363 vtkGetMacro(PlanarLayout, bool);
364
365 vtkSetMacro(BranchDecompositionPlanarLayout, bool);
366 vtkGetMacro(BranchDecompositionPlanarLayout, bool);
367
368 vtkSetMacro(BranchSpacing, double);
369 vtkGetMacro(BranchSpacing, double);
370
371 vtkSetMacro(RescaleTreesIndividually, bool);
372 vtkGetMacro(RescaleTreesIndividually, bool);
373
374 vtkSetMacro(DimensionSpacing, double);
375 vtkGetMacro(DimensionSpacing, double);
376
377 vtkSetMacro(DimensionToShift, int);
378 vtkGetMacro(DimensionToShift, int);
379
380 vtkSetMacro(ImportantPairs, double);
381 vtkGetMacro(ImportantPairs, double);
382
383 vtkSetMacro(MaximumImportantPairs, int);
384 vtkGetMacro(MaximumImportantPairs, int);
385
386 vtkSetMacro(MinimumImportantPairs, int);
387 vtkGetMacro(MinimumImportantPairs, int);
388
389 vtkSetMacro(ImportantPairsSpacing, double);
390 vtkGetMacro(ImportantPairsSpacing, double);
391
392 vtkSetMacro(NonImportantPairsSpacing, double);
393 vtkGetMacro(NonImportantPairsSpacing, double);
394
395 vtkSetMacro(NonImportantPairsProximity, double);
396 vtkGetMacro(NonImportantPairsProximity, double);
397
398 vtkSetMacro(ExcludeImportantPairsLower, const std::string &);
399 vtkGetMacro(ExcludeImportantPairsLower, std::string);
400
401 vtkSetMacro(ExcludeImportantPairsHigher, const std::string &);
402 vtkGetMacro(ExcludeImportantPairsHigher, std::string);
403
410
411protected:
418
423 int FillInputPortInformation(int port, vtkInformation *info) override;
424
429 int FillOutputPortInformation(int port, vtkInformation *info) override;
430
435 int RequestData(vtkInformation *request,
436 vtkInformationVector **inputVector,
437 vtkInformationVector *outputVector) override;
438
439 template <class dataType>
440 int run(vtkInformationVector *outputVector,
441 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
442 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2);
443
444 template <class dataType>
446 vtkInformationVector *outputVector,
447 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
448 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2);
449
450 template <class dataType>
452 vtkInformationVector *outputVector,
453 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
454 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2);
455};
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 wraps the ttk::MergeTreeClustering module.
void SetBarycenterSizeLimitPercent(double percent)
void SetKeepSubtree(bool keepSubtree)
void SetDeterministic(bool deterministic)
void SetComputeBarycenter(bool computeBarycenter)
void SetAssignmentSolver(int assignmentSolver)
void SetEpsilon1UseFarthestSaddle(bool epsilon1UseFarthestSaddle)
void SetJoinSplitMixtureCoefficient(double joinSplitMixtureCoefficient)
void SetEpsilonTree1(double epsilonTree1)
void SetNumberOfBarycenters(unsigned int numberOfBarycenters)
void SetBackend(int newBackend)
int runCompute(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
void SetUseMinMaxPair(bool useMinMaxPair)
static ttkMergeTreeClustering * New()
void SetEpsilon2Tree1(double epsilon2Tree1)
void SetDeleteMultiPersPairs(bool deleteMultiPersPairs)
void SetEpsilon3Tree1(double epsilon3Tree1)
int runOutput(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
void SetBranchDecomposition(bool branchDecomposition)
void SetPersistenceThreshold(double persistenceThreshold)
void SetNonMatchingWeight(double weight)
~ttkMergeTreeClustering() override
void SetNormalizedWasserstein(bool normalizedWasserstein)