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 ? */ 6f84a5eb8SVaclav Hapla 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[]) 7f84a5eb8SVaclav Hapla { 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; 16*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Type_size(dt, &unitsize)); 17*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm(obj, &comm)); 18*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nrranks, &rsize, nrranks, &rarr)); 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nrranks, &rreq, nsranks, &sreq)); 20f84a5eb8SVaclav Hapla /* exchange array size */ 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetNewTag(obj,&tag)); 22f84a5eb8SVaclav Hapla for (r=0; r<nrranks; r++) { 23*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Irecv(&rsize[r], 1, MPIU_INT, rranks[r], tag, comm, &rreq[r])); 24f84a5eb8SVaclav Hapla } 25f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 26*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Isend(&ssize[r], 1, MPIU_INT, sranks[r], tag, comm, &sreq[r])); 27f84a5eb8SVaclav Hapla } 28*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE)); 29*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE)); 30f84a5eb8SVaclav Hapla /* exchange array */ 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetNewTag(obj,&tag)); 32f84a5eb8SVaclav Hapla for (r=0; r<nrranks; r++) { 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc(rsize[r]*unitsize, &rarr[r])); 34*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Irecv(rarr[r], rsize[r], dt, rranks[r], tag, comm, &rreq[r])); 35f84a5eb8SVaclav Hapla } 36f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 37*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Isend(sarr[r], ssize[r], dt, sranks[r], tag, comm, &sreq[r])); 38f84a5eb8SVaclav Hapla } 39*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE)); 40*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(rreq, sreq)); 42f84a5eb8SVaclav Hapla *rsize_out = rsize; 43f84a5eb8SVaclav Hapla *rarr_out = rarr; 44f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 45f84a5eb8SVaclav Hapla } 46f84a5eb8SVaclav Hapla 47f84a5eb8SVaclav Hapla /* TODO VecExchangeBegin/End */ 48f84a5eb8SVaclav Hapla /* TODO move to API ? */ 49f84a5eb8SVaclav Hapla static PetscErrorCode ExchangeVecByRank_Private(PetscObject obj, PetscInt nsranks, const PetscMPIInt sranks[], Vec svecs[], PetscInt nrranks, const PetscMPIInt rranks[], Vec *rvecs[]) 50f84a5eb8SVaclav Hapla { 51f84a5eb8SVaclav Hapla PetscInt r; 52f84a5eb8SVaclav Hapla PetscInt *ssize, *rsize; 53f84a5eb8SVaclav Hapla PetscScalar **rarr; 54f84a5eb8SVaclav Hapla const PetscScalar **sarr; 55f84a5eb8SVaclav Hapla Vec *rvecs_; 56f84a5eb8SVaclav Hapla MPI_Request *sreq, *rreq; 57f84a5eb8SVaclav Hapla 58f84a5eb8SVaclav Hapla PetscFunctionBegin; 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(nsranks, &ssize, nsranks, &sarr, nrranks, &rreq, nsranks, &sreq)); 60f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(svecs[r], &ssize[r])); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(svecs[r], &sarr[r])); 63f84a5eb8SVaclav Hapla } 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(ExchangeArrayByRank_Private(obj, MPIU_SCALAR, nsranks, sranks, ssize, (const void**)sarr, nrranks, rranks, &rsize, (void***)&rarr)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrranks, &rvecs_)); 66f84a5eb8SVaclav Hapla for (r=0; r<nrranks; r++) { 67f84a5eb8SVaclav Hapla /* set array in two steps to mimic PETSC_OWN_POINTER */ 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, rsize[r], NULL, &rvecs_[r])); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecReplaceArray(rvecs_[r], rarr[r])); 70f84a5eb8SVaclav Hapla } 71f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(svecs[r], &sarr[r])); 73f84a5eb8SVaclav Hapla } 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(rsize, rarr)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree4(ssize, sarr, rreq, sreq)); 76f84a5eb8SVaclav Hapla *rvecs = rvecs_; 77f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 78f84a5eb8SVaclav Hapla } 79f84a5eb8SVaclav Hapla 80f84a5eb8SVaclav Hapla static PetscErrorCode SortByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[]) 81f84a5eb8SVaclav Hapla { 82f84a5eb8SVaclav Hapla PetscInt nleaves; 83f84a5eb8SVaclav Hapla PetscInt nranks; 84f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 85f84a5eb8SVaclav Hapla const PetscInt *roffset, *rmine, *rremote; 86f84a5eb8SVaclav Hapla PetscInt n, o, r; 87f84a5eb8SVaclav Hapla 88f84a5eb8SVaclav Hapla PetscFunctionBegin; 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote)); 90f84a5eb8SVaclav Hapla nleaves = roffset[nranks]; 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nleaves, rmine1, nleaves, rremote1)); 92f84a5eb8SVaclav Hapla for (r=0; r<nranks; r++) { 93f84a5eb8SVaclav Hapla /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote 94f84a5eb8SVaclav Hapla - to unify order with the other side */ 95f84a5eb8SVaclav Hapla o = roffset[r]; 96f84a5eb8SVaclav Hapla n = roffset[r+1] - o; 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(&(*rmine1)[o], &rmine[o], n)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(&(*rremote1)[o], &rremote[o], n)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o])); 100f84a5eb8SVaclav Hapla } 101f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 102f84a5eb8SVaclav Hapla } 103f84a5eb8SVaclav Hapla 104f84a5eb8SVaclav Hapla static PetscErrorCode GetRecursiveConeCoordinatesPerRank_Private(DM dm, PetscSF sf, PetscInt rmine[], Vec *coordinatesPerRank[]) 105f84a5eb8SVaclav Hapla { 106f84a5eb8SVaclav Hapla IS pointsPerRank, conesPerRank; 107f84a5eb8SVaclav Hapla PetscInt nranks; 108f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 109f84a5eb8SVaclav Hapla const PetscInt *roffset; 110f84a5eb8SVaclav Hapla PetscInt n, o, r; 111f84a5eb8SVaclav Hapla 112f84a5eb8SVaclav Hapla PetscFunctionBegin; 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocalSetUp(dm)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, NULL, NULL)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nranks, coordinatesPerRank)); 116f84a5eb8SVaclav Hapla for (r=0; r<nranks; r++) { 117f84a5eb8SVaclav Hapla o = roffset[r]; 118f84a5eb8SVaclav Hapla n = roffset[r+1] - o; 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF, n, &rmine[o], PETSC_USE_POINTER, &pointsPerRank)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeRecursiveVertices(dm, pointsPerRank, &conesPerRank)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocalTuple(dm, conesPerRank, NULL, &(*coordinatesPerRank)[r])); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&pointsPerRank)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&conesPerRank)); 124f84a5eb8SVaclav Hapla } 125f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 126f84a5eb8SVaclav Hapla } 127f84a5eb8SVaclav Hapla 128f84a5eb8SVaclav Hapla static PetscErrorCode PetscSFComputeMultiRootOriginalNumberingByRank_Private(PetscSF sf, PetscSF imsf, PetscInt *irmine1[]) 129f84a5eb8SVaclav Hapla { 130f84a5eb8SVaclav Hapla PetscInt *mRootsOrigNumbering; 131f84a5eb8SVaclav Hapla PetscInt nileaves, niranks; 132f84a5eb8SVaclav Hapla const PetscInt *iroffset, *irmine, *degree; 133f84a5eb8SVaclav Hapla PetscInt i, n, o, r; 134f84a5eb8SVaclav Hapla 135f84a5eb8SVaclav Hapla PetscFunctionBegin; 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetRootRanks(imsf, &niranks, NULL, &iroffset, &irmine, NULL)); 1382c71b3e2SJacob Faibussowitsch PetscCheckFalse(nileaves != iroffset[niranks],PETSC_COMM_SELF,PETSC_ERR_PLIB,"nileaves != iroffset[niranks])"); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFComputeDegreeBegin(sf, °ree)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFComputeDegreeEnd(sf, °ree)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFComputeMultiRootOriginalNumbering(sf, degree, NULL, &mRootsOrigNumbering)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nileaves, irmine1)); 143f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 144f84a5eb8SVaclav Hapla o = iroffset[r]; 145f84a5eb8SVaclav Hapla n = iroffset[r+1] - o; 146f84a5eb8SVaclav Hapla for (i=0; i<n; i++) (*irmine1)[o+i] = mRootsOrigNumbering[irmine[o+i]]; 147f84a5eb8SVaclav Hapla } 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mRootsOrigNumbering)); 149f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 150f84a5eb8SVaclav Hapla } 151f84a5eb8SVaclav Hapla 152124f9872SVaclav Hapla /*@ 153a8432d5bSVaclav Hapla DMPlexCheckInterfaceCones - Check that points on inter-partition interfaces have conforming order of cone points. 154124f9872SVaclav Hapla 155124f9872SVaclav Hapla Input Parameters: 156124f9872SVaclav Hapla . dm - The DMPlex object 157124f9872SVaclav Hapla 158a8432d5bSVaclav Hapla Notes: 159a8432d5bSVaclav 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), 160a8432d5bSVaclav 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. 161124f9872SVaclav Hapla 162a8432d5bSVaclav Hapla This is mainly intended for debugging/testing purposes. Does not check cone orientation, for this purpose use DMPlexCheckFaces(). 163a8432d5bSVaclav Hapla 16495eb5ee5SVaclav Hapla For the complete list of DMPlexCheck* functions, see DMSetFromOptions(). 16595eb5ee5SVaclav Hapla 166a8432d5bSVaclav Hapla Developer Note: 167a8432d5bSVaclav Hapla Interface cones are expanded into vertices and then their coordinates are compared. 168124f9872SVaclav Hapla 169124f9872SVaclav Hapla Level: developer 170124f9872SVaclav Hapla 17195eb5ee5SVaclav Hapla .seealso: DMPlexGetCone(), DMPlexGetConeSize(), DMGetPointSF(), DMGetCoordinates(), DMSetFromOptions() 172124f9872SVaclav Hapla @*/ 173a8432d5bSVaclav Hapla PetscErrorCode DMPlexCheckInterfaceCones(DM dm) 174f84a5eb8SVaclav Hapla { 175f84a5eb8SVaclav Hapla PetscSF sf; 176f84a5eb8SVaclav Hapla PetscInt nleaves, nranks, nroots; 177f84a5eb8SVaclav Hapla const PetscInt *mine, *roffset, *rmine, *rremote; 178f84a5eb8SVaclav Hapla const PetscSFNode *remote; 179f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 180f84a5eb8SVaclav Hapla PetscSF msf, imsf; 181f84a5eb8SVaclav Hapla PetscInt nileaves, niranks; 182f84a5eb8SVaclav Hapla const PetscMPIInt *iranks; 183f84a5eb8SVaclav Hapla const PetscInt *iroffset, *irmine, *irremote; 184f84a5eb8SVaclav Hapla PetscInt *rmine1, *rremote1; /* rmine and rremote copies simultaneously sorted by rank and rremote */ 185f84a5eb8SVaclav Hapla PetscInt *mine_orig_numbering; 186f84a5eb8SVaclav Hapla Vec *sntCoordinatesPerRank; 187f84a5eb8SVaclav Hapla Vec *refCoordinatesPerRank; 188ea78f98cSLisandro Dalcin Vec *recCoordinatesPerRank=NULL; 189f84a5eb8SVaclav Hapla PetscInt r; 190f84a5eb8SVaclav Hapla PetscMPIInt commsize, myrank; 191f84a5eb8SVaclav Hapla PetscBool same; 1928f2c89e7SVaclav Hapla PetscBool verbose=PETSC_FALSE; 193f84a5eb8SVaclav Hapla MPI_Comm comm; 194f84a5eb8SVaclav Hapla 195f84a5eb8SVaclav Hapla PetscFunctionBegin; 19610b92ba9SVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm)); 198*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &myrank)); 199*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &commsize)); 200f84a5eb8SVaclav Hapla if (commsize < 2) PetscFunctionReturn(0); 201*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetPointSF(dm, &sf)); 202f84a5eb8SVaclav Hapla if (!sf) PetscFunctionReturn(0); 203*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(sf, &nroots, &nleaves, &mine, &remote)); 204f84a5eb8SVaclav Hapla if (nroots < 0) PetscFunctionReturn(0); 2052c71b3e2SJacob Faibussowitsch PetscCheckFalse(!dm->coordinates && !dm->coordinatesLocal,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM coordinates must be set"); 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetUp(sf)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote)); 208f84a5eb8SVaclav Hapla 209f84a5eb8SVaclav Hapla /* Expand sent cones per rank */ 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(SortByRemote_Private(sf, &rmine1, &rremote1)); 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(GetRecursiveConeCoordinatesPerRank_Private(dm, sf, rmine1, &sntCoordinatesPerRank)); 212f84a5eb8SVaclav Hapla 213f84a5eb8SVaclav Hapla /* Create inverse SF */ 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetMultiSF(sf,&msf)); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFCreateInverseSF(msf,&imsf)); 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetUp(imsf)); 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL)); 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetRootRanks(imsf, &niranks, &iranks, &iroffset, &irmine, &irremote)); 219f84a5eb8SVaclav Hapla 220f84a5eb8SVaclav Hapla /* Compute original numbering of multi-roots (referenced points) */ 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFComputeMultiRootOriginalNumberingByRank_Private(sf, imsf, &mine_orig_numbering)); 222f84a5eb8SVaclav Hapla 223124f9872SVaclav Hapla /* Expand coordinates of the referred cones per rank */ 224*5f80ce2aSJacob Faibussowitsch CHKERRQ(GetRecursiveConeCoordinatesPerRank_Private(dm, imsf, mine_orig_numbering, &refCoordinatesPerRank)); 225f84a5eb8SVaclav Hapla 226f84a5eb8SVaclav Hapla /* Send the coordinates */ 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(ExchangeVecByRank_Private((PetscObject)sf, nranks, ranks, sntCoordinatesPerRank, niranks, iranks, &recCoordinatesPerRank)); 228f84a5eb8SVaclav Hapla 2298f2c89e7SVaclav Hapla /* verbose output */ 230*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_check_cones_conform_on_interfaces_verbose", &verbose, NULL)); 2318f2c89e7SVaclav Hapla if (verbose) { 2329f27a777SBarry Smith PetscViewer sv, v = PETSC_VIEWER_STDOUT_WORLD; 233*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(v, "============\nDMPlexCheckInterfaceCones output\n============\n")); 234*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(v)); 235*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", myrank)); 2368f2c89e7SVaclav Hapla for (r=0; r<nranks; r++) { 237*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, " r=%D ranks[r]=%d sntCoordinatesPerRank[r]:\n", r, ranks[r])); 238*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(v)); 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(v,PETSC_COMM_SELF,&sv)); 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(sntCoordinatesPerRank[r], sv)); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(v,PETSC_COMM_SELF,&sv)); 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(v)); 2438f2c89e7SVaclav Hapla } 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, " ----------\n")); 2458f2c89e7SVaclav Hapla for (r=0; r<niranks; r++) { 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, " r=%D iranks[r]=%d refCoordinatesPerRank[r]:\n", r, iranks[r])); 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(v)); 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(v,PETSC_COMM_SELF,&sv)); 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(refCoordinatesPerRank[r], sv)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(v,PETSC_COMM_SELF,&sv)); 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(v)); 2528f2c89e7SVaclav Hapla } 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, " ----------\n")); 2548f2c89e7SVaclav Hapla for (r=0; r<niranks; r++) { 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(v, " r=%D iranks[r]=%d recCoordinatesPerRank[r]:\n", r, iranks[r])); 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(v)); 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(v,PETSC_COMM_SELF,&sv)); 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(recCoordinatesPerRank[r], sv)); 259*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(v,PETSC_COMM_SELF,&sv)); 260*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(v)); 2618f2c89e7SVaclav Hapla } 262*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(v)); 263*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(v)); 2648f2c89e7SVaclav Hapla } 2658f2c89e7SVaclav Hapla 266f84a5eb8SVaclav Hapla /* Compare recCoordinatesPerRank with refCoordinatesPerRank */ 267f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 268*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(refCoordinatesPerRank[r], recCoordinatesPerRank[r], &same)); 2692c71b3e2SJacob Faibussowitsch PetscCheckFalse(!same,PETSC_COMM_SELF, PETSC_ERR_PLIB, "interface cones do not conform for remote rank %d", iranks[r]); 270f84a5eb8SVaclav Hapla } 271f84a5eb8SVaclav Hapla 272f84a5eb8SVaclav Hapla /* destroy sent stuff */ 273f84a5eb8SVaclav Hapla for (r=0; r<nranks; r++) { 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&sntCoordinatesPerRank[r])); 275f84a5eb8SVaclav Hapla } 276*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(sntCoordinatesPerRank)); 277*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(rmine1, rremote1)); 278*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&imsf)); 279f84a5eb8SVaclav Hapla 280f84a5eb8SVaclav Hapla /* destroy referenced stuff */ 281f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&refCoordinatesPerRank[r])); 283f84a5eb8SVaclav Hapla } 284*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(refCoordinatesPerRank)); 285*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(mine_orig_numbering)); 286f84a5eb8SVaclav Hapla 287f84a5eb8SVaclav Hapla /* destroy received stuff */ 288f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&recCoordinatesPerRank[r])); 290f84a5eb8SVaclav Hapla } 291*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(recCoordinatesPerRank)); 292f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 293f84a5eb8SVaclav Hapla } 294