18#include <unordered_map>
19#include <unordered_set>
25#define OMPI_SKIP_MPICXX 1
31 inline MPI_Datatype getMPIType(
const float ttkNotUsed(val)) {
34 inline MPI_Datatype getMPIType(
const int ttkNotUsed(val)) {
37 inline MPI_Datatype getMPIType(
const unsigned int ttkNotUsed(val)) {
40 inline MPI_Datatype getMPIType(
const double ttkNotUsed(val)) {
43 inline MPI_Datatype getMPIType(
const long double ttkNotUsed(val)) {
44 return MPI_LONG_DOUBLE;
46 inline MPI_Datatype getMPIType(
const long ttkNotUsed(val)) {
49 inline MPI_Datatype getMPIType(
const unsigned long ttkNotUsed(val)) {
50 return MPI_UNSIGNED_LONG;
52 inline MPI_Datatype getMPIType(
const long long ttkNotUsed(val)) {
55 inline MPI_Datatype getMPIType(
const unsigned long long ttkNotUsed(val)) {
56 return MPI_UNSIGNED_LONG_LONG;
58 inline MPI_Datatype getMPIType(
const unsigned char ttkNotUsed(val)) {
59 return MPI_UNSIGNED_CHAR;
62 template <
typename DT,
typename IT>
67 value(DT _scalar, IT _globalId) : scalar(_scalar), globalId(_globalId) {
71 inline bool isRunningWithMPI() {
75 inline bool hasInitializedMPI() {
79 inline int startMPITimer(Timer &t,
int rank,
int size) {
81 MPI_Barrier(ttk::MPIcomm_);
89 inline double endMPITimer(Timer &t,
int rank,
int size) {
90 double elapsedTime = 0;
92 MPI_Barrier(ttk::MPIcomm_);
94 elapsedTime = t.getElapsedTime();
113 template <
typename DT,
typename triangulationType>
114 int getGhostCellScalars(DT *scalarArray,
115 const triangulationType *triangulation,
116 const int rankToSend,
117 MPI_Comm communicator,
118 const int dimensionNumber) {
119 const std::vector<int> &neighbors = triangulation->getNeighborRanks();
120 if(!ttk::isRunningWithMPI()) {
123 MPI_Datatype MPI_DT = getMPIType(
static_cast<DT
>(0));
126 const int tagMultiplier = rankToSend + 1;
127 int valuesTag = 103 * tagMultiplier;
129 int neighborNumber = neighbors.size();
130 const auto &ghostCellsPerOwner = triangulation->getGhostCellsPerOwner();
132 for(
int r = 0; r < neighborNumber; r++) {
134 std::vector<DT> receivedValues(nValues * dimensionNumber);
136 MPI_Recv(receivedValues.data(), nValues * dimensionNumber, MPI_DT,
137 neighbors.at(r), valuesTag, communicator, MPI_STATUS_IGNORE);
140 for(
int j = 0; j < dimensionNumber; j++) {
141 DT receivedVal = receivedValues[i * dimensionNumber + j];
144 scalarArray[localId * dimensionNumber + j] = receivedVal;
152 if(std::find(neighbors.begin(), neighbors.end(), rankToSend)
153 != neighbors.end()) {
155 const auto &ghostCellsForThisRank
156 = triangulation->getRemoteGhostCells()[rankToSend];
160 std::vector<DT> valuesToSend(nValues * dimensionNumber);
162 for(
int j = 0; j < dimensionNumber; j++) {
165 valuesToSend[i * dimensionNumber + j]
166 = scalarArray[localId * dimensionNumber + j];
171 MPI_Send(valuesToSend.data(), nValues * dimensionNumber, MPI_DT,
172 rankToSend, valuesTag, communicator);
180 template <
typename DT,
typename triangulationType>
181 int getGhostVertexScalars(DT *scalarArray,
182 const triangulationType *triangulation,
183 const int rankToSend,
184 MPI_Comm communicator,
185 const int dimensionNumber) {
186 const std::vector<int> &neighbors = triangulation->getNeighborRanks();
187 if(!ttk::isRunningWithMPI()) {
190 MPI_Datatype MPI_DT = getMPIType(
static_cast<DT
>(0));
193 const int tagMultiplier = rankToSend + 1;
194 int valuesTag = 103 * tagMultiplier;
196 int neighborNumber = neighbors.size();
197 const auto &ghostVerticesPerOwner
198 = triangulation->getGhostVerticesPerOwner();
200 for(
int r = 0; r < neighborNumber; r++) {
201 ttk::SimplexId nValues = ghostVerticesPerOwner[neighbors[r]].size();
202 std::vector<DT> receivedValues(nValues * dimensionNumber);
204 MPI_Recv(receivedValues.data(), nValues * dimensionNumber, MPI_DT,
205 neighbors.at(r), valuesTag, communicator, MPI_STATUS_IGNORE);
206#pragma omp parallel for
208 for(
int j = 0; j < dimensionNumber; j++) {
209 DT receivedVal = receivedValues[i * dimensionNumber + j];
212 = triangulation->getVertexLocalId(globalId);
213 scalarArray[localId * dimensionNumber + j] = receivedVal;
221 if(std::find(neighbors.begin(), neighbors.end(), rankToSend)
222 != neighbors.end()) {
224 const auto &ghostVerticesForThisRank
225 = triangulation->getRemoteGhostVertices()[rankToSend];
229 std::vector<DT> valuesToSend(nValues * dimensionNumber);
230#pragma omp parallel for
232 for(
int j = 0; j < dimensionNumber; j++) {
235 = triangulation->getVertexLocalId(globalId);
236 valuesToSend[i * dimensionNumber + j]
237 = scalarArray[localId * dimensionNumber + j];
242 MPI_Send(valuesToSend.data(), nValues * dimensionNumber, MPI_DT,
243 rankToSend, valuesTag, communicator);
272 template <
typename DT,
277 int getGhostDataScalarsWithoutTriangulation(DT *scalarArray,
278 const GVR &getVertexRank,
279 const GVGID &getVertexGlobalId,
280 const GVLID &getVertexLocalId,
281 const std::vector<int> &neighbors,
282 const int rankToSend,
284 MPI_Comm communicator,
285 const int dimensionNumber) {
286 const int neighborNumber = neighbors.size();
287 if(!ttk::isRunningWithMPI()) {
290 MPI_Datatype MPI_DT = getMPIType(
static_cast<DT
>(0));
291 MPI_Datatype MPI_IT = getMPIType(
static_cast<IT
>(0));
292 using globalIdType =
decltype(getVertexGlobalId(0));
293 MPI_Datatype MPI_GIT = getMPIType(
static_cast<globalIdType
>(0));
296 const int tagMultiplier = rankToSend + 1;
297 int amountTag = 101 * tagMultiplier;
298 int idsTag = 102 * tagMultiplier;
299 int valuesTag = 103 * tagMultiplier;
302 std::vector<std::vector<globalIdType>> rankVectors(
306 for(IT i = 0; i < nVerts; i++) {
308 rankVectors[getVertexRank(i)].push_back(getVertexGlobalId(i));
312 for(
int r = 0; r < neighborNumber; r++) {
313 IT nValues = rankVectors[neighbors[r]].size();
314 MPI_Send(&nValues, 1, MPI_IT, neighbors[r], amountTag, communicator);
316 MPI_Send(rankVectors[neighbors[r]].data(), nValues, MPI_GIT,
317 neighbors[r], idsTag, communicator);
321 for(
int r = 0; r < neighborNumber; r++) {
322 IT nValues = rankVectors[neighbors[r]].size();
323 std::vector<DT> receivedValues(nValues * dimensionNumber);
325 MPI_Recv(receivedValues.data(), nValues * dimensionNumber, MPI_DT,
326 neighbors[r], valuesTag, communicator, MPI_STATUS_IGNORE);
327 for(IT i = 0; i < nValues; i++) {
328 for(
int j = 0; j < dimensionNumber; j++) {
329 DT receivedVal = receivedValues[i * dimensionNumber + j];
330 const auto globalId = rankVectors[neighbors[r]][i];
331 IT localId = getVertexLocalId(globalId);
332 scalarArray[localId * dimensionNumber + j] = receivedVal;
340 if(std::find(neighbors.begin(), neighbors.end(), rankToSend)
341 != neighbors.end()) {
345 MPI_Recv(&nValues, 1, MPI_IT, rankToSend, amountTag, communicator,
349 std::vector<globalIdType> receivedIds(nValues);
350 MPI_Recv(receivedIds.data(), nValues, MPI_GIT, rankToSend, idsTag,
351 communicator, MPI_STATUS_IGNORE);
354 std::vector<DT> valuesToSend(nValues * dimensionNumber);
355 for(IT i = 0; i < nValues; i++) {
356 for(
int j = 0; j < dimensionNumber; j++) {
357 IT globalId = receivedIds[i];
358 IT localId = getVertexLocalId(globalId);
359 valuesToSend[i * dimensionNumber + j]
360 = scalarArray[localId * dimensionNumber + j];
365 MPI_Send(valuesToSend.data(), nValues * dimensionNumber, MPI_DT,
366 rankToSend, valuesTag, communicator);
383 template <
typename IT,
typename GVR>
384 int preconditionNeighborsUsingRankArray(std::vector<int> &neighbors,
385 std::map<int, int> &neighborsToId,
386 const GVR &getVertexRank,
388 MPI_Comm communicator) {
389 std::unordered_set<int> neighborSet{};
390 for(IT i = 0; i < nVerts; i++) {
392 neighborSet.emplace(getVertexRank(i));
395 std::vector<int> sendVector(neighborSet.begin(), neighborSet.end());
396 int localSize = neighborSet.size();
400 &localSize, 1, MPI_INT, sizes.data(), 1, MPI_INT, 0, communicator);
404 totalSize += sizes[i];
406 displacements[i] = 0;
408 displacements[i] = displacements[i - 1] + sizes[i - 1];
412 std::vector<int> rootVector(totalSize);
414 MPI_Gatherv(sendVector.data(), sendVector.size(), MPI_INT,
415 rootVector.data(), sizes.data(), displacements.data(), MPI_INT,
417 std::vector<int> scatterVector;
422 std::vector<std::unordered_set<int>> setsFromRanks(
ttk::MPIsize_);
423 auto begin = rootVector.begin();
424 auto end = rootVector.begin();
427 const std::unordered_set<int> s(
begin,
end);
428 setsFromRanks[i] = s;
435 if(setsFromRanks[j].find(i) != setsFromRanks[j].
end()) {
436 setsFromRanks[i].emplace(j);
443 sizes[i] = setsFromRanks[i].size();
445 displacements[i] = 0;
447 displacements[i] = displacements[i - 1] + sizes[i - 1];
449 scatterVector.insert(scatterVector.end(), setsFromRanks[i].begin(),
450 setsFromRanks[i].end());
457 sizes.data(), 1, MPI_INT, &receivedSize, 1, MPI_INT, 0, communicator);
460 std::vector<int> receivedNeighbors(receivedSize);
461 MPI_Scatterv(scatterVector.data(), sizes.data(), displacements.data(),
462 MPI_INT, receivedNeighbors.data(), receivedSize, MPI_INT, 0,
465 const std::unordered_set<int> finalSet(
466 receivedNeighbors.begin(), receivedNeighbors.end());
467 neighborSet = finalSet;
471 for(
const int neighbor : neighborSet) {
472 neighbors.push_back(neighbor);
474 std::sort(neighbors.begin(), neighbors.end());
475 int neighborNumber = neighbors.size();
476 neighborsToId.clear();
477 for(
int i = 0; i < neighborNumber; i++) {
478 neighborsToId[neighbors[i]] = i;
495 template <
typename DT,
typename triangulationType>
496 int exchangeGhostCells(DT *scalarArray,
497 const triangulationType *triangulation,
498 MPI_Comm communicator,
499 const int dimensionNumber = 1) {
500 if(!ttk::isRunningWithMPI()) {
503 if(!triangulation->hasPreconditionedDistributedCells()) {
506 const std::vector<int> &neighbors = triangulation->getNeighborRanks();
507 const int neighborNumber = neighbors.size();
508 const std::map<int, int> &neighborsToId = triangulation->getNeighborsToId();
509 if(!ttk::isRunningWithMPI()) {
513 std::vector<std::vector<DT>> receivedValues(
514 neighborNumber, std::vector<DT>());
515 std::vector<std::vector<DT>> valuesToSend(
516 neighborNumber, std::vector<DT>());
517 MPI_Datatype MPI_DT = getMPIType(
static_cast<DT
>(0));
518 const auto &ghostCellsPerOwner = triangulation->getGhostCellsPerOwner();
519 const auto &remoteGhostCells = triangulation->getRemoteGhostCells();
520 for(
int i = 0; i < neighborNumber; i++) {
521 receivedValues[i].resize(ghostCellsPerOwner[neighbors[i]].size()
523 valuesToSend[i].resize(remoteGhostCells[neighbors[i]].size()
526 for(
int i = 0; i < neighborNumber; i++) {
528#pragma omp parallel for
530 for(
int k = 0; k < dimensionNumber; k++) {
533 valuesToSend[i].at(j * dimensionNumber + k)
534 = scalarArray[localId * dimensionNumber + k];
539 std::vector<MPI_Request> sendRequests(neighborNumber);
540 std::vector<MPI_Request> recvRequests(neighborNumber);
541 for(
int i = 0; i < neighborNumber; i++) {
542 MPI_Isend(valuesToSend[i].data(), valuesToSend[i].size(), MPI_DT,
543 neighbors[i], 0, communicator, &sendRequests[i]);
544 MPI_Irecv(receivedValues[i].data(), receivedValues[i].size(), MPI_DT,
545 neighbors[i], 0, communicator, &recvRequests[i]);
548 std::vector<MPI_Status> recvStatus(neighborNumber);
549 std::vector<int> recvCompleted(neighborNumber, 0);
550 int recvPerformedCount = 0;
551 int recvPerformedCountTotal = 0;
554 while(recvPerformedCountTotal < neighborNumber) {
555 MPI_Waitsome(neighborNumber, recvRequests.data(), &recvPerformedCount,
556 recvCompleted.data(), recvStatus.data());
557 if(recvPerformedCount > 0) {
558 for(
int i = 0; i < recvPerformedCount; i++) {
559 r = recvStatus[i].MPI_SOURCE;
560 rId = neighborsToId.at(r);
562#pragma omp parallel for
564 for(
int k = 0; k < dimensionNumber; k++) {
565 DT receivedVal = receivedValues[rId][j * dimensionNumber + k];
568 scalarArray[localId * dimensionNumber + k] = receivedVal;
572 recvPerformedCountTotal += recvPerformedCount;
575 MPI_Waitall(neighborNumber, sendRequests.data(), MPI_STATUSES_IGNORE);
580 template <
typename DT,
typename triangulationType>
581 int exchangeGhostVertices(DT *scalarArray,
582 const triangulationType *triangulation,
583 MPI_Comm communicator,
584 const int dimensionNumber = 1) {
585 if(!ttk::isRunningWithMPI()) {
588 if(!triangulation->hasPreconditionedDistributedVertices()) {
591 const std::vector<int> &neighbors = triangulation->getNeighborRanks();
592 const int neighborNumber = neighbors.size();
593 const std::map<int, int> &neighborsToId = triangulation->getNeighborsToId();
595 if(!ttk::isRunningWithMPI()) {
599 std::vector<std::vector<DT>> receivedValues(
600 neighborNumber, std::vector<DT>());
601 std::vector<std::vector<DT>> valuesToSend(
602 neighborNumber, std::vector<DT>());
603 MPI_Datatype MPI_DT = getMPIType(
static_cast<DT
>(0));
604 const auto &ghostVerticesPerOwner
605 = triangulation->getGhostVerticesPerOwner();
606 const auto &remoteGhostVertices = triangulation->getRemoteGhostVertices();
607 for(
int i = 0; i < neighborNumber; i++) {
608 receivedValues[i].resize(ghostVerticesPerOwner[neighbors[i]].size()
610 valuesToSend[i].resize(remoteGhostVertices[neighbors[i]].size()
613 for(
int i = 0; i < neighborNumber; i++) {
614 ttk::SimplexId nValues = remoteGhostVertices[neighbors[i]].size();
615#pragma omp parallel for
617 for(
int k = 0; k < dimensionNumber; k++) {
619 ttk::SimplexId localId = triangulation->getVertexLocalId(globalId);
620 valuesToSend[i].at(j * dimensionNumber + k)
621 = scalarArray[localId * dimensionNumber + k];
626 std::vector<MPI_Request> sendRequests(neighborNumber);
627 std::vector<MPI_Request> recvRequests(neighborNumber);
628 for(
int i = 0; i < neighborNumber; i++) {
629 MPI_Isend(valuesToSend[i].data(), valuesToSend[i].size(), MPI_DT,
630 neighbors[i], 0, communicator, &sendRequests[i]);
631 MPI_Irecv(receivedValues[i].data(), receivedValues[i].size(), MPI_DT,
632 neighbors[i], 0, communicator, &recvRequests[i]);
635 std::vector<MPI_Status> recvStatus(neighborNumber);
636 std::vector<int> recvCompleted(neighborNumber, 0);
637 int recvPerformedCount = 0;
638 int recvPerformedCountTotal = 0;
641 while(recvPerformedCountTotal < neighborNumber) {
642 MPI_Waitsome(neighborNumber, recvRequests.data(), &recvPerformedCount,
643 recvCompleted.data(), recvStatus.data());
644 if(recvPerformedCount > 0) {
645 for(
int i = 0; i < recvPerformedCount; i++) {
646 r = recvStatus[i].MPI_SOURCE;
647 rId = neighborsToId.at(r);
649#pragma omp parallel for
651 for(
int k = 0; k < dimensionNumber; k++) {
652 DT receivedVal = receivedValues[rId][j * dimensionNumber + k];
655 = triangulation->getVertexLocalId(globalId);
656 scalarArray[localId * dimensionNumber + k] = receivedVal;
660 recvPerformedCountTotal += recvPerformedCount;
663 MPI_Waitall(neighborNumber, sendRequests.data(), MPI_STATUSES_IGNORE);
689 template <
typename DT,
694 int exchangeGhostDataWithoutTriangulation(DT *scalarArray,
695 const GVR &getVertexRank,
696 const GVGID &getVertexGlobalId,
697 const GVLID &getVertexLocalId,
699 MPI_Comm communicator,
700 const std::vector<int> &neighbors,
701 const int dimensionNumber = 1) {
702 if(!ttk::isRunningWithMPI()) {
706 getGhostDataScalarsWithoutTriangulation(
707 scalarArray, getVertexRank, getVertexGlobalId, getVertexLocalId,
708 neighbors, r, nVerts, communicator, dimensionNumber);
709 MPI_Barrier(communicator);
715 bool inline checkForIntersection(
double *myBB,
double *theirBB) {
718 || myBB[1] < theirBB[0]
719 || myBB[2] > theirBB[3]
720 || myBB[3] < theirBB[2]
721 || myBB[4] > theirBB[5]
722 || myBB[5] < theirBB[4]
726 void inline preconditionNeighborsUsingBoundingBox(
728 std::vector<int> &neighbors,
729 std::map<int, int> &neighborsToId) {
731 std::vector<std::array<double, 6>> rankBoundingBoxes(
ttk::MPIsize_);
735 MPI_Bcast(rankBoundingBoxes[r].data(), 6, MPI_DOUBLE, r, ttk::MPIcomm_);
738 const double epsilon = 0.00001;
740 for(
int i = 0; i < 6; i++) {
742 boundingBox[i] -= epsilon;
744 boundingBox[i] += epsilon;
749 double *theirBoundingBox = rankBoundingBoxes[i].data();
750 if(checkForIntersection(boundingBox, theirBoundingBox)) {
751 neighbors.push_back(i);
755 neighborsToId.clear();
756 int neighborNumber = neighbors.size();
757 for(
int i = 0; i < neighborNumber; i++) {
758 neighborsToId[neighbors[i]] = i;
771 void inline produceRankArray(std::vector<int> &rankArray,
773 const unsigned char *ghostCells,
776 std::vector<int> &neighbors,
777 std::map<int, int> &neighborsToId) {
778 if(neighbors.empty()) {
779 ttk::preconditionNeighborsUsingBoundingBox(
780 boundingBox, neighbors, neighborsToId);
783 std::vector<ttk::SimplexId> currentRankUnknownIds;
784 std::unordered_set<ttk::SimplexId> gIdSet;
785 std::unordered_map<ttk::SimplexId, ttk::SimplexId> gIdToLocalMap;
787 for(
int i = 0; i < nVertices; i++) {
788 const int ghostCellVal = ghostCells[i];
790 if(ghostCellVal == 0) {
794 gIdSet.insert(globalId);
799 currentRankUnknownIds.push_back(globalId);
800 gIdToLocalMap[globalId] = i;
805 std::vector<ttk::SimplexId> gIdsToSend;
806 std::vector<ttk::SimplexId> receivedGlobals;
807 receivedGlobals.resize(sizeOfCurrentRank);
809 std::vector<ttk::SimplexId> neighborUnknownIds;
810 for(
const int neighbor : neighbors) {
812 MPI_Sendrecv(&sizeOfCurrentRank, 1, MIT, neighbor,
ttk::MPIrank_,
813 &sizeOfNeighbor, 1, MIT, neighbor, neighbor, ttk::MPIcomm_,
815 neighborUnknownIds.resize(sizeOfNeighbor);
816 gIdsToSend.reserve(sizeOfNeighbor);
818 MPI_Sendrecv(currentRankUnknownIds.data(), sizeOfCurrentRank, MIT,
820 sizeOfNeighbor, MIT, neighbor, neighbor, ttk::MPIcomm_,
826 if(gIdSet.count(gId)) {
828 gIdsToSend.push_back(gId);
834 MPI_Sendrecv(gIdsToSend.data(), gIdsToSend.size(), MIT, neighbor,
836 currentRankUnknownIds.size(), MIT, neighbor, neighbor,
837 ttk::MPIcomm_, &status);
839 MPI_Get_count(&status, MIT, &amount);
840 receivedGlobals.resize(amount);
844 rankArray[localVal] = neighbor;
848 receivedGlobals.resize(sizeOfCurrentRank);
860 template <
typename DT,
typename IT>
861 void returnVectorForBurstsize(std::vector<value<DT, IT>> &outVector,
862 std::vector<value<DT, IT>> &values,
865 if(burstSize > values.size()) {
866 outVector.resize(values.size(), {0, 0});
867 outVector.assign(values.begin(), values.end());
870 outVector.resize(burstSize, {0, 0});
871 outVector.assign(values.end() - burstSize, values.end());
872 values.erase(values.end() - burstSize, values.end());
885 template <
typename DT,
typename IT>
886 void ReceiveAndAddToVector(
889 std::vector<std::vector<value<DT, IT>>> &unsortedReceivedValues) {
890 std::vector<value<DT, IT>> &receivedValues
891 = unsortedReceivedValues[rankFrom];
895 MPI_Probe(rankFrom, structTag * rankFrom, ttk::MPIcomm_, &status);
896 MPI_Get_count(&status, MPI_CHAR, &amount);
897 receivedValues.resize(amount /
sizeof(value<DT, IT>), {0, 0});
898 MPI_Recv(receivedValues.data(), amount, MPI_CHAR, rankFrom,
899 structTag * rankFrom, ttk::MPIcomm_, MPI_STATUS_IGNORE);
926 template <
typename DT,
typename IT>
927 void getMax(
int intTag,
932 IT &processedValueCounter,
933 std::vector<std::vector<value<DT, IT>>> &unsortedReceivedValues,
934 std::vector<std::vector<IT>> &orderResendValues,
935 std::vector<IT> &orderedValuesForRank,
936 std::vector<value<DT, IT>> &sortingValues) {
938 int rankIdOfMaxScalar = -1;
939 DT maxScalar = std::numeric_limits<DT>::lowest();
941 for(
size_t i = 0; i < unsortedReceivedValues.size(); i++) {
942 if(unsortedReceivedValues[i].size() > 0) {
943 const auto &v = unsortedReceivedValues[i].back();
944 if(v.scalar == maxScalar ? v.globalId > maxGId : v.scalar > maxScalar) {
945 maxScalar = v.scalar;
947 rankIdOfMaxScalar = i;
952 value<DT, IT> currentValue
953 = unsortedReceivedValues[rankIdOfMaxScalar].back();
957 orderResendValues[rankIdOfMaxScalar].push_back(currentValue.globalId);
958 orderResendValues[rankIdOfMaxScalar].push_back(currentOrder);
960 processedValueCounter++;
961 unsortedReceivedValues[rankIdOfMaxScalar].pop_back();
962 if(unsortedReceivedValues[rankIdOfMaxScalar].size() == 0) {
963 if(rankIdOfMaxScalar == 0) {
965 orderedValuesForRank.insert(
966 orderedValuesForRank.end(),
967 orderResendValues[rankIdOfMaxScalar].begin(),
968 orderResendValues[rankIdOfMaxScalar].end());
969 orderResendValues[rankIdOfMaxScalar].clear();
971 if(!sortingValues.empty()) {
972 std::vector<value<DT, IT>> ownValues;
973 returnVectorForBurstsize<DT, IT>(ownValues, sortingValues, burstSize);
975 unsortedReceivedValues[rankIdOfMaxScalar] = ownValues;
980 MPI_Send(orderResendValues[rankIdOfMaxScalar].data(),
981 orderResendValues[rankIdOfMaxScalar].size(), MPI_IT,
982 rankIdOfMaxScalar, intTag * rankIdOfMaxScalar, ttk::MPIcomm_);
983 orderResendValues[rankIdOfMaxScalar].clear();
986 ReceiveAndAddToVector(
987 rankIdOfMaxScalar, structTag, unsortedReceivedValues);
1003 template <
typename IT,
typename GVLID>
1004 void buildArrayForReceivedData(
const size_t nInts,
1005 const IT *
const orderedValuesForRank,
1006 const GVLID &getVertexLocalId,
1008 const int nThreads) {
1012#ifdef TTK_ENABLE_OPENMP
1013#pragma omp parallel for num_threads(nThreads)
1015 for(
size_t i = 0; i < nInts; i += 2) {
1016 order[getVertexLocalId(orderedValuesForRank[i])]
1017 = orderedValuesForRank[i + 1];
1034 template <
typename DT,
typename IT,
typename GVGID,
typename GVR>
1035 void populateVector(std::vector<value<DT, IT>> &valuesToSortVector,
1036 const size_t nVerts,
1037 const DT *
const scalars,
1038 const GVGID &getVertexGlobalId,
1039 const GVR &getVertexRank) {
1040 for(
size_t i = 0; i < nVerts; i++) {
1041 IT globalId = getVertexGlobalId(i);
1043 valuesToSortVector.emplace_back(scalars[i], globalId);
1055 template <
typename DT,
typename IT>
1056 void sortVerticesDistributed(std::vector<value<DT, IT>> &values,
1057 const int nThreads) {
1061 TTK_PSORT(nThreads, values.begin(), values.end(),
1062 [](value<DT, IT> v1, value<DT, IT> v2) {
1063 return (v1.scalar < v2.scalar)
1064 || (v1.scalar == v2.scalar && v1.globalId < v2.globalId);
1079 template <
typename DT,
typename GVGID,
typename GVR,
typename GVLID>
1080 void produceOrdering(
SimplexId *orderArray,
1081 const DT *scalarArray,
1082 const GVGID &getVertexGlobalId,
1083 const GVR &getVertexRank,
1084 const GVLID &getVertexLocalId,
1085 const size_t nVerts,
1086 const int burstSize,
1087 std::vector<int> &neighbors,
1088 std::map<int, int> &neighborsToId) {
1090 int structTag = 102;
1091 if(neighbors.empty()) {
1092 ttk::preconditionNeighborsUsingRankArray(
1093 neighbors, neighborsToId, getVertexRank, nVerts, ttk::MPIcomm_);
1095 MPI_Barrier(ttk::MPIcomm_);
1097 using IT =
decltype(getVertexGlobalId(0));
1098 MPI_Datatype MPI_IT = ttk::getMPIType(
static_cast<IT
>(0));
1100 std::vector<value<DT, IT>> sortingValues;
1102 sortingValues, nVerts, scalarArray, getVertexGlobalId, getVertexRank);
1110 MPI_Barrier(ttk::MPIcomm_);
1112 std::vector<IT> orderedValuesForRank;
1113 IT processedValueCounter = 0;
1114 IT localSize = sortingValues.size();
1117 MPI_Reduce(&localSize, &totalSize, 1, MPI_IT, MPI_SUM, 0, ttk::MPIcomm_);
1119 IT currentOrder = totalSize - 1;
1120 std::vector<std::vector<value<DT, IT>>> unsortedReceivedValues;
1122 std::vector<std::vector<IT>> orderResendValues;
1125 orderResendValues[i].reserve(burstSize);
1131 std::vector<value<DT, IT>> ownValues;
1132 returnVectorForBurstsize<DT, IT>(ownValues, sortingValues, burstSize);
1133 unsortedReceivedValues[i] = ownValues;
1135 ReceiveAndAddToVector<DT, IT>(i, structTag, unsortedReceivedValues);
1138 while(processedValueCounter < totalSize) {
1139 getMax<DT, IT>(intTag, structTag, currentOrder, burstSize, MPI_IT,
1140 processedValueCounter, unsortedReceivedValues,
1141 orderResendValues, orderedValuesForRank, sortingValues);
1147 while(!sortingValues.empty()) {
1148 std::vector<value<DT, IT>> sendValues;
1149 returnVectorForBurstsize<DT, IT>(sendValues, sortingValues, burstSize);
1150 int size = sendValues.size();
1151 MPI_Send(sendValues.data(), size *
sizeof(value<DT, IT>), MPI_CHAR, 0,
1153 std::vector<IT> receivedValues;
1157 receivedValues.resize(burstSize * 2);
1160 MPI_Recv(receivedValues.data(), burstSize * 2, MPI_IT, 0,
1162 MPI_Get_count(&status, MPI_IT, &amount);
1164 receivedValues.resize(amount);
1165 orderedValuesForRank.insert(orderedValuesForRank.end(),
1166 receivedValues.begin(),
1167 receivedValues.end());
1171 MPI_Send(sortingValues.data(), 0, MPI_CHAR, 0, structTag *
ttk::MPIrank_,
1176 MPI_Barrier(ttk::MPIcomm_);
1178 buildArrayForReceivedData<IT>(orderedValuesForRank.size(),
1179 orderedValuesForRank.data(), getVertexLocalId,
1184 ttk::exchangeGhostDataWithoutTriangulation<ttk::SimplexId, IT>(
1185 orderArray, getVertexRank, getVertexGlobalId, getVertexLocalId, nVerts,
1186 ttk::MPIcomm_, neighbors);
1198 template <
typename dataType>
1199 void sendVector(std::vector<dataType> &sendBuffer,
1200 MPI_Datatype messageType,
1203 MPI_Send(&dataSize, 1, getMPIType(dataSize), neighbor,
ttk::MPIrank_,
1205 MPI_Send(sendBuffer.data(), dataSize, messageType, neighbor,
ttk::MPIrank_,
1218 template <
typename dataType>
1219 void recvVector(std::vector<dataType> &receiveBuffer,
1221 MPI_Datatype messageType,
1223 MPI_Recv(&recvMessageSize, 1, getMPIType(recvMessageSize), neighbor,
1224 neighbor, ttk::MPIcomm_, MPI_STATUS_IGNORE);
1225 receiveBuffer.resize(recvMessageSize);
1226 MPI_Recv(receiveBuffer.data(), recvMessageSize, messageType, neighbor,
1227 neighbor, ttk::MPIcomm_, MPI_STATUS_IGNORE);
#define TTK_FORCE_USE(x)
Force the compiler to use the function/method parameter.
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
#define TTK_PSORT(NTHREADS,...)
Parallel sort macro.
COMMON_EXPORTS int MPIsize_
COMMON_EXPORTS int globalThreadNumber_
COMMON_EXPORTS int MPIrank_
long long int LongSimplexId
Identifier type for simplices of any dimension.
int SimplexId
Identifier type for simplices of any dimension.
T end(std::pair< T, T > &p)
T begin(std::pair< T, T > &p)