TTK
Loading...
Searching...
No Matches
Identifiers.h
Go to the documentation of this file.
1
11
12#pragma once
13
14// ttk common includes
15#include <Debug.h>
16#ifdef TTK_ENABLE_MPI
17#include <Geometry.h>
18#include <KDTree.h>
19#include <map>
20#endif
21
22namespace ttk {
23
24#ifdef TTK_ENABLE_MPI
25 struct Point {
26 double x;
27 double y;
28 double z;
29 ttk::SimplexId localId;
30 };
31
32 struct Response {
34 ttk::LongSimplexId globalId;
35 };
36
37#endif
38
45 class Identifiers : virtual public Debug {
46
47 public:
49
50 ~Identifiers() override = default;
51
52 protected:
57#ifdef TTK_ENABLE_MPI
58 std::unordered_map<ttk::SimplexId, ttk::SimplexId> *vertGtoL_;
59 std::vector<int> *neighbors_;
60 std::map<int, int> neighborToId_;
61 int neighborNumber_;
62 double *bounds_;
63 int *dims_;
64 double *spacing_;
65 MPI_Datatype mpiIdType_;
66 MPI_Datatype mpiResponseType_;
67 MPI_Datatype mpiPointType_;
68 int dimension_{};
69 int hasSentData_{0};
70 int *vertexRankArray_{nullptr};
71 int *cellRankArray_{nullptr};
72 unsigned char *vertGhost_{nullptr};
73 unsigned char *cellGhost_{nullptr};
74 std::vector<std::vector<ttk::SimplexId>> pointsToCells_;
75 float *pointSet_;
76 ttk::LongSimplexId *connectivity_;
77 ttk::LongSimplexId *outdatedGlobalPointIds_{nullptr};
78 ttk::LongSimplexId *outdatedGlobalCellIds_{nullptr};
79 std::map<ttk::LongSimplexId, ttk::SimplexId> vertOutdatedGtoL_;
80 std::map<ttk::LongSimplexId, ttk::SimplexId> cellOutdatedGtoL_;
82#endif
83
84 public:
85 void setVertexNumber(const SimplexId &vertexNumber) {
86 vertexNumber_ = vertexNumber;
87 }
88
89 void setCellNumber(const SimplexId &cellNumber) {
90 cellNumber_ = cellNumber;
91 }
92
93 void setVertexIdentifiers(ttk::LongSimplexId *vertexIdentifiers) {
94 this->vertexIdentifiers_ = vertexIdentifiers;
95 }
96
97 void setCellIdentifiers(ttk::LongSimplexId *cellIdentifiers) {
98 this->cellIdentifiers_ = cellIdentifiers;
99 }
100
101#ifdef TTK_ENABLE_MPI
102
103 inline void setDomainDimension(const int &dimension) {
104 dimension_ = dimension;
105 }
106
107 inline void setPointsToCells(
108 const std::vector<std::vector<ttk::SimplexId>> &pointsToCells) {
109 pointsToCells_ = pointsToCells;
110 }
111
112 void setVertexRankArray(int *vertexRankArray) {
113 this->vertexRankArray_ = vertexRankArray;
114 }
115
116 void setCellRankArray(int *cellRankArray) {
117 this->cellRankArray_ = cellRankArray;
118 }
119
120 void setVertGhost(unsigned char *vertGhost) {
121 this->vertGhost_ = vertGhost;
122 }
123
124 void setCellGhost(unsigned char *cellGhost) {
125 this->cellGhost_ = cellGhost;
126 }
127
128 void setBounds(double *bounds) {
129 this->bounds_ = bounds;
130 }
131
132 void setVertGtoL(
133 std::unordered_map<ttk::SimplexId, ttk::SimplexId> *vertGtoL) {
134 this->vertGtoL_ = vertGtoL;
135 }
136
137 void setSpacing(double *spacing) {
138 this->spacing_ = spacing;
139 }
140
141 void setDims(int *dims) {
142 this->dims_ = dims;
143 }
144
145 void setPointSet(float *pointSet) {
146 pointSet_ = pointSet;
147 }
148
149 void setConnectivity(ttk::LongSimplexId *connectivity) {
150 connectivity_ = connectivity;
151 }
152
153 void setOutdatedGlobalPointIds(ttk::LongSimplexId *outdatedGlobalPointIds) {
154 outdatedGlobalPointIds_ = outdatedGlobalPointIds;
155 }
156
157 void buildKDTree() {
159 kdt_.build(pointSet_, vertexNumber_, 3);
160 }
161
162 void setOutdatedGlobalCellIds(ttk::LongSimplexId *outdatedGlobalCellIds) {
163 outdatedGlobalCellIds_ = outdatedGlobalCellIds;
164 }
165
166 void initializeMPITypes() {
167 ttk::SimplexId id{-1};
168 ttk::LongSimplexId longId{-1};
169 // Initialize id type
170 mpiIdType_ = getMPIType(id);
171
172 // Initialize Point Type
173 MPI_Datatype types[] = {MPI_DOUBLE, MPI_DOUBLE, MPI_DOUBLE, mpiIdType_};
174 int lengths[] = {1, 1, 1, 1};
175 const long int mpi_offsets[]
176 = {offsetof(Point, x), offsetof(Point, y), offsetof(Point, z),
177 offsetof(Point, localId)};
178 MPI_Type_create_struct(4, lengths, mpi_offsets, types, &mpiPointType_);
179 MPI_Type_commit(&mpiPointType_);
180
181 // Initialize Response Type
182 MPI_Datatype typesResponse[] = {mpiIdType_, getMPIType(longId)};
183 int lengthsResponse[] = {1, 1};
184 const long int mpi_offsetsResponse[]
185 = {offsetof(Response, id), offsetof(Response, globalId)};
186 MPI_Type_create_struct(2, lengthsResponse, mpi_offsetsResponse,
187 typesResponse, &mpiResponseType_);
188 MPI_Type_commit(&mpiResponseType_);
189 }
190
191 void inline findPoint(ttk::SimplexId &id, float x, float y, float z) {
192 std::array<float, 3> coordinates = {x, y, z};
193 std::vector<ttk::KDTree<float, std::array<float, 3>> *> neighbours;
194 std::vector<float> costs;
195 kdt_.getKClosest(1, coordinates, neighbours, costs);
196 id = neighbours[0]->id_;
197 }
198
199 void initializeNeighbors(double *boundingBox,
200 std::vector<int> &neighborRanks) {
201 if(neighborRanks.empty()) {
202 preconditionNeighborsUsingBoundingBox(boundingBox, neighborRanks);
203 }
204 neighbors_ = &neighborRanks;
205 neighborToId_.clear();
206 neighborNumber_ = neighbors_->size();
207 for(int i = 0; i < neighborNumber_; i++) {
208 neighborToId_[neighbors_->at(i)] = i;
209 }
210 }
211
226 void inline locatePoints(
227#ifdef TTK_ENABLE_OPENMP
228 std::vector<std::vector<Response>> &locatedSimplices,
229#else
230 std::vector<std::vector<Response>> &ttkNotUsed(locatedSimplices),
231#endif
232 std::vector<Point> &receivedPoints,
233 ttk::SimplexId &recvMessageSize,
234 std::vector<Response> &send_buf) {
235 ttk::LongSimplexId globalId{-1};
236 ttk::SimplexId id{-1};
237 send_buf.clear();
238#ifdef TTK_ENABLE_OPENMP
239#pragma omp parallel for num_threads(threadNumber_) firstprivate(globalId, id) \
240 shared(locatedSimplices)
241#endif
242 for(int n = 0; n < recvMessageSize; n++) {
243 if(bounds_[0] <= receivedPoints[n].x
244 && bounds_[1] >= receivedPoints[n].x
245 && bounds_[2] <= receivedPoints[n].y
246 && bounds_[3] >= receivedPoints[n].y
247 && bounds_[4] <= receivedPoints[n].z
248 && bounds_[5] >= receivedPoints[n].z) {
249 this->findPoint(
250 id, receivedPoints[n].x, receivedPoints[n].y, receivedPoints[n].z);
251 if((vertexRankArray_ != nullptr
252 && vertexRankArray_[id] == ttk::MPIrank_)
253 || (vertGhost_ != nullptr && vertGhost_[id] == 0)) {
254 globalId = vertexIdentifiers_[id];
255 if(globalId >= 0) {
256#ifdef TTK_ENABLE_OPENMP
257 locatedSimplices[omp_get_thread_num()].push_back(
258 Response{receivedPoints[n].localId, globalId});
259#else
260 send_buf.push_back(Response{receivedPoints[n].localId, globalId});
261#endif
262 }
263 }
264 }
265 }
266#ifdef TTK_ENABLE_OPENMP
267 for(int n = 0; n < threadNumber_; n++) {
268 send_buf.insert(send_buf.end(), locatedSimplices[n].begin(),
269 locatedSimplices[n].end());
270 locatedSimplices[n].clear();
271 }
272#endif
273 }
287 void inline identifyPoints(
288#ifdef TTK_ENABLE_OPENMP
289 std::vector<std::vector<Response>> &locatedSimplices,
290#else
291 std::vector<std::vector<Response>> &ttkNotUsed(locatedSimplices),
292#endif
293 std::vector<ttk::SimplexId> &receivedOutdatedGlobalIds,
294 ttk::SimplexId &recvMessageSize,
295 std::vector<Response> &send_buf) {
296 send_buf.clear();
297 ttk::SimplexId globalId{-1};
298 std::map<ttk::LongSimplexId, ttk::SimplexId>::iterator search;
299#ifdef TTK_ENABLE_OPENMP
300#pragma omp parallel for num_threads(threadNumber_)
301#endif
302 for(int n = 0; n < recvMessageSize; n++) {
303 search = vertOutdatedGtoL_.find(receivedOutdatedGlobalIds[n]);
304 if(search != vertOutdatedGtoL_.end()) {
305 globalId = vertexIdentifiers_[search->second];
306 if(globalId >= 0) {
307#ifdef TTK_ENABLE_OPENMP
308 locatedSimplices[omp_get_thread_num()].push_back(
309 Response{receivedOutdatedGlobalIds[n], globalId});
310#else
311 send_buf.push_back(
312 Response{receivedOutdatedGlobalIds[n], globalId});
313#endif
314 }
315 }
316 }
317#ifdef TTK_ENABLE_OPENMP
318 for(int n = 0; n < threadNumber_; n++) {
319 send_buf.insert(send_buf.end(), locatedSimplices[n].begin(),
320 locatedSimplices[n].end());
321 locatedSimplices[n].clear();
322 }
323#endif
324 }
338 void inline locateCells(
339#ifdef TTK_ENABLE_OPENMP
340 std::vector<std::vector<Response>> &locatedSimplices,
341#else
342 std::vector<std::vector<Response>> &ttkNotUsed(locatedSimplices),
343#endif
344 std::vector<ttk::SimplexId> &receivedCells,
345 ttk::SimplexId &recvMessageSize,
346 std::vector<Response> &send_buf) {
347 int id{-1};
348 send_buf.clear();
349 std::unordered_map<ttk::SimplexId, ttk::SimplexId>::iterator search;
350 std::vector<ttk::SimplexId> localPointIds;
351 localPointIds.reserve(dimension_ + 1);
352 size_t expectedSize = static_cast<size_t>(dimension_) + 1;
353#ifdef TTK_ENABLE_OPENMP
354#pragma omp parallel for num_threads(threadNumber_) \
355 firstprivate(search, localPointIds)
356#endif
357 for(int n = 0; n < recvMessageSize; n += dimension_ + 2) {
358 localPointIds.clear();
359 for(int k = 1; k < dimension_ + 2; k++) {
360 search = vertGtoL_->find(receivedCells[n + k]);
361 if(search != vertGtoL_->end()) {
362 localPointIds.push_back(search->second);
363 } else {
364 break;
365 }
366 }
367 if(localPointIds.size() == expectedSize) {
368 bool foundIt = false;
369 int k = 0;
370 int l;
371 int m = 0;
372 while(!foundIt && m < dimension_ + 1) {
373 int size = pointsToCells_[localPointIds[m]].size();
374 k = 0;
375 while(!foundIt && k < size) {
376 l = 0;
377 while(l < dimension_ + 1) {
378 id = connectivity_[pointsToCells_[localPointIds[m]][k]
379 * (dimension_ + 1)
380 + l];
381 auto it = find(localPointIds.begin(), localPointIds.end(), id);
382 if(it == localPointIds.end()) {
383 break;
384 }
385 l++;
386 }
387 if(l == dimension_ + 1) {
388 foundIt = true;
389#ifdef TTK_ENABLE_OPENMP
390 locatedSimplices[omp_get_thread_num()].push_back(Response{
391 receivedCells[n],
392 cellIdentifiers_[pointsToCells_[localPointIds[m]][k]]});
393#else
394 send_buf.push_back(Response{
395 receivedCells[n],
396 cellIdentifiers_[pointsToCells_[localPointIds[m]][k]]});
397#endif
398 }
399 k++;
400 }
401 m++;
402 }
403 }
404 }
405#ifdef TTK_ENABLE_OPENMP
406 for(int n = 0; n < threadNumber_; n++) {
407 send_buf.insert(send_buf.end(), locatedSimplices[n].begin(),
408 locatedSimplices[n].end());
409 locatedSimplices[n].clear();
410 }
411#endif
412 }
426 void inline identifyCells(
427#ifdef TTK_ENABLE_OPENMP
428 std::vector<std::vector<Response>> &locatedSimplices,
429#else
430 std::vector<std::vector<Response>> &ttkNotUsed(locatedSimplices),
431#endif
432 std::vector<ttk::SimplexId> &receivedOutdatedGlobalIds,
433 ttk::SimplexId &recvMessageSize,
434 std::vector<Response> &send_buf) {
435 ttk::LongSimplexId globalId{-1};
436 std::map<ttk::LongSimplexId, ttk::SimplexId>::iterator search;
437 send_buf.clear();
438#ifdef TTK_ENABLE_OPENMP
439#pragma omp parallel for num_threads(threadNumber_)
440#endif
441 for(int n = 0; n < recvMessageSize; n++) {
442 search = cellOutdatedGtoL_.find(receivedOutdatedGlobalIds[n]);
443 if(search != cellOutdatedGtoL_.end()) {
444 globalId = cellIdentifiers_[search->second];
445 if(globalId >= 0) {
446#ifdef TTK_ENABLE_OPENMP
447 locatedSimplices[omp_get_thread_num()].push_back(
448 Response{receivedOutdatedGlobalIds[n], globalId});
449#else
450 send_buf.push_back(
451 Response{receivedOutdatedGlobalIds[n], globalId});
452#endif
453 }
454 }
455 }
456#ifdef TTK_ENABLE_OPENMP
457 for(int n = 0; n < threadNumber_; n++) {
458 send_buf.insert(send_buf.end(), locatedSimplices[n].begin(),
459 locatedSimplices[n].end());
460 locatedSimplices[n].clear();
461 }
462#endif
463 }
464
465 template <typename dataType>
466 void sendToAllNeighbors(std::vector<dataType> &vectorToSend,
467 MPI_Datatype messageType) const {
468 for(int j = 0; j < neighborNumber_; j++) {
469 sendVector<dataType>(vectorToSend, messageType, neighbors_->at(j));
470 }
471 }
472
473 template <typename dataType>
474 void sendToAllNeighborsUsingRankArray(
475 std::vector<std::vector<dataType>> &vectorToSend,
476 MPI_Datatype messageType) {
477
478 for(int j = 0; j < neighborNumber_; j++) {
479 sendVector<dataType>(vectorToSend[neighborToId_[neighbors_->at(j)]],
480 messageType, neighbors_->at(j));
481 }
482 }
483
502 void generateGlobalIds(
503 std::vector<std::vector<Point>> &vertGhostCoordinatesPerRank,
504 std::vector<Point> &vertGhostCoordinates,
505 std::vector<std::vector<ttk::SimplexId>> &vertGhostGlobalIdsPerRank,
506 std::vector<ttk::SimplexId> &vertGhostGlobalIds) {
507 ttk::SimplexId realVertexNumber = vertexNumber_;
508 ttk::SimplexId realCellNumber = cellNumber_;
509 float p[3];
510
511 // Computes the number of vertices owned by the current process
512 // If the vertex is not owned, it will be added to the vector of
513 // ghosts.
514
515 if(vertexRankArray_ != nullptr) {
516 for(ttk::SimplexId i = 0; i < vertexNumber_; i++) {
517 if(vertexRankArray_[i] != ttk::MPIrank_) {
518 realVertexNumber--;
519 if(outdatedGlobalPointIds_ == nullptr) {
520 p[0] = pointSet_[i * 3];
521 p[1] = pointSet_[i * 3 + 1];
522 p[2] = pointSet_[i * 3 + 2];
523 vertGhostCoordinatesPerRank[neighborToId_[vertexRankArray_[i]]]
524 .push_back(Point{p[0], p[1], p[2], i});
525 } else {
526 vertGhostGlobalIdsPerRank[neighborToId_[vertexRankArray_[i]]]
527 .push_back(
528 static_cast<ttk::SimplexId>(outdatedGlobalPointIds_[i]));
529 }
530 }
531 }
532 } else {
533 for(ttk::SimplexId i = 0; i < vertexNumber_; i++) {
534 if(vertGhost_[i] != 0) {
535 realVertexNumber--;
536 if(outdatedGlobalPointIds_ == nullptr) {
537
538 p[0] = pointSet_[i * 3];
539 p[1] = pointSet_[i * 3 + 1];
540 p[2] = pointSet_[i * 3 + 2];
541 vertGhostCoordinates.push_back(Point{p[0], p[1], p[2], i});
542 } else {
543 vertGhostGlobalIds.push_back(
544 static_cast<ttk::SimplexId>(outdatedGlobalPointIds_[i]));
545 }
546 }
547 }
548 }
549
550 // The number of cells owned by the process is computed
551 if(cellRankArray_ != nullptr) {
552 for(ttk::SimplexId i = 0; i < cellNumber_; i++) {
553 if(cellRankArray_[i] != ttk::MPIrank_) {
554 realCellNumber--;
555 }
556 }
557 } else {
558 for(ttk::SimplexId i = 0; i < cellNumber_; i++) {
559 if(cellGhost_[i] != 0) {
560 realCellNumber--;
561 }
562 }
563 }
564
565 ttk::SimplexId vertIndex;
566 ttk::SimplexId cellIndex;
567
568 // Perform exclusive prefix sum to find local offset for vertices and
569 // cells
570 MPI_Exscan(
571 &realVertexNumber, &vertIndex, 1, mpiIdType_, MPI_SUM, ttk::MPIcomm_);
572 MPI_Exscan(
573 &realCellNumber, &cellIndex, 1, mpiIdType_, MPI_SUM, ttk::MPIcomm_);
574
575 // Rank 0 received garbage values, it is replaced by the correct offset
576 // (always 0)
577 if(ttk::MPIrank_ == 0) {
578 vertIndex = 0;
579 cellIndex = 0;
580 }
581
582 // Generate global ids for vertices
583 if(vertexRankArray_ != nullptr) {
584 for(ttk::SimplexId i = 0; i < vertexNumber_; i++) {
585 if(vertexRankArray_[i] == ttk::MPIrank_) {
586 vertexIdentifiers_[i] = vertIndex;
587 (*vertGtoL_)[vertIndex] = i;
588 vertIndex++;
589 }
590 if(outdatedGlobalPointIds_ != nullptr) {
591 vertOutdatedGtoL_[outdatedGlobalPointIds_[i]] = i;
592 }
593 }
594 } else {
595 for(ttk::SimplexId i = 0; i < vertexNumber_; i++) {
596 if(vertGhost_[i] == 0) {
597 vertexIdentifiers_[i] = vertIndex;
598 (*vertGtoL_)[vertIndex] = i;
599 vertIndex++;
600 }
601 if(outdatedGlobalPointIds_ != nullptr) {
602 vertOutdatedGtoL_[outdatedGlobalPointIds_[i]] = i;
603 }
604 }
605 }
606
607 // Generate global ids for cells
608 if(cellRankArray_ != nullptr) {
609 for(ttk::SimplexId i = 0; i < cellNumber_; i++) {
610 if(cellRankArray_[i] == ttk::MPIrank_) {
611 cellIdentifiers_[i] = cellIndex;
612 cellIndex++;
613 }
614 if(outdatedGlobalCellIds_ != nullptr) {
615 cellOutdatedGtoL_[outdatedGlobalCellIds_[i]] = i;
616 }
617 }
618 } else {
619 for(ttk::SimplexId i = 0; i < cellNumber_; i++) {
620 if(cellGhost_[i] == 0) {
621 cellIdentifiers_[i] = cellIndex;
622 cellIndex++;
623 }
624 if(outdatedGlobalCellIds_ != nullptr) {
625 cellOutdatedGtoL_[outdatedGlobalCellIds_[i]] = i;
626 }
627 }
628 }
629 }
630
637 int executePolyData() {
638 vertGtoL_->clear();
639 std::vector<Point> vertGhostCoordinates;
640 std::vector<ttk::SimplexId> vertGhostGlobalIds;
641 std::vector<ttk::SimplexId> cellGhostGlobalVertexIds;
642 std::vector<ttk::SimplexId> cellGhostGlobalIds;
643 std::vector<std::vector<Point>> vertGhostCoordinatesPerRank;
644 std::vector<std::vector<ttk::SimplexId>> vertGhostGlobalIdsPerRank;
645 std::vector<std::vector<ttk::SimplexId>> cellGhostGlobalVertexIdsPerRank;
646 std::vector<std::vector<ttk::SimplexId>> cellGhostGlobalIdsPerRank;
647 // In case RankArray is defined for vertices, a vector needs to be resized
648 // In case outdated global point ids are defined, then it is
649 // vertGhostGlobalIdsPerRank that need to be resized. This vector will
650 // store the outdated global ids of ghost points and their local id.
651 // Otherwise, vertGhostCoordinatesPerRank needs to be resized. This vector
652 // will store the coordinates of a ghost point and its local id.
653 if(vertexRankArray_ != nullptr) {
654 if(outdatedGlobalPointIds_ == nullptr) {
655 vertGhostCoordinatesPerRank.resize(neighborNumber_);
656 } else {
657 vertGhostGlobalIdsPerRank.resize(neighborNumber_);
658 }
659 }
660
661 // Similar as above, but with cells
662 if(cellRankArray_ != nullptr) {
663 if(outdatedGlobalCellIds_ == nullptr) {
664 cellGhostGlobalVertexIdsPerRank.resize(neighborNumber_);
665 } else {
666 cellGhostGlobalIdsPerRank.resize(neighborNumber_);
667 }
668 }
669
670 this->generateGlobalIds(vertGhostCoordinatesPerRank, vertGhostCoordinates,
671 vertGhostGlobalIdsPerRank, vertGhostGlobalIds);
672
673 // Start exchange of data to identify the global ids of ghost simplices
674 // records whether the process has sent its ghost simplices
675 hasSentData_ = 0;
676 ttk::SimplexId recvMessageSize = 0;
677 std::vector<Point> receivedPoints;
678 std::vector<ttk::SimplexId> receivedIds;
679 std::vector<Response> receivedResponse;
680 std::vector<Response> send_buf;
681 std::vector<std::vector<Response>> locatedSimplices(threadNumber_);
682 // For each neighbor, a process will receive the ghost points of all its
683 // neighbors or, if it is its turn, will receive the ghost points of all
684 // its neighbors
685 for(int i = 0; i < neighborNumber_ + 1; i++) {
686 if((i == neighborNumber_ && !hasSentData_)
687 || (!hasSentData_ && ttk::MPIrank_ < neighbors_->at(i))) {
688 // it is the turn of the current process to send its ghosts
689 if(outdatedGlobalPointIds_ == nullptr) {
690 if(vertexRankArray_ == nullptr) {
691 sendToAllNeighbors<Point>(vertGhostCoordinates, mpiPointType_);
692 } else {
693 sendToAllNeighborsUsingRankArray<Point>(
694 vertGhostCoordinatesPerRank, mpiPointType_);
695 }
696 } else {
697 if(vertexRankArray_ == nullptr) {
698 sendToAllNeighbors<ttk::SimplexId>(
699 vertGhostGlobalIds, mpiIdType_);
700 } else {
701 sendToAllNeighborsUsingRankArray<ttk::SimplexId>(
702 vertGhostGlobalIdsPerRank, mpiIdType_);
703 }
704 }
705 // The current process receives the responses of all the processes
706 // and associates its ghosts with the right global identifier
707 for(int j = 0; j < neighborNumber_; j++) {
708 recvVector<Response>(receivedResponse, recvMessageSize,
709 mpiResponseType_, neighbors_->at(j));
710 if(outdatedGlobalPointIds_ == nullptr) {
711 for(int n = 0; n < recvMessageSize; n++) {
712 vertexIdentifiers_[receivedResponse[n].id]
713 = receivedResponse[n].globalId;
714 (*vertGtoL_)[receivedResponse[n].globalId]
715 = receivedResponse[n].id;
716 }
717 } else {
718 for(int n = 0; n < recvMessageSize; n++) {
719 vertexIdentifiers_[vertOutdatedGtoL_[receivedResponse[n].id]]
720 = receivedResponse[n].globalId;
721 (*vertGtoL_)[receivedResponse[n].globalId]
722 = receivedResponse[n].id;
723 }
724 }
725 }
726 hasSentData_ = 1;
727 } else {
728 // It is not the turn of the current process to send its data.
729 if(outdatedGlobalPointIds_ == nullptr) {
730 recvVector<Point>(receivedPoints, recvMessageSize, mpiPointType_,
731 neighbors_->at(i - hasSentData_));
732 // Point coordinates are matched to a local point using a kd-tree
733 // and its global id is added to the vector send_buf
734 locatePoints(
735 locatedSimplices, receivedPoints, recvMessageSize, send_buf);
736 } else {
737 recvVector<ttk::SimplexId>(receivedIds, recvMessageSize, mpiIdType_,
738 neighbors_->at(i - hasSentData_));
739 // Outdated global ids are matched to a local point
740 // and its global id is added to the vector send_buf
741 identifyPoints(
742 locatedSimplices, receivedIds, recvMessageSize, send_buf);
743 }
744 // The points founds are sent back to the neighbor with their global
745 // ids
746 sendVector<Response>(
747 send_buf, mpiResponseType_, neighbors_->at(i - hasSentData_));
748 }
749 }
750 // Start of computation of ghost information for cells
751 int id{-1};
752 // If outdated global ids exist, for each ghost cell is added its local id
753 // and its outdated global id
754 if(cellRankArray_ != nullptr) {
755 for(ttk::SimplexId i = 0; i < cellNumber_; i++) {
756 if(cellRankArray_[i] != ttk::MPIrank_) {
757 if(outdatedGlobalCellIds_ == nullptr) {
758 cellGhostGlobalVertexIdsPerRank[neighborToId_[cellRankArray_[i]]]
759 .push_back(i);
760 for(int k = 0; k < dimension_ + 1; k++) {
761 id = connectivity_[i * (dimension_ + 1) + k];
762 cellGhostGlobalVertexIdsPerRank
763 [neighborToId_[cellRankArray_[i]]]
764 .push_back(vertexIdentifiers_[id]);
765 }
766 } else {
767 cellGhostGlobalIdsPerRank[neighborToId_[cellRankArray_[i]]]
768 .push_back(
769 static_cast<ttk::SimplexId>(outdatedGlobalCellIds_[i]));
770 }
771 }
772 }
773 } else {
774 // If no outdated global ids exist, for each ghost cell is added its
775 // local id and the global id of all of its vertices
776 for(ttk::SimplexId i = 0; i < cellNumber_; i++) {
777 if(cellGhost_[i] != 0) {
778 if(outdatedGlobalCellIds_ == nullptr) {
779 cellGhostGlobalVertexIds.push_back(i);
780 for(int k = 0; k < dimension_ + 1; k++) {
781 id = connectivity_[i * (dimension_ + 1) + k];
782 cellGhostGlobalVertexIds.push_back(vertexIdentifiers_[id]);
783 }
784 } else {
785 cellGhostGlobalIds.push_back(
786 static_cast<ttk::SimplexId>(outdatedGlobalCellIds_[i]));
787 }
788 }
789 }
790 }
791 hasSentData_ = 0;
792
793 // Exchange cells similarly to what is done for vertices
794 for(int i = 0; i < neighborNumber_ + 1; i++) {
795 if((i == neighborNumber_ && !hasSentData_)
796 || (!hasSentData_ && ttk::MPIrank_ < neighbors_->at(i))) {
797 // For each neighbor, a process will receive the ghost cells of all
798 // its neighbors or, if it is its turn, will receive the ghost cells
799 // of all its neighbors
800 if(outdatedGlobalCellIds_ == nullptr) {
801 if(cellRankArray_ == nullptr) {
802 sendToAllNeighbors<ttk::SimplexId>(
803 cellGhostGlobalVertexIds, mpiIdType_);
804 } else {
805 sendToAllNeighborsUsingRankArray<ttk::SimplexId>(
806 cellGhostGlobalVertexIdsPerRank, mpiIdType_);
807 }
808 } else {
809 if(cellRankArray_ == nullptr) {
810 sendToAllNeighbors<ttk::SimplexId>(
811 cellGhostGlobalIds, mpiIdType_);
812 } else {
813 sendToAllNeighborsUsingRankArray<ttk::SimplexId>(
814 cellGhostGlobalIdsPerRank, mpiIdType_);
815 }
816 }
817 // The current process receives the responses of all the processes
818 // and associates its ghosts with the right global identifier
819 for(int j = 0; j < neighborNumber_; j++) {
820 recvVector<Response>(receivedResponse, recvMessageSize,
821 mpiResponseType_, neighbors_->at(j));
822 if(outdatedGlobalCellIds_ == nullptr) {
823 for(int n = 0; n < recvMessageSize; n++) {
824 cellIdentifiers_[receivedResponse[n].id]
825 = receivedResponse[n].globalId;
826 }
827 } else {
828 for(int n = 0; n < recvMessageSize; n++) {
829 cellIdentifiers_[cellOutdatedGtoL_[receivedResponse[n].id]]
830 = receivedResponse[n].globalId;
831 }
832 }
833 }
834 hasSentData_ = 1;
835 } else {
836 // It is not the turn of the current process to send its data.
837 recvVector<ttk::SimplexId>(receivedIds, recvMessageSize, mpiIdType_,
838 neighbors_->at(i - hasSentData_));
839 if(outdatedGlobalCellIds_ == nullptr) {
840 // Point global ids are matched to a local cell and the global
841 // id of the cell is added to the vector send_buf
842 locateCells(
843 locatedSimplices, receivedIds, recvMessageSize, send_buf);
844 } else {
845 // Outdated global ids are matched to a local cell
846 // and its global id is added to the vector send_buf
847 identifyCells(
848 locatedSimplices, receivedIds, recvMessageSize, send_buf);
849 }
850 // The cells founds are sent back to the neighbor with their global
851 // ids
852 sendVector<Response>(
853 send_buf, mpiResponseType_, neighbors_->at(i - hasSentData_));
854 }
855 }
856
857 return 1; // return success
858 }
859
865 int executeImageData() {
866
867 // print horizontal separator
868 this->printMsg(ttk::debug::Separator::L1); // L1 is the '=' separator
869
870 // Reorganize bounds to only execute Allreduce twice
871 double tempBounds[6] = {
872 bounds_[0], bounds_[2], bounds_[4], bounds_[1], bounds_[3], bounds_[5]};
873 double tempGlobalBounds[6];
874 // Compute and send to all processes the lower bounds of the data set
875 MPI_Allreduce(
876 tempBounds, tempGlobalBounds, 3, MPI_DOUBLE, MPI_MIN, ttk::MPIcomm_);
877
878 // Compute and send to all processes the higher bounds of the data set
879 MPI_Allreduce(tempBounds + 3, tempGlobalBounds + 3, 3, MPI_DOUBLE,
880 MPI_MAX, ttk::MPIcomm_);
881 // Global bounds
882 double globalBounds[6]
883 = {tempGlobalBounds[0], tempGlobalBounds[3], tempGlobalBounds[1],
884 tempGlobalBounds[4], tempGlobalBounds[2], tempGlobalBounds[5]};
885 // Compute global width and height of the data set
886 int width
887 = static_cast<int>((globalBounds[1] - globalBounds[0]) / spacing_[0])
888 + 1;
889 int height
890 = static_cast<int>((globalBounds[3] - globalBounds[2]) / spacing_[1])
891 + 1;
892
893 // Compute offset of the current process for each direction
894 int offsetWidth = static_cast<int>(
895 std::round((bounds_[0] - globalBounds[0]) / spacing_[0]));
896 int offsetHeight = static_cast<int>(
897 std::round((bounds_[2] - globalBounds[2]) / spacing_[1]));
898 int offsetLength = static_cast<int>(
899 std::round((bounds_[4] - globalBounds[4]) / spacing_[2]));
900 // Generate global ids for vertices
901#ifdef TTK_ENABLE_OPENMP
902#pragma omp parallel for num_threads(threadNumber_)
903#endif
904 for(int k = 0; k < dims_[2]; k++) {
905 for(int j = 0; j < dims_[1]; j++) {
906 for(int i = 0; i < dims_[0]; i++) {
907 vertexIdentifiers_[k * dims_[0] * dims_[1] + j * dims_[0] + i]
908 = i + offsetWidth + (j + offsetHeight) * width
909 + (k + offsetLength) * width * height;
910 }
911 }
912 }
913
914 // Generate global ids for cells
915 dims_[0] -= 1;
916 dims_[1] -= 1;
917 dims_[2] -= 1;
918 width -= 1;
919 height -= 1;
920#ifdef TTK_ENABLE_OPENMP
921#pragma omp parallel for num_threads(threadNumber_)
922#endif
923 for(int k = 0; k < dims_[2]; k++) {
924 for(int j = 0; j < dims_[1]; j++) {
925 for(int i = 0; i < dims_[0]; i++) {
926 cellIdentifiers_[k * dims_[0] * dims_[1] + j * dims_[0] + i]
927 = i + offsetWidth + (j + offsetHeight) * width
928 + (k + offsetLength) * width * height;
929 }
930 }
931 }
932
933 return 1; // return success
934 }
935
936#endif
937
944 // print horizontal separator
945 this->printMsg(ttk::debug::Separator::L1); // L1 is the '=' separator
946
947#ifdef TTK_ENABLE_OPENMP
948#pragma omp parallel for num_threads(threadNumber_)
949#endif
950 for(SimplexId i = 0; i < vertexNumber_; i++) {
951 // avoid any processing if the abort signal is sent
952 vertexIdentifiers_[i] = i;
953 }
954
955#ifdef TTK_ENABLE_OPENMP
956#pragma omp parallel for num_threads(threadNumber_)
957#endif
958 for(SimplexId i = 0; i < cellNumber_; i++) {
959 // avoid any processing if the abort signal is sent
960 cellIdentifiers_[i] = i;
961 }
962
963 return 1; // return success
964 }
965
966 }; // Identifiers class
967
968} // namespace ttk
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition: BaseClass.h:47
int threadNumber_
Definition: BaseClass.h:95
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
void setCellIdentifiers(ttk::LongSimplexId *cellIdentifiers)
Definition: Identifiers.h:97
ttk::SimplexId cellNumber_
Definition: Identifiers.h:54
void setVertexNumber(const SimplexId &vertexNumber)
Definition: Identifiers.h:85
void setCellNumber(const SimplexId &cellNumber)
Definition: Identifiers.h:89
ttk::LongSimplexId * vertexIdentifiers_
Definition: Identifiers.h:55
int executeSequential()
Generates global ids for all data set type in sequential.
Definition: Identifiers.h:943
void setVertexIdentifiers(ttk::LongSimplexId *vertexIdentifiers)
Definition: Identifiers.h:93
~Identifiers() override=default
ttk::SimplexId vertexNumber_
Definition: Identifiers.h:53
ttk::LongSimplexId * cellIdentifiers_
Definition: Identifiers.h:56
TTK KD-Tree.
Definition: KDTree.h:21
void getKClosest(const unsigned int k, const Container &coordinates, KDTreeMap &neighbours, std::vector< dataType > &costs, const int weight_index=0)
Definition: KDTree.h:335
KDTreeMap build(dataType *data, const int &ptNumber, const int &dimension, const std::vector< std::vector< dataType > > &weights={}, const int weight_number=1)
Definition: KDTree.h:144
The Topology ToolKit.
COMMON_EXPORTS int MPIrank_
Definition: BaseClass.cpp:9
long long int LongSimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:15
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22