xref: /petsc/src/dm/impls/plex/plexcheckinterface.c (revision f84a5eb86a678fd1c5282af6f4d0ed72cbd71325)
1*f84a5eb8SVaclav Hapla #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*f84a5eb8SVaclav Hapla 
3*f84a5eb8SVaclav Hapla /* TODO PetscArrayExchangeBegin/End */
4*f84a5eb8SVaclav Hapla /* TODO blocksize */
5*f84a5eb8SVaclav Hapla /* TODO move to API ? */
6*f84a5eb8SVaclav 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[])
7*f84a5eb8SVaclav Hapla {
8*f84a5eb8SVaclav Hapla   PetscInt r;
9*f84a5eb8SVaclav Hapla   PetscInt *rsize;
10*f84a5eb8SVaclav Hapla   void **rarr;
11*f84a5eb8SVaclav Hapla   MPI_Request *sreq, *rreq;
12*f84a5eb8SVaclav Hapla   PetscMPIInt tag, unitsize;
13*f84a5eb8SVaclav Hapla   MPI_Comm comm;
14*f84a5eb8SVaclav Hapla   PetscErrorCode ierr;
15*f84a5eb8SVaclav Hapla 
16*f84a5eb8SVaclav Hapla   PetscFunctionBegin;
17*f84a5eb8SVaclav Hapla   ierr = MPI_Type_size(dt, &unitsize);CHKERRQ(ierr);
18*f84a5eb8SVaclav Hapla   ierr = PetscObjectGetComm(obj, &comm);CHKERRQ(ierr);
19*f84a5eb8SVaclav Hapla   ierr = PetscMalloc2(nrranks, &rsize, nrranks, &rarr);CHKERRQ(ierr);
20*f84a5eb8SVaclav Hapla   ierr = PetscMalloc2(nrranks, &rreq, nsranks, &sreq);CHKERRQ(ierr);
21*f84a5eb8SVaclav Hapla   /* exchange array size */
22*f84a5eb8SVaclav Hapla   ierr = PetscObjectGetNewTag(obj,&tag);CHKERRQ(ierr);
23*f84a5eb8SVaclav Hapla   for (r=0; r<nrranks; r++) {
24*f84a5eb8SVaclav Hapla     ierr = MPI_Irecv(&rsize[r], 1, MPIU_INT, rranks[r], tag, comm, &rreq[r]);CHKERRQ(ierr);
25*f84a5eb8SVaclav Hapla   }
26*f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
27*f84a5eb8SVaclav Hapla     ierr = MPI_Isend(&ssize[r], 1, MPIU_INT, sranks[r], tag, comm, &sreq[r]);CHKERRQ(ierr);
28*f84a5eb8SVaclav Hapla   }
29*f84a5eb8SVaclav Hapla   ierr = MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE);CHKERRQ(ierr);
30*f84a5eb8SVaclav Hapla   /* exchange array */
31*f84a5eb8SVaclav Hapla   ierr = PetscObjectGetNewTag(obj,&tag);CHKERRQ(ierr);
32*f84a5eb8SVaclav Hapla   for (r=0; r<nrranks; r++) {
33*f84a5eb8SVaclav Hapla     ierr = PetscMalloc(rsize[r]*unitsize, &rarr[r]);CHKERRQ(ierr);
34*f84a5eb8SVaclav Hapla     ierr = MPI_Irecv(rarr[r], rsize[r], dt, rranks[r], tag, comm, &rreq[r]);CHKERRQ(ierr);
35*f84a5eb8SVaclav Hapla   }
36*f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
37*f84a5eb8SVaclav Hapla     ierr = MPI_Isend(sarr[r], ssize[r], dt, sranks[r], tag, comm, &sreq[r]);CHKERRQ(ierr);
38*f84a5eb8SVaclav Hapla   }
39*f84a5eb8SVaclav Hapla   ierr = MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE);CHKERRQ(ierr);
40*f84a5eb8SVaclav Hapla   ierr = MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE);CHKERRQ(ierr);
41*f84a5eb8SVaclav Hapla   ierr = PetscFree2(rreq, sreq);CHKERRQ(ierr);
42*f84a5eb8SVaclav Hapla   *rsize_out = rsize;
43*f84a5eb8SVaclav Hapla   *rarr_out = rarr;
44*f84a5eb8SVaclav Hapla   PetscFunctionReturn(0);
45*f84a5eb8SVaclav Hapla }
46*f84a5eb8SVaclav Hapla 
47*f84a5eb8SVaclav Hapla /* TODO VecExchangeBegin/End */
48*f84a5eb8SVaclav Hapla /* TODO move to API ? */
49*f84a5eb8SVaclav Hapla static PetscErrorCode ExchangeVecByRank_Private(PetscObject obj, PetscInt nsranks, const PetscMPIInt sranks[], Vec svecs[], PetscInt nrranks, const PetscMPIInt rranks[], Vec *rvecs[])
50*f84a5eb8SVaclav Hapla {
51*f84a5eb8SVaclav Hapla   PetscInt r;
52*f84a5eb8SVaclav Hapla   PetscInt *ssize, *rsize;
53*f84a5eb8SVaclav Hapla   PetscScalar **rarr;
54*f84a5eb8SVaclav Hapla   const PetscScalar **sarr;
55*f84a5eb8SVaclav Hapla   Vec *rvecs_;
56*f84a5eb8SVaclav Hapla   MPI_Request *sreq, *rreq;
57*f84a5eb8SVaclav Hapla   PetscErrorCode ierr;
58*f84a5eb8SVaclav Hapla 
59*f84a5eb8SVaclav Hapla   PetscFunctionBegin;
60*f84a5eb8SVaclav Hapla   ierr = PetscMalloc4(nsranks, &ssize, nsranks, &sarr, nrranks, &rreq, nsranks, &sreq);CHKERRQ(ierr);
61*f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
62*f84a5eb8SVaclav Hapla     ierr = VecGetLocalSize(svecs[r], &ssize[r]);CHKERRQ(ierr);
63*f84a5eb8SVaclav Hapla     ierr = VecGetArrayRead(svecs[r], &sarr[r]);CHKERRQ(ierr);
64*f84a5eb8SVaclav Hapla   }
65*f84a5eb8SVaclav Hapla   ierr = ExchangeArrayByRank_Private(obj, MPIU_SCALAR, nsranks, sranks, ssize, (const void**)sarr, nrranks, rranks, &rsize, (void***)&rarr);CHKERRQ(ierr);
66*f84a5eb8SVaclav Hapla   ierr = PetscMalloc1(nrranks, &rvecs_);CHKERRQ(ierr);
67*f84a5eb8SVaclav Hapla   for (r=0; r<nrranks; r++) {
68*f84a5eb8SVaclav Hapla     /* set array in two steps to mimic PETSC_OWN_POINTER */
69*f84a5eb8SVaclav Hapla     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, 1, rsize[r], NULL, &rvecs_[r]);CHKERRQ(ierr);
70*f84a5eb8SVaclav Hapla     ierr = VecReplaceArray(rvecs_[r], rarr[r]);CHKERRQ(ierr);
71*f84a5eb8SVaclav Hapla   }
72*f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
73*f84a5eb8SVaclav Hapla     ierr = VecRestoreArrayRead(svecs[r], &sarr[r]);CHKERRQ(ierr);
74*f84a5eb8SVaclav Hapla   }
75*f84a5eb8SVaclav Hapla   ierr = PetscFree2(rsize, rarr);CHKERRQ(ierr);
76*f84a5eb8SVaclav Hapla   ierr = PetscFree4(ssize, sarr, rreq, sreq);CHKERRQ(ierr);
77*f84a5eb8SVaclav Hapla   *rvecs = rvecs_;
78*f84a5eb8SVaclav Hapla   PetscFunctionReturn(0);
79*f84a5eb8SVaclav Hapla }
80*f84a5eb8SVaclav Hapla 
81*f84a5eb8SVaclav Hapla static PetscErrorCode SortByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[])
82*f84a5eb8SVaclav Hapla {
83*f84a5eb8SVaclav Hapla   PetscInt            nleaves;
84*f84a5eb8SVaclav Hapla   PetscInt            nranks;
85*f84a5eb8SVaclav Hapla   const PetscMPIInt   *ranks;
86*f84a5eb8SVaclav Hapla   const PetscInt      *roffset, *rmine, *rremote;
87*f84a5eb8SVaclav Hapla   PetscInt            n, o, r;
88*f84a5eb8SVaclav Hapla   PetscErrorCode      ierr;
89*f84a5eb8SVaclav Hapla 
90*f84a5eb8SVaclav Hapla   PetscFunctionBegin;
91*f84a5eb8SVaclav Hapla   ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr);
92*f84a5eb8SVaclav Hapla   nleaves = roffset[nranks];
93*f84a5eb8SVaclav Hapla   ierr = PetscMalloc2(nleaves, rmine1, nleaves, rremote1);CHKERRQ(ierr);
94*f84a5eb8SVaclav Hapla   for (r=0; r<nranks; r++) {
95*f84a5eb8SVaclav Hapla     /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote
96*f84a5eb8SVaclav Hapla        - to unify order with the other side */
97*f84a5eb8SVaclav Hapla     o = roffset[r];
98*f84a5eb8SVaclav Hapla     n = roffset[r+1] - o;
99*f84a5eb8SVaclav Hapla     ierr = PetscMemcpy(&(*rmine1)[o], &rmine[o], n*sizeof(PetscInt));CHKERRQ(ierr);
100*f84a5eb8SVaclav Hapla     ierr = PetscMemcpy(&(*rremote1)[o], &rremote[o], n*sizeof(PetscInt));CHKERRQ(ierr);
101*f84a5eb8SVaclav Hapla     ierr = PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);CHKERRQ(ierr);
102*f84a5eb8SVaclav Hapla   }
103*f84a5eb8SVaclav Hapla   PetscFunctionReturn(0);
104*f84a5eb8SVaclav Hapla }
105*f84a5eb8SVaclav Hapla 
106*f84a5eb8SVaclav Hapla static PetscErrorCode GetRecursiveConeCoordinatesPerRank_Private(DM dm, PetscSF sf, PetscInt rmine[], Vec *coordinatesPerRank[])
107*f84a5eb8SVaclav Hapla {
108*f84a5eb8SVaclav Hapla   IS                  pointsPerRank, conesPerRank;
109*f84a5eb8SVaclav Hapla   PetscInt            nranks;
110*f84a5eb8SVaclav Hapla   const PetscMPIInt   *ranks;
111*f84a5eb8SVaclav Hapla   const PetscInt      *roffset;
112*f84a5eb8SVaclav Hapla   PetscInt            n, o, r;
113*f84a5eb8SVaclav Hapla   PetscErrorCode      ierr;
114*f84a5eb8SVaclav Hapla 
115*f84a5eb8SVaclav Hapla   PetscFunctionBegin;
116*f84a5eb8SVaclav Hapla   ierr = DMGetCoordinatesLocalSetUp(dm);CHKERRQ(ierr);
117*f84a5eb8SVaclav Hapla   ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);CHKERRQ(ierr);
118*f84a5eb8SVaclav Hapla   ierr = PetscMalloc1(nranks, coordinatesPerRank);CHKERRQ(ierr);
119*f84a5eb8SVaclav Hapla   for (r=0; r<nranks; r++) {
120*f84a5eb8SVaclav Hapla     o = roffset[r];
121*f84a5eb8SVaclav Hapla     n = roffset[r+1] - o;
122*f84a5eb8SVaclav Hapla     ierr = ISCreateGeneral(PETSC_COMM_SELF, n, &rmine[o], PETSC_USE_POINTER, &pointsPerRank);CHKERRQ(ierr);
123*f84a5eb8SVaclav Hapla     ierr = DMPlexGetConeRecursive(dm, pointsPerRank, &conesPerRank);CHKERRQ(ierr);
124*f84a5eb8SVaclav Hapla     ierr = DMGetCoordinatesLocalTuple(dm, conesPerRank, NULL, &(*coordinatesPerRank)[r]);CHKERRQ(ierr);
125*f84a5eb8SVaclav Hapla     ierr = ISDestroy(&pointsPerRank);CHKERRQ(ierr);
126*f84a5eb8SVaclav Hapla     ierr = ISDestroy(&conesPerRank);CHKERRQ(ierr);
127*f84a5eb8SVaclav Hapla   }
128*f84a5eb8SVaclav Hapla   PetscFunctionReturn(0);
129*f84a5eb8SVaclav Hapla }
130*f84a5eb8SVaclav Hapla 
131*f84a5eb8SVaclav Hapla static PetscErrorCode PetscSFComputeMultiRootOriginalNumberingByRank_Private(PetscSF sf, PetscSF imsf, PetscInt *irmine1[])
132*f84a5eb8SVaclav Hapla {
133*f84a5eb8SVaclav Hapla   PetscInt            *mRootsOrigNumbering;
134*f84a5eb8SVaclav Hapla   PetscInt            nileaves, niranks;
135*f84a5eb8SVaclav Hapla   const PetscInt      *iroffset, *irmine, *degree;
136*f84a5eb8SVaclav Hapla   PetscInt            i, n, o, r;
137*f84a5eb8SVaclav Hapla   PetscErrorCode      ierr;
138*f84a5eb8SVaclav Hapla 
139*f84a5eb8SVaclav Hapla   PetscFunctionBegin;
140*f84a5eb8SVaclav Hapla   ierr = PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL);CHKERRQ(ierr);
141*f84a5eb8SVaclav Hapla   ierr = PetscSFGetRanks(imsf, &niranks, NULL, &iroffset, &irmine, NULL);CHKERRQ(ierr);
142*f84a5eb8SVaclav Hapla #if defined(PETSC_USE_DEBUG)
143*f84a5eb8SVaclav Hapla   if (PetscUnlikely(nileaves != iroffset[niranks])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nileaves != iroffset[niranks])");
144*f84a5eb8SVaclav Hapla #endif
145*f84a5eb8SVaclav Hapla   ierr = PetscSFComputeDegreeBegin(sf, &degree);CHKERRQ(ierr);
146*f84a5eb8SVaclav Hapla   ierr = PetscSFComputeDegreeEnd(sf, &degree);CHKERRQ(ierr);
147*f84a5eb8SVaclav Hapla   ierr = PetscSFComputeMultiRootOriginalNumbering(sf, degree, &mRootsOrigNumbering);CHKERRQ(ierr);
148*f84a5eb8SVaclav Hapla   ierr = PetscMalloc1(nileaves, irmine1);CHKERRQ(ierr);
149*f84a5eb8SVaclav Hapla   for (r=0; r<niranks; r++) {
150*f84a5eb8SVaclav Hapla     o = iroffset[r];
151*f84a5eb8SVaclav Hapla     n = iroffset[r+1] - o;
152*f84a5eb8SVaclav Hapla     for (i=0; i<n; i++) (*irmine1)[o+i] = mRootsOrigNumbering[irmine[o+i]];
153*f84a5eb8SVaclav Hapla   }
154*f84a5eb8SVaclav Hapla   ierr = PetscFree(mRootsOrigNumbering);CHKERRQ(ierr);
155*f84a5eb8SVaclav Hapla   PetscFunctionReturn(0);
156*f84a5eb8SVaclav Hapla }
157*f84a5eb8SVaclav Hapla 
158*f84a5eb8SVaclav Hapla /* TODO documentation */
159*f84a5eb8SVaclav Hapla PetscErrorCode DMPlexCheckConesConformOnInterfaces(DM dm)
160*f84a5eb8SVaclav Hapla {
161*f84a5eb8SVaclav Hapla   PetscSF             sf;
162*f84a5eb8SVaclav Hapla   PetscInt            nleaves, nranks, nroots;
163*f84a5eb8SVaclav Hapla   const PetscInt      *mine, *roffset, *rmine, *rremote;
164*f84a5eb8SVaclav Hapla   const PetscSFNode   *remote;
165*f84a5eb8SVaclav Hapla   const PetscMPIInt   *ranks;
166*f84a5eb8SVaclav Hapla   PetscSF             msf, imsf;
167*f84a5eb8SVaclav Hapla   PetscInt            nileaves, niranks;
168*f84a5eb8SVaclav Hapla   const PetscMPIInt   *iranks;
169*f84a5eb8SVaclav Hapla   const PetscInt      *iroffset, *irmine, *irremote;
170*f84a5eb8SVaclav Hapla   PetscInt            *rmine1, *rremote1; /* rmine and rremote copies simultaneously sorted by rank and rremote */
171*f84a5eb8SVaclav Hapla   PetscInt            *mine_orig_numbering;
172*f84a5eb8SVaclav Hapla   Vec                 *sntCoordinatesPerRank;
173*f84a5eb8SVaclav Hapla   Vec                 *refCoordinatesPerRank;
174*f84a5eb8SVaclav Hapla   Vec                 *recCoordinatesPerRank;
175*f84a5eb8SVaclav Hapla   PetscInt            r;
176*f84a5eb8SVaclav Hapla   PetscMPIInt         commsize, myrank;
177*f84a5eb8SVaclav Hapla   PetscBool           same;
178*f84a5eb8SVaclav Hapla   MPI_Comm            comm;
179*f84a5eb8SVaclav Hapla   PetscErrorCode      ierr;
180*f84a5eb8SVaclav Hapla 
181*f84a5eb8SVaclav Hapla   PetscFunctionBegin;
182*f84a5eb8SVaclav Hapla   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
183*f84a5eb8SVaclav Hapla   ierr = MPI_Comm_rank(comm, &myrank);CHKERRQ(ierr);
184*f84a5eb8SVaclav Hapla   ierr = MPI_Comm_size(comm, &commsize);CHKERRQ(ierr);
185*f84a5eb8SVaclav Hapla   if (commsize < 2) PetscFunctionReturn(0);
186*f84a5eb8SVaclav Hapla   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
187*f84a5eb8SVaclav Hapla   if (!sf) PetscFunctionReturn(0);
188*f84a5eb8SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, &mine, &remote);CHKERRQ(ierr);
189*f84a5eb8SVaclav Hapla   if (nroots < 0) PetscFunctionReturn(0);
190*f84a5eb8SVaclav Hapla   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
191*f84a5eb8SVaclav Hapla   ierr = PetscSFGetRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);CHKERRQ(ierr);
192*f84a5eb8SVaclav Hapla 
193*f84a5eb8SVaclav Hapla   /* Expand sent cones per rank */
194*f84a5eb8SVaclav Hapla   ierr = SortByRemote_Private(sf, &rmine1, &rremote1);CHKERRQ(ierr);
195*f84a5eb8SVaclav Hapla   ierr = GetRecursiveConeCoordinatesPerRank_Private(dm, sf, rmine1, &sntCoordinatesPerRank);CHKERRQ(ierr);
196*f84a5eb8SVaclav Hapla 
197*f84a5eb8SVaclav Hapla   /* Create inverse SF */
198*f84a5eb8SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
199*f84a5eb8SVaclav Hapla   ierr = PetscSFCreateInverseSF(msf,&imsf);CHKERRQ(ierr);
200*f84a5eb8SVaclav Hapla   ierr = PetscSFSetUp(imsf);CHKERRQ(ierr);
201*f84a5eb8SVaclav Hapla   ierr = PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL);CHKERRQ(ierr);
202*f84a5eb8SVaclav Hapla   ierr = PetscSFGetRanks(imsf, &niranks, &iranks, &iroffset, &irmine, &irremote);CHKERRQ(ierr);
203*f84a5eb8SVaclav Hapla 
204*f84a5eb8SVaclav Hapla   /* Compute original numbering of multi-roots (referenced points) */
205*f84a5eb8SVaclav Hapla   ierr = PetscSFComputeMultiRootOriginalNumberingByRank_Private(sf, imsf, &mine_orig_numbering);CHKERRQ(ierr);
206*f84a5eb8SVaclav Hapla 
207*f84a5eb8SVaclav Hapla   /* Expand referred cones per rank */
208*f84a5eb8SVaclav Hapla   ierr = GetRecursiveConeCoordinatesPerRank_Private(dm, imsf, mine_orig_numbering, &refCoordinatesPerRank);CHKERRQ(ierr);
209*f84a5eb8SVaclav Hapla 
210*f84a5eb8SVaclav Hapla   /* Send the coordinates */
211*f84a5eb8SVaclav Hapla   ierr = ExchangeVecByRank_Private((PetscObject)sf, nranks, ranks, sntCoordinatesPerRank, niranks, iranks, &recCoordinatesPerRank);CHKERRQ(ierr);
212*f84a5eb8SVaclav Hapla 
213*f84a5eb8SVaclav Hapla   /* Compare recCoordinatesPerRank with refCoordinatesPerRank */
214*f84a5eb8SVaclav Hapla   for (r=0; r<niranks; r++) {
215*f84a5eb8SVaclav Hapla     ierr = VecEqual(refCoordinatesPerRank[r], recCoordinatesPerRank[r], &same);CHKERRQ(ierr);
216*f84a5eb8SVaclav Hapla     if (!same) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "interface cones do not conform for remote rank %d", iranks[r]);
217*f84a5eb8SVaclav Hapla   }
218*f84a5eb8SVaclav Hapla 
219*f84a5eb8SVaclav Hapla   /* destroy sent stuff */
220*f84a5eb8SVaclav Hapla   for (r=0; r<nranks; r++) {
221*f84a5eb8SVaclav Hapla     ierr = VecDestroy(&sntCoordinatesPerRank[r]);CHKERRQ(ierr);
222*f84a5eb8SVaclav Hapla   }
223*f84a5eb8SVaclav Hapla   ierr = PetscFree(sntCoordinatesPerRank);CHKERRQ(ierr);
224*f84a5eb8SVaclav Hapla   ierr = PetscFree2(rmine1, rremote1);CHKERRQ(ierr);
225*f84a5eb8SVaclav Hapla   ierr = PetscSFDestroy(&imsf);CHKERRQ(ierr);
226*f84a5eb8SVaclav Hapla 
227*f84a5eb8SVaclav Hapla   /* destroy referenced stuff */
228*f84a5eb8SVaclav Hapla   for (r=0; r<niranks; r++) {
229*f84a5eb8SVaclav Hapla     ierr = VecDestroy(&refCoordinatesPerRank[r]);CHKERRQ(ierr);
230*f84a5eb8SVaclav Hapla   }
231*f84a5eb8SVaclav Hapla   ierr = PetscFree(refCoordinatesPerRank);CHKERRQ(ierr);
232*f84a5eb8SVaclav Hapla   ierr = PetscFree(mine_orig_numbering);CHKERRQ(ierr);
233*f84a5eb8SVaclav Hapla 
234*f84a5eb8SVaclav Hapla   /* destroy received stuff */
235*f84a5eb8SVaclav Hapla   for (r=0; r<niranks; r++) {
236*f84a5eb8SVaclav Hapla     ierr = VecDestroy(&recCoordinatesPerRank[r]);CHKERRQ(ierr);
237*f84a5eb8SVaclav Hapla   }
238*f84a5eb8SVaclav Hapla   ierr = PetscFree(recCoordinatesPerRank);CHKERRQ(ierr);
239*f84a5eb8SVaclav Hapla   PetscFunctionReturn(0);
240*f84a5eb8SVaclav Hapla }
241