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 PetscErrorCode ierr; 15f84a5eb8SVaclav Hapla 16f84a5eb8SVaclav Hapla PetscFunctionBegin; 17f84a5eb8SVaclav Hapla ierr = MPI_Type_size(dt, &unitsize);CHKERRQ(ierr); 18f84a5eb8SVaclav Hapla ierr = PetscObjectGetComm(obj, &comm);CHKERRQ(ierr); 19f84a5eb8SVaclav Hapla ierr = PetscMalloc2(nrranks, &rsize, nrranks, &rarr);CHKERRQ(ierr); 20f84a5eb8SVaclav Hapla ierr = PetscMalloc2(nrranks, &rreq, nsranks, &sreq);CHKERRQ(ierr); 21f84a5eb8SVaclav Hapla /* exchange array size */ 22f84a5eb8SVaclav Hapla ierr = PetscObjectGetNewTag(obj,&tag);CHKERRQ(ierr); 23f84a5eb8SVaclav Hapla for (r=0; r<nrranks; r++) { 24f84a5eb8SVaclav Hapla ierr = MPI_Irecv(&rsize[r], 1, MPIU_INT, rranks[r], tag, comm, &rreq[r]);CHKERRQ(ierr); 25f84a5eb8SVaclav Hapla } 26f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 27f84a5eb8SVaclav Hapla ierr = MPI_Isend(&ssize[r], 1, MPIU_INT, sranks[r], tag, comm, &sreq[r]);CHKERRQ(ierr); 28f84a5eb8SVaclav Hapla } 29f84a5eb8SVaclav Hapla ierr = MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE);CHKERRQ(ierr); 30f84a5eb8SVaclav Hapla /* exchange array */ 31f84a5eb8SVaclav Hapla ierr = PetscObjectGetNewTag(obj,&tag);CHKERRQ(ierr); 32f84a5eb8SVaclav Hapla for (r=0; r<nrranks; r++) { 33f84a5eb8SVaclav Hapla ierr = PetscMalloc(rsize[r]*unitsize, &rarr[r]);CHKERRQ(ierr); 34f84a5eb8SVaclav Hapla ierr = MPI_Irecv(rarr[r], rsize[r], dt, rranks[r], tag, comm, &rreq[r]);CHKERRQ(ierr); 35f84a5eb8SVaclav Hapla } 36f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 37f84a5eb8SVaclav Hapla ierr = MPI_Isend(sarr[r], ssize[r], dt, sranks[r], tag, comm, &sreq[r]);CHKERRQ(ierr); 38f84a5eb8SVaclav Hapla } 39f84a5eb8SVaclav Hapla ierr = MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE);CHKERRQ(ierr); 40f84a5eb8SVaclav Hapla ierr = MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE);CHKERRQ(ierr); 41f84a5eb8SVaclav Hapla ierr = PetscFree2(rreq, sreq);CHKERRQ(ierr); 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 PetscErrorCode ierr; 58f84a5eb8SVaclav Hapla 59f84a5eb8SVaclav Hapla PetscFunctionBegin; 60f84a5eb8SVaclav Hapla ierr = PetscMalloc4(nsranks, &ssize, nsranks, &sarr, nrranks, &rreq, nsranks, &sreq);CHKERRQ(ierr); 61f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 62f84a5eb8SVaclav Hapla ierr = VecGetLocalSize(svecs[r], &ssize[r]);CHKERRQ(ierr); 63f84a5eb8SVaclav Hapla ierr = VecGetArrayRead(svecs[r], &sarr[r]);CHKERRQ(ierr); 64f84a5eb8SVaclav Hapla } 65f84a5eb8SVaclav Hapla ierr = ExchangeArrayByRank_Private(obj, MPIU_SCALAR, nsranks, sranks, ssize, (const void**)sarr, nrranks, rranks, &rsize, (void***)&rarr);CHKERRQ(ierr); 66f84a5eb8SVaclav Hapla ierr = PetscMalloc1(nrranks, &rvecs_);CHKERRQ(ierr); 67f84a5eb8SVaclav Hapla for (r=0; r<nrranks; r++) { 68f84a5eb8SVaclav Hapla /* set array in two steps to mimic PETSC_OWN_POINTER */ 69f84a5eb8SVaclav Hapla ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, 1, rsize[r], NULL, &rvecs_[r]);CHKERRQ(ierr); 70f84a5eb8SVaclav Hapla ierr = VecReplaceArray(rvecs_[r], rarr[r]);CHKERRQ(ierr); 71f84a5eb8SVaclav Hapla } 72f84a5eb8SVaclav Hapla for (r=0; r<nsranks; r++) { 73f84a5eb8SVaclav Hapla ierr = VecRestoreArrayRead(svecs[r], &sarr[r]);CHKERRQ(ierr); 74f84a5eb8SVaclav Hapla } 75f84a5eb8SVaclav Hapla ierr = PetscFree2(rsize, rarr);CHKERRQ(ierr); 76f84a5eb8SVaclav Hapla ierr = PetscFree4(ssize, sarr, rreq, sreq);CHKERRQ(ierr); 77f84a5eb8SVaclav Hapla *rvecs = rvecs_; 78f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 79f84a5eb8SVaclav Hapla } 80f84a5eb8SVaclav Hapla 81f84a5eb8SVaclav Hapla static PetscErrorCode SortByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[]) 82f84a5eb8SVaclav Hapla { 83f84a5eb8SVaclav Hapla PetscInt nleaves; 84f84a5eb8SVaclav Hapla PetscInt nranks; 85f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 86f84a5eb8SVaclav Hapla const PetscInt *roffset, *rmine, *rremote; 87f84a5eb8SVaclav Hapla PetscInt n, o, r; 88f84a5eb8SVaclav Hapla PetscErrorCode ierr; 89f84a5eb8SVaclav Hapla 90f84a5eb8SVaclav Hapla PetscFunctionBegin; 91f84a5eb8SVaclav Hapla ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr); 92f84a5eb8SVaclav Hapla nleaves = roffset[nranks]; 93f84a5eb8SVaclav Hapla ierr = PetscMalloc2(nleaves, rmine1, nleaves, rremote1);CHKERRQ(ierr); 94f84a5eb8SVaclav Hapla for (r=0; r<nranks; r++) { 95f84a5eb8SVaclav Hapla /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote 96f84a5eb8SVaclav Hapla - to unify order with the other side */ 97f84a5eb8SVaclav Hapla o = roffset[r]; 98f84a5eb8SVaclav Hapla n = roffset[r+1] - o; 99f84a5eb8SVaclav Hapla ierr = PetscMemcpy(&(*rmine1)[o], &rmine[o], n*sizeof(PetscInt));CHKERRQ(ierr); 100f84a5eb8SVaclav Hapla ierr = PetscMemcpy(&(*rremote1)[o], &rremote[o], n*sizeof(PetscInt));CHKERRQ(ierr); 101f84a5eb8SVaclav Hapla ierr = PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);CHKERRQ(ierr); 102f84a5eb8SVaclav Hapla } 103f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 104f84a5eb8SVaclav Hapla } 105f84a5eb8SVaclav Hapla 106f84a5eb8SVaclav Hapla static PetscErrorCode GetRecursiveConeCoordinatesPerRank_Private(DM dm, PetscSF sf, PetscInt rmine[], Vec *coordinatesPerRank[]) 107f84a5eb8SVaclav Hapla { 108f84a5eb8SVaclav Hapla IS pointsPerRank, conesPerRank; 109f84a5eb8SVaclav Hapla PetscInt nranks; 110f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 111f84a5eb8SVaclav Hapla const PetscInt *roffset; 112f84a5eb8SVaclav Hapla PetscInt n, o, r; 113f84a5eb8SVaclav Hapla PetscErrorCode ierr; 114f84a5eb8SVaclav Hapla 115f84a5eb8SVaclav Hapla PetscFunctionBegin; 116f84a5eb8SVaclav Hapla ierr = DMGetCoordinatesLocalSetUp(dm);CHKERRQ(ierr); 117f84a5eb8SVaclav Hapla ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);CHKERRQ(ierr); 118f84a5eb8SVaclav Hapla ierr = PetscMalloc1(nranks, coordinatesPerRank);CHKERRQ(ierr); 119f84a5eb8SVaclav Hapla for (r=0; r<nranks; r++) { 120f84a5eb8SVaclav Hapla o = roffset[r]; 121f84a5eb8SVaclav Hapla n = roffset[r+1] - o; 122f84a5eb8SVaclav Hapla ierr = ISCreateGeneral(PETSC_COMM_SELF, n, &rmine[o], PETSC_USE_POINTER, &pointsPerRank);CHKERRQ(ierr); 123f84a5eb8SVaclav Hapla ierr = DMPlexGetConeRecursive(dm, pointsPerRank, &conesPerRank);CHKERRQ(ierr); 124f84a5eb8SVaclav Hapla ierr = DMGetCoordinatesLocalTuple(dm, conesPerRank, NULL, &(*coordinatesPerRank)[r]);CHKERRQ(ierr); 125f84a5eb8SVaclav Hapla ierr = ISDestroy(&pointsPerRank);CHKERRQ(ierr); 126f84a5eb8SVaclav Hapla ierr = ISDestroy(&conesPerRank);CHKERRQ(ierr); 127f84a5eb8SVaclav Hapla } 128f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 129f84a5eb8SVaclav Hapla } 130f84a5eb8SVaclav Hapla 131f84a5eb8SVaclav Hapla static PetscErrorCode PetscSFComputeMultiRootOriginalNumberingByRank_Private(PetscSF sf, PetscSF imsf, PetscInt *irmine1[]) 132f84a5eb8SVaclav Hapla { 133f84a5eb8SVaclav Hapla PetscInt *mRootsOrigNumbering; 134f84a5eb8SVaclav Hapla PetscInt nileaves, niranks; 135f84a5eb8SVaclav Hapla const PetscInt *iroffset, *irmine, *degree; 136f84a5eb8SVaclav Hapla PetscInt i, n, o, r; 137f84a5eb8SVaclav Hapla PetscErrorCode ierr; 138f84a5eb8SVaclav Hapla 139f84a5eb8SVaclav Hapla PetscFunctionBegin; 140f84a5eb8SVaclav Hapla ierr = PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL);CHKERRQ(ierr); 141f84a5eb8SVaclav Hapla ierr = PetscSFGetRanks(imsf, &niranks, NULL, &iroffset, &irmine, NULL);CHKERRQ(ierr); 142f84a5eb8SVaclav Hapla #if defined(PETSC_USE_DEBUG) 143f84a5eb8SVaclav Hapla if (PetscUnlikely(nileaves != iroffset[niranks])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nileaves != iroffset[niranks])"); 144f84a5eb8SVaclav Hapla #endif 145f84a5eb8SVaclav Hapla ierr = PetscSFComputeDegreeBegin(sf, °ree);CHKERRQ(ierr); 146f84a5eb8SVaclav Hapla ierr = PetscSFComputeDegreeEnd(sf, °ree);CHKERRQ(ierr); 147f84a5eb8SVaclav Hapla ierr = PetscSFComputeMultiRootOriginalNumbering(sf, degree, &mRootsOrigNumbering);CHKERRQ(ierr); 148f84a5eb8SVaclav Hapla ierr = PetscMalloc1(nileaves, irmine1);CHKERRQ(ierr); 149f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 150f84a5eb8SVaclav Hapla o = iroffset[r]; 151f84a5eb8SVaclav Hapla n = iroffset[r+1] - o; 152f84a5eb8SVaclav Hapla for (i=0; i<n; i++) (*irmine1)[o+i] = mRootsOrigNumbering[irmine[o+i]]; 153f84a5eb8SVaclav Hapla } 154f84a5eb8SVaclav Hapla ierr = PetscFree(mRootsOrigNumbering);CHKERRQ(ierr); 155f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 156f84a5eb8SVaclav Hapla } 157f84a5eb8SVaclav Hapla 158124f9872SVaclav Hapla /*@ 159124f9872SVaclav Hapla DMPlexCheckConesConformOnInterfaces - Check that points on inter-partition interfaces have conforming order of cone points. 160124f9872SVaclav 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), 161124f9872SVaclav 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. 162124f9872SVaclav Hapla 163124f9872SVaclav Hapla Input Parameters: 164124f9872SVaclav Hapla . dm - The DMPlex object 165124f9872SVaclav Hapla 166124f9872SVaclav Hapla Note: This is mainly intended for debugging/testing purposes. Does not check cone orientation, for this purpose use DMPlexCheckFaces(). 167124f9872SVaclav Hapla 168124f9872SVaclav Hapla Developer Note: Interface cones are expanded into vertices and then their coordinates are compared. 169124f9872SVaclav Hapla 170124f9872SVaclav Hapla Level: developer 171124f9872SVaclav Hapla 172124f9872SVaclav Hapla .seealso: DMPlexGetCone(), DMPlexGetConeSize(), DMGetPointSF(), DMGetCoordinates(), DMPlexCheckFaces(), DMPlexCheckPointSF(), DMPlexCheckSymmetry(), DMPlexCheckSkeleton() 173124f9872SVaclav Hapla @*/ 174f84a5eb8SVaclav Hapla PetscErrorCode DMPlexCheckConesConformOnInterfaces(DM dm) 175f84a5eb8SVaclav Hapla { 176f84a5eb8SVaclav Hapla PetscSF sf; 177f84a5eb8SVaclav Hapla PetscInt nleaves, nranks, nroots; 178f84a5eb8SVaclav Hapla const PetscInt *mine, *roffset, *rmine, *rremote; 179f84a5eb8SVaclav Hapla const PetscSFNode *remote; 180f84a5eb8SVaclav Hapla const PetscMPIInt *ranks; 181f84a5eb8SVaclav Hapla PetscSF msf, imsf; 182f84a5eb8SVaclav Hapla PetscInt nileaves, niranks; 183f84a5eb8SVaclav Hapla const PetscMPIInt *iranks; 184f84a5eb8SVaclav Hapla const PetscInt *iroffset, *irmine, *irremote; 185f84a5eb8SVaclav Hapla PetscInt *rmine1, *rremote1; /* rmine and rremote copies simultaneously sorted by rank and rremote */ 186f84a5eb8SVaclav Hapla PetscInt *mine_orig_numbering; 187f84a5eb8SVaclav Hapla Vec *sntCoordinatesPerRank; 188f84a5eb8SVaclav Hapla Vec *refCoordinatesPerRank; 189f84a5eb8SVaclav Hapla Vec *recCoordinatesPerRank; 190f84a5eb8SVaclav Hapla PetscInt r; 191f84a5eb8SVaclav Hapla PetscMPIInt commsize, myrank; 192f84a5eb8SVaclav Hapla PetscBool same; 193f84a5eb8SVaclav Hapla MPI_Comm comm; 194f84a5eb8SVaclav Hapla PetscErrorCode ierr; 195f84a5eb8SVaclav Hapla 196f84a5eb8SVaclav Hapla PetscFunctionBegin; 197*10b92ba9SVaclav Hapla PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 198f84a5eb8SVaclav Hapla ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 199f84a5eb8SVaclav Hapla ierr = MPI_Comm_rank(comm, &myrank);CHKERRQ(ierr); 200f84a5eb8SVaclav Hapla ierr = MPI_Comm_size(comm, &commsize);CHKERRQ(ierr); 201f84a5eb8SVaclav Hapla if (commsize < 2) PetscFunctionReturn(0); 202f84a5eb8SVaclav Hapla ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 203f84a5eb8SVaclav Hapla if (!sf) PetscFunctionReturn(0); 204f84a5eb8SVaclav Hapla ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &mine, &remote);CHKERRQ(ierr); 205f84a5eb8SVaclav Hapla if (nroots < 0) PetscFunctionReturn(0); 206*10b92ba9SVaclav Hapla if (!dm->coordinates && !dm->coordinatesLocal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM coordinates must be set"); 207f84a5eb8SVaclav Hapla ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 208f84a5eb8SVaclav Hapla ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr); 209f84a5eb8SVaclav Hapla 210f84a5eb8SVaclav Hapla /* Expand sent cones per rank */ 211f84a5eb8SVaclav Hapla ierr = SortByRemote_Private(sf, &rmine1, &rremote1);CHKERRQ(ierr); 212f84a5eb8SVaclav Hapla ierr = GetRecursiveConeCoordinatesPerRank_Private(dm, sf, rmine1, &sntCoordinatesPerRank);CHKERRQ(ierr); 213f84a5eb8SVaclav Hapla 214f84a5eb8SVaclav Hapla /* Create inverse SF */ 215f84a5eb8SVaclav Hapla ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr); 216f84a5eb8SVaclav Hapla ierr = PetscSFCreateInverseSF(msf,&imsf);CHKERRQ(ierr); 217f84a5eb8SVaclav Hapla ierr = PetscSFSetUp(imsf);CHKERRQ(ierr); 218f84a5eb8SVaclav Hapla ierr = PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL);CHKERRQ(ierr); 219f84a5eb8SVaclav Hapla ierr = PetscSFGetRanks(imsf, &niranks, &iranks, &iroffset, &irmine, &irremote);CHKERRQ(ierr); 220f84a5eb8SVaclav Hapla 221f84a5eb8SVaclav Hapla /* Compute original numbering of multi-roots (referenced points) */ 222f84a5eb8SVaclav Hapla ierr = PetscSFComputeMultiRootOriginalNumberingByRank_Private(sf, imsf, &mine_orig_numbering);CHKERRQ(ierr); 223f84a5eb8SVaclav Hapla 224124f9872SVaclav Hapla /* Expand coordinates of the referred cones per rank */ 225f84a5eb8SVaclav Hapla ierr = GetRecursiveConeCoordinatesPerRank_Private(dm, imsf, mine_orig_numbering, &refCoordinatesPerRank);CHKERRQ(ierr); 226f84a5eb8SVaclav Hapla 227f84a5eb8SVaclav Hapla /* Send the coordinates */ 228f84a5eb8SVaclav Hapla ierr = ExchangeVecByRank_Private((PetscObject)sf, nranks, ranks, sntCoordinatesPerRank, niranks, iranks, &recCoordinatesPerRank);CHKERRQ(ierr); 229f84a5eb8SVaclav Hapla 230f84a5eb8SVaclav Hapla /* Compare recCoordinatesPerRank with refCoordinatesPerRank */ 231f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 232f84a5eb8SVaclav Hapla ierr = VecEqual(refCoordinatesPerRank[r], recCoordinatesPerRank[r], &same);CHKERRQ(ierr); 233f84a5eb8SVaclav Hapla if (!same) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "interface cones do not conform for remote rank %d", iranks[r]); 234f84a5eb8SVaclav Hapla } 235f84a5eb8SVaclav Hapla 236f84a5eb8SVaclav Hapla /* destroy sent stuff */ 237f84a5eb8SVaclav Hapla for (r=0; r<nranks; r++) { 238f84a5eb8SVaclav Hapla ierr = VecDestroy(&sntCoordinatesPerRank[r]);CHKERRQ(ierr); 239f84a5eb8SVaclav Hapla } 240f84a5eb8SVaclav Hapla ierr = PetscFree(sntCoordinatesPerRank);CHKERRQ(ierr); 241f84a5eb8SVaclav Hapla ierr = PetscFree2(rmine1, rremote1);CHKERRQ(ierr); 242f84a5eb8SVaclav Hapla ierr = PetscSFDestroy(&imsf);CHKERRQ(ierr); 243f84a5eb8SVaclav Hapla 244f84a5eb8SVaclav Hapla /* destroy referenced stuff */ 245f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 246f84a5eb8SVaclav Hapla ierr = VecDestroy(&refCoordinatesPerRank[r]);CHKERRQ(ierr); 247f84a5eb8SVaclav Hapla } 248f84a5eb8SVaclav Hapla ierr = PetscFree(refCoordinatesPerRank);CHKERRQ(ierr); 249f84a5eb8SVaclav Hapla ierr = PetscFree(mine_orig_numbering);CHKERRQ(ierr); 250f84a5eb8SVaclav Hapla 251f84a5eb8SVaclav Hapla /* destroy received stuff */ 252f84a5eb8SVaclav Hapla for (r=0; r<niranks; r++) { 253f84a5eb8SVaclav Hapla ierr = VecDestroy(&recCoordinatesPerRank[r]);CHKERRQ(ierr); 254f84a5eb8SVaclav Hapla } 255f84a5eb8SVaclav Hapla ierr = PetscFree(recCoordinatesPerRank);CHKERRQ(ierr); 256f84a5eb8SVaclav Hapla PetscFunctionReturn(0); 257f84a5eb8SVaclav Hapla } 258