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