TTK
Loading...
Searching...
No Matches
JacobiSet_Template.h
Go to the documentation of this file.
1#pragma once
2
3#include <JacobiSet.h>
4
5template <class dataTypeU, class dataTypeV, typename triangulationType>
6int ttk::JacobiSet::execute(std::vector<std::pair<SimplexId, char>> &jacobiSet,
7 const dataTypeU *const uField,
8 const dataTypeV *const vField,
9 const triangulationType &triangulation,
10 std::vector<char> *isPareto) {
11
12 Timer t;
13
14 // check the consistency of the variables -- to adapt
15#ifndef TTK_ENABLE_KAMIKAZE
16 if((triangulation.isEmpty())) {
17 if(vertexNumber_) {
18 return executeLegacy(jacobiSet, uField, vField);
19 }
20 return -1;
21 }
22 if(!uField)
23 return -2;
24 if(!vField)
25 return -3;
26#endif
27
28 jacobiSet.clear();
29
30 SimplexId edgeNumber = triangulation.getNumberOfEdges();
31
32 std::vector<std::vector<std::pair<SimplexId, char>>> threadedCriticalTypes(
34
35#ifdef TTK_ENABLE_OPENMP
36#pragma omp parallel for num_threads(threadNumber_)
37#endif
38 for(SimplexId i = 0; i < edgeNumber; i++) {
39
40 char type = getCriticalType(i, uField, vField, triangulation);
41
42 if(type != -2) {
43 // -2: regular vertex
44#ifdef TTK_ENABLE_OPENMP
45 const auto tid = omp_get_thread_num();
46#else
47 const auto tid = 0;
48#endif // TTK_ENABLE_OPENMP
49 threadedCriticalTypes[tid].emplace_back(i, type);
50 }
51 }
52
53 // now merge the threaded lists
54 size_t jacobiSetSize{};
55 for(const auto &vec : threadedCriticalTypes) {
56 jacobiSetSize += vec.size();
57 }
58 jacobiSet.reserve(jacobiSetSize);
59
60 for(SimplexId i = 0; i < threadNumber_; i++) {
61 for(size_t j = 0; j < threadedCriticalTypes[i].size(); j++) {
62 jacobiSet.push_back(threadedCriticalTypes[i][j]);
63 }
64 }
65
66 SimplexId minimumNumber = 0, saddleNumber = 0, maximumNumber = 0,
67 monkeySaddleNumber = 0;
68
69 for(size_t i = 0; i < jacobiSet.size(); i++) {
70 switch(jacobiSet[i].second) {
71 case 0:
72 minimumNumber++;
73 break;
74 case 1:
75 saddleNumber++;
76 break;
77 case 2:
78 maximumNumber++;
79 break;
80 case -1:
81 monkeySaddleNumber++;
82 break;
83 }
84 }
85
86 if(isPareto) {
87 isPareto->resize(jacobiSet.size(), 0);
88#ifdef TTK_ENABLE_OPENMP
89#pragma omp parallel for num_threads(threadNumber_)
90#endif
91 for(int i = 0; i < (int)jacobiSet.size(); i++) {
92 int edgeId = jacobiSet[i].first;
93 SimplexId vertexId0 = -1, vertexId1 = -1;
94 triangulation.getEdgeVertex(edgeId, 0, vertexId0);
95 triangulation.getEdgeVertex(edgeId, 1, vertexId1);
96 double denominator = vField[vertexId1] - vField[vertexId0];
97 if(fabs(denominator) < Geometry::powIntTen(-DBL_DIG))
98 denominator = 1;
99 if((uField[vertexId1] - uField[vertexId0]) / denominator < 0)
100 (*isPareto)[i] = 1;
101 }
102 }
103
104 this->printMsg(
105 std::vector<std::vector<std::string>>{
106 {"#Minimum edges", std::to_string(minimumNumber)},
107 {"#Saddle edges", std::to_string(saddleNumber)},
108 {"#Maximum edges", std::to_string(maximumNumber)},
109 {"#Multi-saddle edges", std::to_string(monkeySaddleNumber)}},
111
112 this->printMsg(
113 "Data-set processed", 1.0, t.getElapsedTime(), this->threadNumber_);
114 this->printMsg("Jacobi edge rate: "
115 + std::to_string(100.0 * jacobiSet.size() / (double)edgeNumber)
116 + "%");
117
118 return 0;
119}
120
121template <class dataTypeU, class dataTypeV>
123 std::vector<std::pair<SimplexId, char>> &jacobiSet,
124 const dataTypeU *const uField,
125 const dataTypeV *const vField) {
126
127 Timer t;
128
129 this->printWrn("Using legacy implementation...");
130
131 // check the consistency of the variables -- to adapt
132#ifndef TTK_ENABLE_KAMIKAZE
133 if(!vertexNumber_)
134 return -1;
135 if(!uField)
136 return -2;
137 if(!vField)
138 return -3;
139 if(!edgeList_)
140 return -4;
141 if(!edgeFanLinkEdgeLists_)
142 return -5;
143 if(!edgeFans_)
144 return -6;
145 if(!sosOffsetsU_)
146 return -7;
147#endif
148
149 SimplexId count = 0;
150
151 jacobiSet.clear();
152
153 // distance fields (not really memory efficient)
154 // for each thread
155 // for each vertex: distance field map
156 std::vector<std::vector<double>> threadedDistanceField(threadNumber_);
157 for(size_t i = 0; i < threadedDistanceField.size(); i++) {
158 threadedDistanceField[i].resize(vertexNumber_);
159 }
160
161 std::vector<ScalarFieldCriticalPoints> threadedCriticalPoints(threadNumber_);
162 for(ThreadId i = 0; i < threadNumber_; i++) {
163 threadedCriticalPoints[i].setDomainDimension(2);
164 threadedCriticalPoints[i].setVertexNumber(vertexNumber_);
165 }
166
167 std::vector<std::vector<std::pair<SimplexId, char>>> threadedCriticalTypes(
168 threadNumber_);
169
170#ifdef TTK_ENABLE_OPENMP
171#pragma omp parallel for num_threads(threadNumber_)
172#endif
173 for(size_t i = 0; i < edgeList_->size(); i++) {
174
175 // avoid any processing if the abort signal is sent
176 if((!wrapper_) || ((wrapper_) && (!wrapper_->needsToAbort()))) {
177
178 ThreadId threadId = 0;
179#ifdef TTK_ENABLE_OPENMP
180 threadId = omp_get_thread_num();
181#endif
182
183 // processing here!
184 SimplexId pivotVertexId = (*edgeList_)[i].first;
185 SimplexId otherExtremityId = (*edgeList_)[i].second;
186
187 // A) compute the distance field
188 double projectedPivotVertex[2];
189 projectedPivotVertex[0] = uField[pivotVertexId];
190 projectedPivotVertex[1] = vField[pivotVertexId];
191
192 double projectedOtherVertex[2];
193 projectedOtherVertex[0] = uField[otherExtremityId];
194 projectedOtherVertex[1] = vField[otherExtremityId];
195
196 double rangeEdge[2];
197 rangeEdge[0] = projectedOtherVertex[0] - projectedPivotVertex[0];
198 rangeEdge[1] = projectedOtherVertex[1] - projectedPivotVertex[1];
199
200 double rangeNormal[2];
201 rangeNormal[0] = -rangeEdge[1];
202 rangeNormal[1] = rangeEdge[0];
203
204 for(size_t j = 0; j < (*edgeFans_)[i].size() / 4; j++) {
205 for(int k = 0; k < 3; k++) {
206
207 SimplexId vertexId = (*edgeFans_)[i][j * 4 + 1 + k];
208
209 // we can compute the distance field (in the rage)
210 double projectedVertex[2];
211 projectedVertex[0] = uField[vertexId];
212 projectedVertex[1] = vField[vertexId];
213
214 double vertexRangeEdge[2];
215 vertexRangeEdge[0] = projectedVertex[0] - projectedPivotVertex[0];
216 vertexRangeEdge[1] = projectedVertex[1] - projectedPivotVertex[1];
217
218 // signed distance: linear function of the dot product
219 threadedDistanceField[threadId][vertexId]
220 = vertexRangeEdge[0] * rangeNormal[0]
221 + vertexRangeEdge[1] * rangeNormal[1];
222 }
223 }
224
225 // B) compute critical points
226 // watch out between local and global Ids
227 // what I could do is to translate the ids from global to local
228 // also, lots of things in there can be done out of the loop
229
230 // in the loop
231 char type = threadedCriticalPoints[threadId].getCriticalType(
232 pivotVertexId, sosOffsetsU_, (*edgeFanLinkEdgeLists_)[i]);
233
234 if(type != -2) {
235 // -2: regular vertex
236 threadedCriticalTypes[threadId].push_back(
237 std::pair<SimplexId, char>(i, type));
238 }
239
240 // update the progress bar of the wrapping code -- to adapt
241 if(debugLevel_ > static_cast<int>(debug::Priority::DETAIL)) {
242#ifdef TTK_ENABLE_OPENMP
243#pragma omp critical
244#endif
245 {
246 if((wrapper_) && (!(count % ((vertexNumber_) / 10)))) {
247 wrapper_->updateProgress((count + 1.0) / vertexNumber_);
248 }
249
250 count++;
251 }
252 }
253 }
254 }
255
256 // now merge the threaded lists
257 for(ThreadId i = 0; i < threadNumber_; i++) {
258 for(size_t j = 0; j < threadedCriticalTypes[i].size(); j++) {
259 jacobiSet.push_back(threadedCriticalTypes[i][j]);
260 }
261 }
262
263 SimplexId minimumNumber = 0, saddleNumber = 0, maximumNumber = 0,
264 monkeySaddleNumber = 0;
265
266 for(size_t i = 0; i < jacobiSet.size(); i++) {
267 switch(jacobiSet[i].second) {
268 case 0:
269 minimumNumber++;
270 break;
271 case 1:
272 saddleNumber++;
273 break;
274 case 2:
275 maximumNumber++;
276 break;
277 case -1:
278 monkeySaddleNumber++;
279 break;
280 }
281 }
282
283 this->printMsg(
284 std::vector<std::vector<std::string>>{
285 {"#Minimum edges", std::to_string(minimumNumber)},
286 {"#Saddle edges", std::to_string(saddleNumber)},
287 {"#Maximum edges", std::to_string(maximumNumber)},
288 {"#Multi-saddle edges", std::to_string(monkeySaddleNumber)}},
290
291 this->printMsg(
292 "Data-set processed", 1.0, t.getElapsedTime(), this->threadNumber_);
293 this->printMsg(
294 "Jacobi edge rate: "
295 + std::to_string(100.0 * (jacobiSet.size() / ((double)edgeList_->size())))
296 + "%");
297
298 return 0;
299}
300
301template <class dataTypeU, class dataTypeV, typename triangulationType>
303 const dataTypeU *const uField,
304 const dataTypeV *const vField,
305 const triangulationType &triangulation) {
306
307 SimplexId vertexId0 = -1, vertexId1 = -1;
308 triangulation.getEdgeVertex(edgeId, 0, vertexId0);
309 triangulation.getEdgeVertex(edgeId, 1, vertexId1);
310
311 double projectedPivotVertex[2];
312 projectedPivotVertex[0] = uField[vertexId0];
313 projectedPivotVertex[1] = vField[vertexId0];
314
315 double projectedOtherVertex[2];
316 projectedOtherVertex[0] = uField[vertexId1];
317 projectedOtherVertex[1] = vField[vertexId1];
318
319 double rangeEdge[2];
320 rangeEdge[0] = projectedOtherVertex[0] - projectedPivotVertex[0];
321 rangeEdge[1] = projectedOtherVertex[1] - projectedPivotVertex[1];
322
323 double rangeNormal[2];
324 rangeNormal[0] = -rangeEdge[1];
325 rangeNormal[1] = rangeEdge[0];
326
327 SimplexId starNumber = triangulation.getEdgeStarNumber(edgeId);
328 std::vector<SimplexId> lowerNeighbors, upperNeighbors;
329
330 SimplexId neighborNumber = 0;
331
332 for(SimplexId i = 0; i < starNumber; i++) {
333
334 SimplexId tetId = -1;
335 triangulation.getEdgeStar(edgeId, i, tetId);
336
337 SimplexId vertexNumber = triangulation.getCellVertexNumber(tetId);
338 for(SimplexId j = 0; j < vertexNumber; j++) {
339 SimplexId vertexId = -1;
340 triangulation.getCellVertex(tetId, j, vertexId);
341
342 if((vertexId != -1) && (vertexId != vertexId0)
343 && (vertexId != vertexId1)) {
344 // new neighbor
345 bool isIn = false;
346 for(size_t k = 0; k < lowerNeighbors.size(); k++) {
347 if(vertexId == lowerNeighbors[k]) {
348 isIn = true;
349 break;
350 }
351 }
352
353 if(!isIn) {
354 for(size_t k = 0; k < upperNeighbors.size(); k++) {
355 if(vertexId == upperNeighbors[k]) {
356 isIn = true;
357 break;
358 }
359 }
360 }
361
362 if(!isIn) {
363 // compute the actual distance field
364 // A) compute the distance field
365 double projectedVertex[2];
366 projectedVertex[0] = uField[vertexId];
367 projectedVertex[1] = vField[vertexId];
368
369 double vertexRangeEdge[2];
370 vertexRangeEdge[0] = projectedVertex[0] - projectedPivotVertex[0];
371 vertexRangeEdge[1] = projectedVertex[1] - projectedPivotVertex[1];
372
373 // signed distance: linear function of the dot product
374 double distance = vertexRangeEdge[0] * rangeNormal[0]
375 + vertexRangeEdge[1] * rangeNormal[1];
376
377 neighborNumber++;
378
379 if(distance < 0) {
380 lowerNeighbors.push_back(vertexId);
381 } else if(distance > 0) {
382 upperNeighbors.push_back(vertexId);
383 } else {
384 // degenerate
385 // compute the distance field out of the offset positions
386 double offsetProjectedPivotVertex[2];
387 offsetProjectedPivotVertex[0] = sosOffsetsU_[vertexId0];
388 offsetProjectedPivotVertex[1]
389 = sosOffsetsV_[vertexId0] * sosOffsetsV_[vertexId0];
390
391 double offsetProjectedOtherVertex[2];
392 offsetProjectedOtherVertex[0] = sosOffsetsU_[vertexId1];
393 offsetProjectedOtherVertex[1]
394 = sosOffsetsV_[vertexId1] * sosOffsetsV_[vertexId1];
395
396 double offsetRangeEdge[2];
397 offsetRangeEdge[0]
398 = offsetProjectedOtherVertex[0] - offsetProjectedPivotVertex[0];
399 offsetRangeEdge[1]
400 = offsetProjectedOtherVertex[1] - offsetProjectedPivotVertex[1];
401
402 double offsetRangeNormal[2];
403 offsetRangeNormal[0] = -offsetRangeEdge[1];
404 offsetRangeNormal[1] = offsetRangeEdge[0];
405
406 projectedVertex[0] = sosOffsetsU_[vertexId];
407 projectedVertex[1]
408 = sosOffsetsV_[vertexId] * sosOffsetsV_[vertexId];
409
410 vertexRangeEdge[0]
411 = projectedVertex[0] - offsetProjectedPivotVertex[0];
412 vertexRangeEdge[1]
413 = projectedVertex[1] - offsetProjectedPivotVertex[1];
414
415 distance = vertexRangeEdge[0] * offsetRangeNormal[0]
416 + vertexRangeEdge[1] * offsetRangeNormal[1];
417
418 if(distance < 0) {
419 lowerNeighbors.push_back(vertexId);
420 } else if(distance > 0) {
421 upperNeighbors.push_back(vertexId);
422 } else {
423 this->printWrn(
424 "Inconsistent (non-bijective?) offsets for vertex #"
425 + std::to_string(vertexId));
426 }
427 }
428 }
429 }
430 }
431 }
432
433 // at this point, we know if each vertex of the edge link is higher or not.
434 if((SimplexId)(lowerNeighbors.size() + upperNeighbors.size())
435 != neighborNumber) {
436 // Inconsistent offsets (cf above error message)
437 return -2;
438 }
439
440 if(lowerNeighbors.empty()) {
441 if(rangeNormal[0] + rangeNormal[1] > 0)
442 // minimum
443 return 0;
444 // else maximum
445 return triangulation.getDimensionality() - 1;
446 }
447 if(upperNeighbors.empty()) {
448 if(rangeNormal[0] + rangeNormal[1] > 0)
449 // maximum
450 return triangulation.getDimensionality() - 1;
451 // else minimum
452 return 0;
453 }
454
455 // let's check the connectivity now
456 std::vector<UnionFind> lowerSeeds(lowerNeighbors.size());
457 std::vector<UnionFind *> lowerList(lowerNeighbors.size());
458 std::vector<UnionFind> upperSeeds(upperNeighbors.size());
459 std::vector<UnionFind *> upperList(upperNeighbors.size());
460
461 for(size_t i = 0; i < lowerSeeds.size(); i++) {
462 lowerList[i] = &(lowerSeeds[i]);
463 }
464 for(size_t i = 0; i < upperSeeds.size(); i++) {
465 upperList[i] = &(upperSeeds[i]);
466 }
467
468 for(SimplexId i = 0; i < starNumber; i++) {
469
470 SimplexId tetId = -1;
471 triangulation.getEdgeStar(edgeId, i, tetId);
472
473 SimplexId vertexNumber = triangulation.getCellVertexNumber(tetId);
474 for(SimplexId j = 0; j < vertexNumber; j++) {
475 SimplexId edgeVertexId0 = -1;
476 triangulation.getCellVertex(tetId, j, edgeVertexId0);
477 if((edgeVertexId0 != vertexId0) && (edgeVertexId0 != vertexId1)) {
478 for(SimplexId k = j + 1; k < vertexNumber; k++) {
479 SimplexId edgeVertexId1 = -1;
480 triangulation.getCellVertex(tetId, k, edgeVertexId1);
481 if((edgeVertexId1 != vertexId0) && (edgeVertexId1 != vertexId1)) {
482 // processing the edge (edgeVertexId0, edgeVertexId1)
483
484 // we need to find out if they're lower or not
485 bool lower0 = false;
486 for(size_t l = 0; l < lowerNeighbors.size(); l++) {
487 if(lowerNeighbors[l] == edgeVertexId0) {
488 lower0 = true;
489 break;
490 }
491 }
492 bool lower1 = false;
493 for(size_t l = 0; l < lowerNeighbors.size(); l++) {
494 if(lowerNeighbors[l] == edgeVertexId1) {
495 lower1 = true;
496 break;
497 }
498 }
499
500 auto *neighbors = &lowerNeighbors;
501 auto *seeds = &lowerList;
502
503 if(!lower0) {
504 neighbors = &upperNeighbors;
505 seeds = &upperList;
506 }
507
508 if(lower0 == lower1) {
509 // connect their union-find sets!
510 SimplexId lowerId0 = -1, lowerId1 = -1;
511 for(size_t l = 0; l < neighbors->size(); l++) {
512 if((*neighbors)[l] == edgeVertexId0) {
513 lowerId0 = l;
514 }
515 if((*neighbors)[l] == edgeVertexId1) {
516 lowerId1 = l;
517 }
518 }
519
520 if((lowerId0 != -1) && (lowerId1 != -1)) {
521 (*seeds)[lowerId0] = UnionFind::makeUnion(
522 (*seeds)[lowerId0], (*seeds)[lowerId1]);
523 (*seeds)[lowerId1] = (*seeds)[lowerId0];
524 }
525 }
526
527 break;
528 }
529 }
530 }
531 }
532 }
533
534 // update the UF if necessary
535 for(size_t i = 0; i < lowerList.size(); i++)
536 lowerList[i] = lowerList[i]->find();
537 for(size_t i = 0; i < upperList.size(); i++)
538 upperList[i] = upperList[i]->find();
539
540 {
541 std::sort(lowerList.begin(), lowerList.end());
542 const auto it = std::unique(lowerList.begin(), lowerList.end());
543 lowerList.erase(it, lowerList.end());
544 }
545 {
546 std::sort(upperList.begin(), upperList.end());
547 const auto it = std::unique(upperList.begin(), upperList.end());
548 upperList.erase(it, upperList.end());
549 }
550
551 if((upperList.size() == 1) && (lowerList.size() == 1))
552 return -2;
553
554 return 1;
555}
556
557template <class dataTypeU, class dataTypeV>
558int ttk::JacobiSet::perturb(const dataTypeU *const uField,
559 const dataTypeV *const vField,
560 const dataTypeU uEpsilon,
561 const dataTypeV vEpsilon) const {
562
563#ifndef TTK_ENABLE_KAMIKAZE
564 if(!uField)
565 return -1;
566 if(!vField)
567 return -2;
568 if(!vertexNumber_)
569 return -3;
570#endif
571
572#ifdef TTK_ENABLE_OPENMP
573#pragma omp parallel for num_threads(threadNumber_)
574#endif
575 for(SimplexId i = 0; i < vertexNumber_; i++) {
576 // simulation of simplicity in 2 dimensions, need to use degree 2 polynoms
577
578 uField[i] += i * uEpsilon;
579 vField[i] += (i * vEpsilon) * (i * vEpsilon);
580 }
581
582 return 0;
583}
int threadNumber_
Definition: BaseClass.h:95
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 executeLegacy(std::vector< std::pair< SimplexId, char > > &jacobiSet, const dataTypeU *const uField, const dataTypeV *const vField)
SimplexId vertexNumber_
Definition: JacobiSet.h:131
char getCriticalType(const SimplexId &edgeId, const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation)
int perturb(const dataTypeU *const uField, const dataTypeV *const vField, const dataTypeU uEpsilon=Geometry::powIntTen(-DBL_DIG), const dataTypeV vEpsilon=Geometry::powIntTen(-DBL_DIG)) const
int execute(std::vector< std::pair< SimplexId, char > > &jacobiSet, const dataTypeU *const uField, const dataTypeV *const vField, const triangulationType &triangulation, std::vector< char > *isPareto=nullptr)
double getElapsedTime()
Definition: Timer.h:15
static UnionFind * makeUnion(UnionFind *uf0, UnionFind *uf1)
Definition: UnionFind.h:40
T powIntTen(const int n)
Compute the nth power of ten.
Definition: Geometry.h:360
int ThreadId
Identifier type for threads (i.e. with OpenMP).
Definition: DataTypes.h:26
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22
printMsg(debug::output::GREEN+"                           "+debug::output::ENDCOLOR+debug::output::GREEN+"▒"+debug::output::ENDCOLOR+debug::output::GREEN+"▒▒▒▒▒▒▒▒▒▒▒▒▒░"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)