1f84a5eb8SVaclav Hapla #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2f84a5eb8SVaclav Hapla 3f84a5eb8SVaclav Hapla /* TODO PetscArrayExchangeBegin/End */ 4f84a5eb8SVaclav Hapla /* TODO blocksize */ 5f84a5eb8SVaclav Hapla /* TODO move to API ? */ 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode ExchangeArrayByRank_Private(PetscObject obj, MPI_Datatype dt, PetscInt nsranks, const PetscMPIInt sranks[], PetscInt ssize[], const void *sarr[], PetscInt nrranks, const PetscMPIInt rranks[], PetscInt *rsize_out[], void **rarr_out[]) 7d71ae5a4SJacob Faibussowitsch { 8f84a5eb8SVaclav Hapla PetscInt r; 9f84a5eb8SVaclav Hapla PetscInt *rsize; 10f84a5eb8SVaclav Hapla void **rarr; 11f84a5eb8SVaclav Hapla MPI_Request *sreq, *rreq; 12f84a5eb8SVaclav Hapla PetscMPIInt tag, unitsize; 13f84a5eb8SVaclav Hapla MPI_Comm comm; 14f84a5eb8SVaclav Hapla 15f84a5eb8SVaclav Hapla PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_size(dt, &unitsize)); 179566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm(obj, &comm)); 189566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrranks, &rsize, nrranks, &rarr)); 199566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrranks, &rreq, nsranks, &sreq)); 20f84a5eb8SVaclav Hapla /* exchange array size */ 219566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag(obj, &tag)); 2248a46eb9SPierre Jolivet for (r = 0; r < nrranks; r++) PetscCallMPI(MPI_Irecv(&rsize[r], 1, MPIU_INT, rranks[r], tag, comm, &rreq[r])); 2348a46eb9SPierre Jolivet for (r = 0; r < nsranks; r++) PetscCallMPI(MPI_Isend(&ssize[r], 1, MPIU_INT, sranks[r], tag, comm, &sreq[r])); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE)); 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE)); 26f84a5eb8SVaclav Hapla /* exchange array */ 279566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag(obj, &tag)); 28f84a5eb8SVaclav Hapla for (r = 0; r < nrranks; r++) { 299566063dSJacob Faibussowitsch PetscCall(PetscMalloc(rsize[r] * unitsize, &rarr[r])); 309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(rarr[r], rsize[r], dt, rranks[r], tag, comm, &rreq[r])); 31f84a5eb8SVaclav Hapla } 3248a46eb9SPierre Jolivet for (r = 0; r < nsranks; r++) PetscCallMPI(MPI_Isend(sarr[r], ssize[r], dt, sranks[r], tag, comm, &sreq[r])); 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE)); 349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE)); 359566063dSJacob Faibussowitsch PetscCall(PetscFree2(rreq, sreq)); 36f84a5eb8SVaclav Hapla *rsize_out = rsize; 37f84a5eb8SVaclav Hapla *rarr_out = rarr; 38f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 39f84a5eb8SVaclav Hapla } 40f84a5eb8SVaclav Hapla 41f84a5eb8SVaclav Hapla /* TODO VecExchangeBegin/End */ 42f84a5eb8SVaclav Hapla /* TODO move to API ? */ 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode ExchangeVecByRank_Private(PetscObject obj, PetscInt nsranks, const PetscMPIInt sranks[], Vec svecs[], PetscInt nrranks, const PetscMPIInt rranks[], Vec *rvecs[]) 44d71ae5a4SJacob Faibussowitsch { 45f84a5eb8SVaclav Hapla PetscInt r; 46f84a5eb8SVaclav Hapla PetscInt *ssize, *rsize; 47f84a5eb8SVaclav Hapla PetscScalar **rarr; 48f84a5eb8SVaclav Hapla const PetscScalar **sarr; 49f84a5eb8SVaclav Hapla Vec *rvecs_; 50f84a5eb8SVaclav Hapla MPI_Request *sreq, *rreq; 51f84a5eb8SVaclav Hapla 52f84a5eb8SVaclav Hapla PetscFunctionBegin; 539566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(nsranks, &ssize, nsranks, &sarr, nrranks, &rreq, nsranks, &sreq)); 54f84a5eb8SVaclav Hapla for (r = 0; r < nsranks; r++) { 559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(svecs[r], &ssize[r])); 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(svecs[r], &sarr[r])); 57f84a5eb8SVaclav Hapla } 589566063dSJacob Faibussowitsch PetscCall(ExchangeArrayByRank_Private(obj, MPIU_SCALAR, nsranks, sranks, ssize, (const void **)sarr, nrranks, rranks, &rsize, (void ***)&rarr)); 599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrranks, &rvecs_)); 60f84a5eb8SVaclav Hapla for (r = 0; r < nrranks; r++) { 61f84a5eb8SVaclav Hapla /* set array in two steps to mimic PETSC_OWN_POINTER */ 629566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, rsize[r], NULL, &rvecs_[r])); 639566063dSJacob Faibussowitsch PetscCall(VecReplaceArray(rvecs_[r], rarr[r])); 64f84a5eb8SVaclav Hapla } 6548a46eb9SPierre Jolivet for (r = 0; r < nsranks; r++) PetscCall(VecRestoreArrayRead(svecs[r], &sarr[r])); 669566063dSJacob Faibussowitsch PetscCall(PetscFree2(rsize, rarr)); 679566063dSJacob Faibussowitsch PetscCall(PetscFree4(ssize, sarr, rreq, sreq)); 68f84a5eb8SVaclav Hapla *rvecs = rvecs_; 69f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 70f84a5eb8SVaclav Hapla } 71f84a5eb8SVaclav Hapla 72d71ae5a4SJacob Faibussowitsch static PetscErrorCode SortByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[]) 73d71ae5a4SJacob Faibussowitsch { 74f84a5eb8SVaclav Hapla PetscInt nleaves; 75f84a5eb8SVaclav Hapla PetscInt nranks; 76f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 77f84a5eb8SVaclav Hapla const PetscInt *roffset, *rmine, *rremote; 78f84a5eb8SVaclav Hapla PetscInt n, o, r; 79f84a5eb8SVaclav Hapla 80f84a5eb8SVaclav Hapla PetscFunctionBegin; 819566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote)); 82f84a5eb8SVaclav Hapla nleaves = roffset[nranks]; 839566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nleaves, rmine1, nleaves, rremote1)); 84f84a5eb8SVaclav Hapla for (r = 0; r < nranks; r++) { 85f84a5eb8SVaclav Hapla /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote 86f84a5eb8SVaclav Hapla - to unify order with the other side */ 87f84a5eb8SVaclav Hapla o = roffset[r]; 88f84a5eb8SVaclav Hapla n = roffset[r + 1] - o; 899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&(*rmine1)[o], &rmine[o], n)); 909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&(*rremote1)[o], &rremote[o], n)); 919566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o])); 92f84a5eb8SVaclav Hapla } 93f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 94f84a5eb8SVaclav Hapla } 95f84a5eb8SVaclav Hapla 96d71ae5a4SJacob Faibussowitsch static PetscErrorCode GetRecursiveConeCoordinatesPerRank_Private(DM dm, PetscSF sf, PetscInt rmine[], Vec *coordinatesPerRank[]) 97d71ae5a4SJacob Faibussowitsch { 98f84a5eb8SVaclav Hapla IS pointsPerRank, conesPerRank; 99f84a5eb8SVaclav Hapla PetscInt nranks; 100f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 101f84a5eb8SVaclav Hapla const PetscInt *roffset; 102f84a5eb8SVaclav Hapla PetscInt n, o, r; 103f84a5eb8SVaclav Hapla 104f84a5eb8SVaclav Hapla PetscFunctionBegin; 1059566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalSetUp(dm)); 1069566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, NULL, NULL)); 1079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nranks, coordinatesPerRank)); 108f84a5eb8SVaclav Hapla for (r = 0; r < nranks; r++) { 109f84a5eb8SVaclav Hapla o = roffset[r]; 110f84a5eb8SVaclav Hapla n = roffset[r + 1] - o; 1119566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, &rmine[o], PETSC_USE_POINTER, &pointsPerRank)); 1129566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeRecursiveVertices(dm, pointsPerRank, &conesPerRank)); 1139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalTuple(dm, conesPerRank, NULL, &(*coordinatesPerRank)[r])); 1149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointsPerRank)); 1159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&conesPerRank)); 116f84a5eb8SVaclav Hapla } 117f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 118f84a5eb8SVaclav Hapla } 119f84a5eb8SVaclav Hapla 120d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFComputeMultiRootOriginalNumberingByRank_Private(PetscSF sf, PetscSF imsf, PetscInt *irmine1[]) 121d71ae5a4SJacob Faibussowitsch { 122f84a5eb8SVaclav Hapla PetscInt *mRootsOrigNumbering; 123f84a5eb8SVaclav Hapla PetscInt nileaves, niranks; 124f84a5eb8SVaclav Hapla const PetscInt *iroffset, *irmine, *degree; 125f84a5eb8SVaclav Hapla PetscInt i, n, o, r; 126f84a5eb8SVaclav Hapla 127f84a5eb8SVaclav Hapla PetscFunctionBegin; 1289566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL)); 1299566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks(imsf, &niranks, NULL, &iroffset, &irmine, NULL)); 13008401ef6SPierre Jolivet PetscCheck(nileaves == iroffset[niranks], PETSC_COMM_SELF, PETSC_ERR_PLIB, "nileaves != iroffset[niranks])"); 1319566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, °ree)); 1329566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, °ree)); 1339566063dSJacob Faibussowitsch PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, degree, NULL, &mRootsOrigNumbering)); 1349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nileaves, irmine1)); 135f84a5eb8SVaclav Hapla for (r = 0; r < niranks; r++) { 136f84a5eb8SVaclav Hapla o = iroffset[r]; 137f84a5eb8SVaclav Hapla n = iroffset[r + 1] - o; 138f84a5eb8SVaclav Hapla for (i = 0; i < n; i++) (*irmine1)[o + i] = mRootsOrigNumbering[irmine[o + i]]; 139f84a5eb8SVaclav Hapla } 1409566063dSJacob Faibussowitsch PetscCall(PetscFree(mRootsOrigNumbering)); 141f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 142f84a5eb8SVaclav Hapla } 143f84a5eb8SVaclav Hapla 144124f9872SVaclav Hapla /*@ 145a8432d5bSVaclav Hapla DMPlexCheckInterfaceCones - Check that points on inter-partition interfaces have conforming order of cone points. 146124f9872SVaclav Hapla 147124f9872SVaclav Hapla Input Parameters: 148*a1cb98faSBarry Smith . dm - The `DMPLEX` object 149*a1cb98faSBarry Smith 150*a1cb98faSBarry Smith Level: developer 151124f9872SVaclav Hapla 152a8432d5bSVaclav Hapla Notes: 153a8432d5bSVaclav Hapla For example, if there is an edge (rank,index)=(0,2) connecting points cone(0,2)=[(0,0),(0,1)] in this order, and the point SF containts connections 0 <- (1,0), 1 <- (1,1) and 2 <- (1,2), 154a8432d5bSVaclav Hapla then this check would pass if the edge (1,2) has cone(1,2)=[(1,0),(1,1)]. By contrast, if cone(1,2)=[(1,1),(1,0)], then this check would fail. 155124f9872SVaclav Hapla 156*a1cb98faSBarry Smith This is mainly intended for debugging/testing purposes. Does not check cone orientation, for this purpose use `DMPlexCheckFaces()`. 157a8432d5bSVaclav Hapla 158*a1cb98faSBarry Smith For the complete list of DMPlexCheck* functions, see `DMSetFromOptions()`. 15995eb5ee5SVaclav Hapla 160a8432d5bSVaclav Hapla Developer Note: 161a8432d5bSVaclav Hapla Interface cones are expanded into vertices and then their coordinates are compared. 162124f9872SVaclav Hapla 163*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `DMPlexGetCone()`, `DMPlexGetConeSize()`, `DMGetPointSF()`, `DMGetCoordinates()`, `DMSetFromOptions()` 164124f9872SVaclav Hapla @*/ 165d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCheckInterfaceCones(DM dm) 166d71ae5a4SJacob Faibussowitsch { 167f84a5eb8SVaclav Hapla PetscSF sf; 168f84a5eb8SVaclav Hapla PetscInt nleaves, nranks, nroots; 169f84a5eb8SVaclav Hapla const PetscInt *mine, *roffset, *rmine, *rremote; 170f84a5eb8SVaclav Hapla const PetscSFNode *remote; 171f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 172f84a5eb8SVaclav Hapla PetscSF msf, imsf; 173f84a5eb8SVaclav Hapla PetscInt nileaves, niranks; 174f84a5eb8SVaclav Hapla const PetscMPIInt *iranks; 175f84a5eb8SVaclav Hapla const PetscInt *iroffset, *irmine, *irremote; 176f84a5eb8SVaclav Hapla PetscInt *rmine1, *rremote1; /* rmine and rremote copies simultaneously sorted by rank and rremote */ 177f84a5eb8SVaclav Hapla PetscInt *mine_orig_numbering; 178f84a5eb8SVaclav Hapla Vec *sntCoordinatesPerRank; 179f84a5eb8SVaclav Hapla Vec *refCoordinatesPerRank; 180ea78f98cSLisandro Dalcin Vec *recCoordinatesPerRank = NULL; 181f84a5eb8SVaclav Hapla PetscInt r; 182f84a5eb8SVaclav Hapla PetscMPIInt commsize, myrank; 183f84a5eb8SVaclav Hapla PetscBool same; 1848f2c89e7SVaclav Hapla PetscBool verbose = PETSC_FALSE; 185f84a5eb8SVaclav Hapla MPI_Comm comm; 186f84a5eb8SVaclav Hapla 187f84a5eb8SVaclav Hapla PetscFunctionBegin; 18810b92ba9SVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1899566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &myrank)); 1919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &commsize)); 192f84a5eb8SVaclav Hapla if (commsize < 2) PetscFunctionReturn(0); 1939566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 194f84a5eb8SVaclav Hapla if (!sf) PetscFunctionReturn(0); 1959566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &mine, &remote)); 196f84a5eb8SVaclav Hapla if (nroots < 0) PetscFunctionReturn(0); 1976858538eSMatthew G. Knepley PetscCheck(dm->coordinates[0].x || dm->coordinates[0].xl, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM coordinates must be set"); 1989566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 1999566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote)); 200f84a5eb8SVaclav Hapla 201f84a5eb8SVaclav Hapla /* Expand sent cones per rank */ 2029566063dSJacob Faibussowitsch PetscCall(SortByRemote_Private(sf, &rmine1, &rremote1)); 2039566063dSJacob Faibussowitsch PetscCall(GetRecursiveConeCoordinatesPerRank_Private(dm, sf, rmine1, &sntCoordinatesPerRank)); 204f84a5eb8SVaclav Hapla 205f84a5eb8SVaclav Hapla /* Create inverse SF */ 2069566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &msf)); 2079566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(msf, &imsf)); 2089566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(imsf)); 2099566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL)); 2109566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks(imsf, &niranks, &iranks, &iroffset, &irmine, &irremote)); 211f84a5eb8SVaclav Hapla 212f84a5eb8SVaclav Hapla /* Compute original numbering of multi-roots (referenced points) */ 2139566063dSJacob Faibussowitsch PetscCall(PetscSFComputeMultiRootOriginalNumberingByRank_Private(sf, imsf, &mine_orig_numbering)); 214f84a5eb8SVaclav Hapla 215124f9872SVaclav Hapla /* Expand coordinates of the referred cones per rank */ 2169566063dSJacob Faibussowitsch PetscCall(GetRecursiveConeCoordinatesPerRank_Private(dm, imsf, mine_orig_numbering, &refCoordinatesPerRank)); 217f84a5eb8SVaclav Hapla 218f84a5eb8SVaclav Hapla /* Send the coordinates */ 2199566063dSJacob Faibussowitsch PetscCall(ExchangeVecByRank_Private((PetscObject)sf, nranks, ranks, sntCoordinatesPerRank, niranks, iranks, &recCoordinatesPerRank)); 220f84a5eb8SVaclav Hapla 2218f2c89e7SVaclav Hapla /* verbose output */ 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_check_cones_conform_on_interfaces_verbose", &verbose, NULL)); 2238f2c89e7SVaclav Hapla if (verbose) { 2249f27a777SBarry Smith PetscViewer sv, v = PETSC_VIEWER_STDOUT_WORLD; 2259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "============\nDMPlexCheckInterfaceCones output\n============\n")); 2269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(v)); 2279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", myrank)); 2288f2c89e7SVaclav Hapla for (r = 0; r < nranks; r++) { 22963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, " r=%" PetscInt_FMT " ranks[r]=%d sntCoordinatesPerRank[r]:\n", r, ranks[r])); 2309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(v, PETSC_COMM_SELF, &sv)); 2329566063dSJacob Faibussowitsch PetscCall(VecView(sntCoordinatesPerRank[r], sv)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(v, PETSC_COMM_SELF, &sv)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 2358f2c89e7SVaclav Hapla } 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, " ----------\n")); 2378f2c89e7SVaclav Hapla for (r = 0; r < niranks; r++) { 23863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, " r=%" PetscInt_FMT " iranks[r]=%d refCoordinatesPerRank[r]:\n", r, iranks[r])); 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(v, PETSC_COMM_SELF, &sv)); 2419566063dSJacob Faibussowitsch PetscCall(VecView(refCoordinatesPerRank[r], sv)); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(v, PETSC_COMM_SELF, &sv)); 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 2448f2c89e7SVaclav Hapla } 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, " ----------\n")); 2468f2c89e7SVaclav Hapla for (r = 0; r < niranks; r++) { 24763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(v, " r=%" PetscInt_FMT " iranks[r]=%d recCoordinatesPerRank[r]:\n", r, iranks[r])); 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(v, PETSC_COMM_SELF, &sv)); 2509566063dSJacob Faibussowitsch PetscCall(VecView(recCoordinatesPerRank[r], sv)); 2519566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(v, PETSC_COMM_SELF, &sv)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 2538f2c89e7SVaclav Hapla } 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(v)); 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(v)); 2568f2c89e7SVaclav Hapla } 2578f2c89e7SVaclav Hapla 258f84a5eb8SVaclav Hapla /* Compare recCoordinatesPerRank with refCoordinatesPerRank */ 259f84a5eb8SVaclav Hapla for (r = 0; r < niranks; r++) { 2609566063dSJacob Faibussowitsch PetscCall(VecEqual(refCoordinatesPerRank[r], recCoordinatesPerRank[r], &same)); 26128b400f6SJacob Faibussowitsch PetscCheck(same, PETSC_COMM_SELF, PETSC_ERR_PLIB, "interface cones do not conform for remote rank %d", iranks[r]); 262f84a5eb8SVaclav Hapla } 263f84a5eb8SVaclav Hapla 264f84a5eb8SVaclav Hapla /* destroy sent stuff */ 26548a46eb9SPierre Jolivet for (r = 0; r < nranks; r++) PetscCall(VecDestroy(&sntCoordinatesPerRank[r])); 2669566063dSJacob Faibussowitsch PetscCall(PetscFree(sntCoordinatesPerRank)); 2679566063dSJacob Faibussowitsch PetscCall(PetscFree2(rmine1, rremote1)); 2689566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&imsf)); 269f84a5eb8SVaclav Hapla 270f84a5eb8SVaclav Hapla /* destroy referenced stuff */ 27148a46eb9SPierre Jolivet for (r = 0; r < niranks; r++) PetscCall(VecDestroy(&refCoordinatesPerRank[r])); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(refCoordinatesPerRank)); 2739566063dSJacob Faibussowitsch PetscCall(PetscFree(mine_orig_numbering)); 274f84a5eb8SVaclav Hapla 275f84a5eb8SVaclav Hapla /* destroy received stuff */ 27648a46eb9SPierre Jolivet for (r = 0; r < niranks; r++) PetscCall(VecDestroy(&recCoordinatesPerRank[r])); 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(recCoordinatesPerRank)); 278f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 279f84a5eb8SVaclav Hapla } 280