TTK
Loading...
Searching...
No Matches
ttkTrackingFromPersistenceDiagrams.cpp
Go to the documentation of this file.
2#include <ttkMacros.h>
4
6
8 SetNumberOfInputPorts(1);
9 SetNumberOfOutputPorts(1);
10}
11
13 int port, vtkInformation *info) {
14 if(port == 0) {
15 info->Set(vtkAlgorithm::INPUT_IS_REPEATABLE(), 1);
16 return 1;
17 }
18 return 0;
19}
20
22 int port, vtkInformation *info) {
23 if(port == 0) {
24 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
25 return 1;
26 }
27 return 0;
28}
29
31 const std::vector<ttk::trackingTuple> &trackings,
32 const std::vector<std::vector<ttk::MatchingType>> &outputMatchings,
33 const std::vector<ttk::DiagramType> &inputPersistenceDiagrams,
34 const bool useGeometricSpacing,
35 const double spacing,
36 const bool doPostProc,
37 const std::vector<std::set<int>> &trackingTupleToMerged,
38 vtkPoints *points,
39 vtkUnstructuredGrid *persistenceDiagram,
40 vtkDoubleArray *persistenceScalars,
41 vtkDoubleArray *valueScalars,
42 vtkIntArray *matchingIdScalars,
43 vtkIntArray *lengthScalars,
44 vtkIntArray *timeScalars,
45 vtkIntArray *componentIds,
46 vtkIntArray *pointTypeScalars,
47 const ttk::Debug &dbg) {
49 int currentVertex = 0;
50 for(size_t k = 0; k < trackings.size(); ++k) {
51 const ttk::trackingTuple &tt = trackings[k];
52
53 int const numStart = std::get<0>(tt);
54 // int numEnd = std::get<1>(tt);
55 const std::vector<ttk::SimplexId> &chain = std::get<2>(tt);
56 int const chainLength = chain.size();
57
58 if(chainLength <= 1) {
59 dbg.printErr("Got an unexpected 0-size chain.");
60 return 0;
61 }
62
63 for(int c = 0; c < chainLength - 1; ++c) {
64 const auto &matchings1 = outputMatchings[numStart + c];
65 int const d1id = numStart + c;
66 int const d2id = d1id + 1; // c % 2 == 0 ? d1id + 1 : d1id;
67 const auto &diagram1 = inputPersistenceDiagrams[d1id];
68 const auto &diagram2 = inputPersistenceDiagrams[d2id];
69
70 // Insert segments
71 vtkIdType ids[2];
72 const auto n1 = chain[c];
73 const auto n2 = chain[c + 1];
74
75 // Search for right matching.
76 double cost = 0.0;
77 for(const auto &tup : matchings1) {
78 auto d1id1 = (int)std::get<0>(tup);
79 if(d1id1 == n1) {
80 cost = std::get<2>(tup);
81 break;
82 }
83 }
84
85 const auto &pair0 = diagram1[n1];
86 const auto &pair1 = diagram2[n2];
87
88 double x1, y1, z1, x2, y2, z2;
89
90 auto point1Type1 = pair0.birth.type;
91 auto point1Type2 = pair0.death.type;
92 auto point1Type = point1Type1 == CriticalType::Local_maximum
93 || point1Type2 == CriticalType::Local_maximum
94 ? CriticalType::Local_maximum
95 : point1Type1 == CriticalType::Local_minimum
96 || point1Type2 == CriticalType::Local_minimum
97 ? CriticalType::Local_minimum
98 : point1Type1 == CriticalType::Saddle2
99 || point1Type2 == CriticalType::Saddle2
100 ? CriticalType::Saddle2
101 : CriticalType::Saddle1;
102 bool t11Min = point1Type1 == CriticalType::Local_minimum;
103 bool t11Max = point1Type1 == CriticalType::Local_maximum;
104 bool t12Min = point1Type2 == CriticalType::Local_minimum;
105 bool t12Max = point1Type2 == CriticalType::Local_maximum;
106 bool bothEx1 = (t11Min && t12Max) || (t11Max && t12Min);
107 if(bothEx1) {
108 x1 = t12Max ? pair0.death.coords[0] : pair0.birth.coords[0];
109 y1 = t12Max ? pair0.death.coords[1] : pair0.birth.coords[1];
110 z1 = t12Max ? pair0.death.coords[2] : pair0.birth.coords[2];
111 if(useGeometricSpacing)
112 z1 += spacing * (numStart + c);
113 } else {
114 x1 = t12Max ? pair0.death.coords[0]
115 : t11Min ? pair0.birth.coords[0]
116 : (pair0.death.coords[0] + pair0.birth.coords[0]) / 2.0;
117 y1 = t12Max ? pair0.death.coords[1]
118 : t11Min ? pair0.birth.coords[1]
119 : (pair0.death.coords[1] + pair0.birth.coords[1]) / 2.0;
120 z1 = t12Max ? pair0.death.coords[2]
121 : t11Min ? pair0.birth.coords[2]
122 : (pair0.death.coords[2] + pair0.birth.coords[2]) / 2.0;
123 if(useGeometricSpacing)
124 z1 += spacing * (numStart + c);
125 }
126
127 // Postproc component ids.
128 int cid = k;
129 bool hasMergedFirst = false;
130 if(doPostProc) {
131 const auto &connected = trackingTupleToMerged[k];
132 if(!connected.empty()) {
133 int const min = *(connected.begin());
134 const ttk::trackingTuple &ttt = trackings[min];
135 // int numStart2 = std::get<0>(ttt);
136 int const numEnd2 = std::get<1>(ttt);
137 if((numEnd2 > 0 && numStart + c > numEnd2 + 1) && min < (int)k) {
138 dbg.printMsg("Switched " + std::to_string(k) + " for "
139 + std::to_string(min));
140 cid = min;
141 hasMergedFirst = numStart + c <= numEnd2 + 3;
142 }
143
144 if(hasMergedFirst) {
145 // std::cout << "Has merged first " << std::endl;
146
147 // Replace former first end of the segment with previous ending
148 // segment.
149 const std::vector<ttk::SimplexId> &chain3 = std::get<2>(ttt);
150 const auto nn = chain3.back();
151 const auto &diagramRematch = inputPersistenceDiagrams[numEnd2];
152 const auto &pairN = diagramRematch[nn];
153
154 point1Type1 = pairN.birth.type;
155 point1Type2 = pairN.death.type;
156 point1Type = point1Type1 == CriticalType::Local_maximum
157 || point1Type2 == CriticalType::Local_maximum
158 ? CriticalType::Local_maximum
159 : point1Type1 == CriticalType::Local_minimum
160 || point1Type2 == CriticalType::Local_minimum
161 ? CriticalType::Local_minimum
162 : point1Type1 == CriticalType::Saddle2
163 || point1Type2 == CriticalType::Saddle2
164 ? CriticalType::Saddle2
165 : CriticalType::Saddle1;
166 t11Min = point1Type1 == CriticalType::Local_minimum;
167 t11Max = point1Type1 == CriticalType::Local_maximum;
168 t12Min = point1Type2 == CriticalType::Local_minimum;
169 t12Max = point1Type2 == CriticalType::Local_maximum;
170 bothEx1 = (t11Min && t12Max) || (t11Max && t12Min);
171 // std::cout << "xyz " << x1 << ", " << y1 << ", " << z1 <<
172 // std::endl;
173 if(bothEx1) {
174 x1 = t12Max ? pairN.death.coords[0] : pairN.birth.coords[0];
175 y1 = t12Max ? pairN.death.coords[1] : pairN.birth.coords[1];
176 z1 = t12Max ? pairN.death.coords[2] : pairN.birth.coords[2];
177 if(useGeometricSpacing)
178 z1 += spacing * (numStart + c);
179 } else {
180 x1 = t12Max ? pairN.death.coords[0]
181 : t11Min
182 ? pairN.birth.coords[0]
183 : (pairN.birth.coords[0] + pairN.death.coords[0]) / 2.0;
184 y1 = t12Max ? pairN.death.coords[1]
185 : t11Min
186 ? pairN.birth.coords[1]
187 : (pairN.birth.coords[1] + pairN.death.coords[1]) / 2.0;
188 z1 = t12Max ? pairN.death.coords[2]
189 : t11Min
190 ? pairN.birth.coords[2]
191 : (pairN.birth.coords[2] + pairN.death.coords[2]) / 2.0;
192 if(useGeometricSpacing)
193 z1 += spacing * (numStart + c);
194 }
195 // std::cout << "xyz " << x1 << ", " << y1 << ", " << z1 <<
196 // std::endl;
197 }
198 }
199 }
200
201 points->InsertNextPoint(x1, y1, z1);
202 ids[0] = 2 * currentVertex;
203 pointTypeScalars->InsertTuple1(ids[0], (double)(int)point1Type);
204 timeScalars->InsertTuple1(ids[0], (double)numStart + c);
205 componentIds->InsertTuple1(ids[0], cid);
206
207 const auto point2Type1 = pair1.birth.type;
208 const auto point2Type2 = pair1.death.type;
209 const auto point2Type = point2Type1 == CriticalType::Local_maximum
210 || point2Type2 == CriticalType::Local_maximum
211 ? CriticalType::Local_maximum
212 : point2Type1 == CriticalType::Local_minimum
213 || point2Type2 == CriticalType::Local_minimum
214 ? CriticalType::Local_minimum
215 : point2Type1 == CriticalType::Saddle2
216 || point2Type2 == CriticalType::Saddle2
217 ? CriticalType::Saddle2
218 : CriticalType::Saddle1;
219 bool const t21Ex = point2Type1 == CriticalType::Local_minimum
220 || point2Type1 == CriticalType::Local_maximum;
221 bool const t22Ex = point2Type2 == CriticalType::Local_minimum
222 || point2Type2 == CriticalType::Local_maximum;
223 bool const bothEx2 = t21Ex && t22Ex;
224 if(bothEx2) {
225 x2 = point2Type2 == CriticalType::Local_maximum ? pair1.death.coords[0]
226 : pair1.birth.coords[0];
227 y2 = point2Type2 == CriticalType::Local_maximum ? pair1.death.coords[1]
228 : pair1.birth.coords[1];
229 z2 = point2Type2 == CriticalType::Local_maximum ? pair1.death.coords[2]
230 : pair1.birth.coords[2];
231 if(useGeometricSpacing)
232 z2 += spacing * (numStart + c + 1);
233 } else {
234 x2 = t22Ex ? pair1.death.coords[0]
235 : t21Ex ? pair1.birth.coords[0]
236 : (pair1.birth.coords[0] + pair1.death.coords[0]) / 2;
237 y2 = t22Ex ? pair1.death.coords[1]
238 : t21Ex ? pair1.birth.coords[1]
239 : (pair1.birth.coords[1] + pair1.death.coords[1]) / 2;
240 z2 = t22Ex ? pair1.death.coords[2]
241 : t21Ex ? pair1.birth.coords[2]
242 : (pair1.birth.coords[2] + pair1.death.coords[2]) / 2;
243 if(useGeometricSpacing)
244 z2 += spacing * (numStart + c + 1);
245 }
246 points->InsertNextPoint(x2, y2, z2);
247 ids[1] = 2 * currentVertex + 1;
248 pointTypeScalars->InsertTuple1(ids[1], (double)(int)point2Type);
249 timeScalars->InsertTuple1(ids[1], (double)numStart + c);
250
251 // Postproc component ids.
252 componentIds->InsertTuple1(ids[1], cid);
253
254 persistenceDiagram->InsertNextCell(VTK_LINE, 2, ids);
255
256 persistenceScalars->InsertTuple1(currentVertex, cost);
257 valueScalars->InsertTuple1(
258 currentVertex, (pair0.death.sfValue + pair1.death.sfValue) / 2);
259 matchingIdScalars->InsertTuple1(currentVertex, currentVertex);
260 lengthScalars->InsertTuple1(currentVertex, chainLength);
261
262 currentVertex++;
263 }
264 }
265
266 persistenceDiagram->SetPoints(points);
267 persistenceDiagram->GetCellData()->AddArray(persistenceScalars);
268 persistenceDiagram->GetCellData()->AddArray(valueScalars);
269 persistenceDiagram->GetCellData()->AddArray(matchingIdScalars);
270 persistenceDiagram->GetCellData()->AddArray(lengthScalars);
271 persistenceDiagram->GetPointData()->AddArray(timeScalars);
272 persistenceDiagram->GetPointData()->AddArray(componentIds);
273 persistenceDiagram->GetPointData()->AddArray(pointTypeScalars);
274
275 return 0;
276}
277
279 vtkInformation *ttkNotUsed(request),
280 vtkInformationVector **inputVector,
281 vtkInformationVector *outputVector) {
282
283 // Unified bound fields
284 vtkUnstructuredGrid *mesh = vtkUnstructuredGrid::GetData(outputVector, 0);
285
286 // Number of input files
287 int const numInputs = inputVector[0]->GetNumberOfInformationObjects();
288 this->printMsg("Number of inputs: " + std::to_string(numInputs));
289
290 // Get input data
291 std::vector<vtkUnstructuredGrid *> inputVTUs(numInputs);
292 for(int i = 0; i < numInputs; i++) {
293 inputVTUs[i] = vtkUnstructuredGrid::GetData(inputVector[0], i);
294 }
295
296 std::vector<ttk::DiagramType> inputPersistenceDiagrams(numInputs);
297 std::vector<vtkNew<vtkUnstructuredGrid>> outputDiags(2 * numInputs - 2);
298 std::vector<std::vector<ttk::MatchingType>> outputMatchings(numInputs - 1);
299
300 // Input parameters.
301 double const spacing = Spacing;
302 std::string const algorithm = DistanceAlgorithm;
303 double const tolerance = Tolerance;
304 std::string const wasserstein = WassersteinMetric;
305
306 // Transform inputs into the right structure.
307 for(int i = 0; i < numInputs; ++i) {
308 VTUToDiagram(inputPersistenceDiagrams[i], inputVTUs[i], *this);
309 }
310
311 this->performMatchings(
312 numInputs, inputPersistenceDiagrams, outputMatchings,
313 algorithm, // Not from paraview, from enclosing tracking plugin
314 wasserstein, tolerance, PX, PY, PZ, PS, PE // Coefficients
315 );
316
317 // Get back meshes.
318 // #pragma omp parallel for num_threads(ThreadNumber)
319 for(int i = 0; i < numInputs - 1; ++i) {
320 outputDiags[2 * i + 0]->ShallowCopy(inputVTUs[i]);
321 outputDiags[2 * i + 1]->ShallowCopy(inputVTUs[i + 1]);
322
323 int const status = augmentDiagrams(
324 outputMatchings[i], outputDiags[2 * i + 0], outputDiags[2 * i + 1]);
325 if(status < 0)
326 return -2;
327 }
328
329 for(int i = 0; i < numInputs; ++i) {
330 const auto &grid = outputDiags[i];
331 if(!grid || !grid->GetCellData()
332 || !grid->GetCellData()->GetArray("Persistence")) {
333 this->printErr("Inputs are not persistence diagrams");
334 return 0;
335 }
336
337 // Check if inputs have the same data type and the same number of points
338 if(grid->GetCellData()->GetArray("Persistence")->GetDataType()
339 != outputDiags[0]
340 ->GetCellData()
341 ->GetArray("Persistence")
342 ->GetDataType()) {
343 this->printErr("Inputs of different data types");
344 return 0;
345 }
346 }
347
348 for(int i = 0; i < numInputs - 1; ++i) {
349 const auto &grid1 = outputDiags[i];
350 const auto &grid2 = outputDiags[i + 1];
351 if(i % 2 == 1 && i < numInputs - 1
352 && grid1->GetCellData()->GetNumberOfTuples()
353 != grid2->GetCellData()->GetNumberOfTuples()) {
354 this->printErr("Inconsistent length or order of input diagrams.");
355 return 0;
356 }
357 }
358
359 vtkUnstructuredGrid *outputMesh = vtkUnstructuredGrid::SafeDownCast(mesh);
360
361 vtkNew<vtkPoints> const points{};
362 vtkNew<vtkUnstructuredGrid> const persistenceDiagram{};
363
364 vtkNew<vtkDoubleArray> persistenceScalars{};
365 vtkNew<vtkDoubleArray> valueScalars{};
366 vtkNew<vtkIntArray> matchingIdScalars{};
367 vtkNew<vtkIntArray> lengthScalars{};
368 vtkNew<vtkIntArray> timeScalars{};
369 vtkNew<vtkIntArray> componentIds{};
370 vtkNew<vtkIntArray> pointTypeScalars{};
371 persistenceScalars->SetName("Cost");
372 valueScalars->SetName("Scalar");
373 matchingIdScalars->SetName("MatchingIdentifier");
374 lengthScalars->SetName("ComponentLength");
375 timeScalars->SetName("TimeStep");
376 componentIds->SetName("ConnectedComponentId");
377 pointTypeScalars->SetName("CriticalType");
378
379 // (+ vertex id)
380 std::vector<ttk::trackingTuple>
381 trackingsBase; // structure containing all trajectories
382 this->performTracking(
383 inputPersistenceDiagrams, outputMatchings, trackingsBase);
384
385 std::vector<std::set<int>> trackingTupleToMerged(trackingsBase.size());
386 if(DoPostProc)
387 this->performPostProcess(inputPersistenceDiagrams, trackingsBase,
388 trackingTupleToMerged, PostProcThresh);
389
390 // bool Is3D = true;
391 bool const useGeometricSpacing = UseGeometricSpacing;
392 // auto spacing = (float) Spacing;
393
394 // std::vector<trackingTuple> trackings = *trackingsBase;
395 // Row = iteration number
396 // Col = (id in pd 2, arrival point id)
397
398 // Build mesh.
399 buildMesh(trackingsBase, outputMatchings, inputPersistenceDiagrams,
400 useGeometricSpacing, spacing, DoPostProc, trackingTupleToMerged,
401 points, persistenceDiagram, persistenceScalars, valueScalars,
402 matchingIdScalars, lengthScalars, timeScalars, componentIds,
403 pointTypeScalars, *this);
404
405 outputMesh->ShallowCopy(persistenceDiagram);
406
407 return 1;
408}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
static int buildMesh(const std::vector< ttk::trackingTuple > &trackings, const std::vector< std::vector< ttk::MatchingType > > &outputMatchings, const std::vector< ttk::DiagramType > &inputPersistenceDiagrams, const bool useGeometricSpacing, const double spacing, const bool doPostProc, const std::vector< std::set< int > > &trackingTupleToMerged, vtkPoints *points, vtkUnstructuredGrid *persistenceDiagram, vtkDoubleArray *persistenceScalars, vtkDoubleArray *valueScalars, vtkIntArray *matchingIdScalars, vtkIntArray *lengthScalars, vtkIntArray *timeScalars, vtkIntArray *componentIds, vtkIntArray *pointTypeScalars, const ttk::Debug &dbg)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Minimalist debugging class.
Definition Debug.h:88
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition Debug.h:118
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
int performMatchings(int numInputs, std::vector< ttk::DiagramType > &inputPersistenceDiagrams, std::vector< std::vector< MatchingType > > &outputMatchings, const std::string &algorithm, const std::string &wasserstein, double tolerance, double px, double py, double pz, double ps, double pe)
int performTracking(std::vector< ttk::DiagramType > &allDiagrams, std::vector< std::vector< MatchingType > > &allMatchings, std::vector< trackingTuple > &trackings)
int performPostProcess(const std::vector< ttk::DiagramType > &allDiagrams, const std::vector< trackingTuple > &trackings, std::vector< std::set< int > > &trackingTupleToMerged, const double postProcThresh)
CriticalType
default value for critical index
Definition DataTypes.h:88
std::tuple< int, int, std::vector< SimplexId > > trackingTuple
int augmentDiagrams(const std::vector< ttk::MatchingType > &matchings, vtkUnstructuredGrid *const vtu0, vtkUnstructuredGrid *const vtu1)
int VTUToDiagram(ttk::DiagramType &diagram, vtkUnstructuredGrid *vtu, const ttk::Debug &dbg)
Converts a Persistence Diagram in the VTK Unstructured Grid format (as generated by the ttkPersistenc...
vtkStandardNewMacro(ttkTrackingFromPersistenceDiagrams)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)