xref: /petsc/src/dm/impls/plex/plexcheckinterface.c (revision 6858538e710279fe46cd8279ab34c98b10293bbd)
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;
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));
22f84a5eb8SVaclav Hapla   for (r=0; r<nrranks; r++) {
239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(&rsize[r], 1, MPIU_INT, rranks[r], tag, comm, &rreq[r]));
24f84a5eb8SVaclav Hapla   }
25f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
269566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(&ssize[r], 1, MPIU_INT, sranks[r], tag, comm, &sreq[r]));
27f84a5eb8SVaclav Hapla   }
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE));
299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE));
30f84a5eb8SVaclav Hapla   /* exchange array */
319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag(obj,&tag));
32f84a5eb8SVaclav Hapla   for (r=0; r<nrranks; r++) {
339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(rsize[r]*unitsize, &rarr[r]));
349566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(rarr[r], rsize[r], dt, rranks[r], tag, comm, &rreq[r]));
35f84a5eb8SVaclav Hapla   }
36f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
379566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(sarr[r], ssize[r], dt, sranks[r], tag, comm, &sreq[r]));
38f84a5eb8SVaclav Hapla   }
399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nrranks, rreq, MPI_STATUSES_IGNORE));
409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(nsranks, sreq, MPI_STATUSES_IGNORE));
419566063dSJacob Faibussowitsch   PetscCall(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;
599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(nsranks, &ssize, nsranks, &sarr, nrranks, &rreq, nsranks, &sreq));
60f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
619566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(svecs[r], &ssize[r]));
629566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(svecs[r], &sarr[r]));
63f84a5eb8SVaclav Hapla   }
649566063dSJacob Faibussowitsch   PetscCall(ExchangeArrayByRank_Private(obj, MPIU_SCALAR, nsranks, sranks, ssize, (const void**)sarr, nrranks, rranks, &rsize, (void***)&rarr));
659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrranks, &rvecs_));
66f84a5eb8SVaclav Hapla   for (r=0; r<nrranks; r++) {
67f84a5eb8SVaclav Hapla     /* set array in two steps to mimic PETSC_OWN_POINTER */
689566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, rsize[r], NULL, &rvecs_[r]));
699566063dSJacob Faibussowitsch     PetscCall(VecReplaceArray(rvecs_[r], rarr[r]));
70f84a5eb8SVaclav Hapla   }
71f84a5eb8SVaclav Hapla   for (r=0; r<nsranks; r++) {
729566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(svecs[r], &sarr[r]));
73f84a5eb8SVaclav Hapla   }
749566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rsize, rarr));
759566063dSJacob Faibussowitsch   PetscCall(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;
899566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote));
90f84a5eb8SVaclav Hapla   nleaves = roffset[nranks];
919566063dSJacob Faibussowitsch   PetscCall(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;
979566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&(*rmine1)[o], &rmine[o], n));
989566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&(*rremote1)[o], &rremote[o], n));
999566063dSJacob Faibussowitsch     PetscCall(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;
1139566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(dm));
1149566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, NULL, NULL));
1159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nranks, coordinatesPerRank));
116f84a5eb8SVaclav Hapla   for (r=0; r<nranks; r++) {
117f84a5eb8SVaclav Hapla     o = roffset[r];
118f84a5eb8SVaclav Hapla     n = roffset[r+1] - o;
1199566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, &rmine[o], PETSC_USE_POINTER, &pointsPerRank));
1209566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeRecursiveVertices(dm, pointsPerRank, &conesPerRank));
1219566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocalTuple(dm, conesPerRank, NULL, &(*coordinatesPerRank)[r]));
1229566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&pointsPerRank));
1239566063dSJacob Faibussowitsch     PetscCall(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;
1369566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL));
1379566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks(imsf, &niranks, NULL, &iroffset, &irmine, NULL));
13808401ef6SPierre Jolivet   PetscCheck(nileaves == iroffset[niranks],PETSC_COMM_SELF,PETSC_ERR_PLIB,"nileaves != iroffset[niranks])");
1399566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
1409566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
1419566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, degree, NULL, &mRootsOrigNumbering));
1429566063dSJacob Faibussowitsch   PetscCall(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   }
1489566063dSJacob Faibussowitsch   PetscCall(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 
171db781477SPatrick Sanan .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);
1979566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &myrank));
1999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &commsize));
200f84a5eb8SVaclav Hapla   if (commsize < 2) PetscFunctionReturn(0);
2019566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
202f84a5eb8SVaclav Hapla   if (!sf) PetscFunctionReturn(0);
2039566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &mine, &remote));
204f84a5eb8SVaclav Hapla   if (nroots < 0) PetscFunctionReturn(0);
205*6858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[0].x || dm->coordinates[0].xl, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM coordinates must be set");
2069566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
2079566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote));
208f84a5eb8SVaclav Hapla 
209f84a5eb8SVaclav Hapla   /* Expand sent cones per rank */
2109566063dSJacob Faibussowitsch   PetscCall(SortByRemote_Private(sf, &rmine1, &rremote1));
2119566063dSJacob Faibussowitsch   PetscCall(GetRecursiveConeCoordinatesPerRank_Private(dm, sf, rmine1, &sntCoordinatesPerRank));
212f84a5eb8SVaclav Hapla 
213f84a5eb8SVaclav Hapla   /* Create inverse SF */
2149566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf,&msf));
2159566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateInverseSF(msf,&imsf));
2169566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(imsf));
2179566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(imsf, NULL, &nileaves, NULL, NULL));
2189566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks(imsf, &niranks, &iranks, &iroffset, &irmine, &irremote));
219f84a5eb8SVaclav Hapla 
220f84a5eb8SVaclav Hapla   /* Compute original numbering of multi-roots (referenced points) */
2219566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeMultiRootOriginalNumberingByRank_Private(sf, imsf, &mine_orig_numbering));
222f84a5eb8SVaclav Hapla 
223124f9872SVaclav Hapla   /* Expand coordinates of the referred cones per rank */
2249566063dSJacob Faibussowitsch   PetscCall(GetRecursiveConeCoordinatesPerRank_Private(dm, imsf, mine_orig_numbering, &refCoordinatesPerRank));
225f84a5eb8SVaclav Hapla 
226f84a5eb8SVaclav Hapla   /* Send the coordinates */
2279566063dSJacob Faibussowitsch   PetscCall(ExchangeVecByRank_Private((PetscObject)sf, nranks, ranks, sntCoordinatesPerRank, niranks, iranks, &recCoordinatesPerRank));
228f84a5eb8SVaclav Hapla 
2298f2c89e7SVaclav Hapla   /* verbose output */
2309566063dSJacob Faibussowitsch   PetscCall(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;
2339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(v, "============\nDMPlexCheckInterfaceCones output\n============\n"));
2349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(v));
2359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", myrank));
2368f2c89e7SVaclav Hapla     for (r=0; r<nranks; r++) {
23763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(v, "  r=%" PetscInt_FMT " ranks[r]=%d sntCoordinatesPerRank[r]:\n", r, ranks[r]));
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(v,PETSC_COMM_SELF,&sv));
2409566063dSJacob Faibussowitsch       PetscCall(VecView(sntCoordinatesPerRank[r], sv));
2419566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(v,PETSC_COMM_SELF,&sv));
2429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
2438f2c89e7SVaclav Hapla     }
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "  ----------\n"));
2458f2c89e7SVaclav Hapla     for (r=0; r<niranks; r++) {
24663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(v, "  r=%" PetscInt_FMT " iranks[r]=%d refCoordinatesPerRank[r]:\n", r, iranks[r]));
2479566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
2489566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(v,PETSC_COMM_SELF,&sv));
2499566063dSJacob Faibussowitsch       PetscCall(VecView(refCoordinatesPerRank[r], sv));
2509566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(v,PETSC_COMM_SELF,&sv));
2519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
2528f2c89e7SVaclav Hapla     }
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "  ----------\n"));
2548f2c89e7SVaclav Hapla     for (r=0; r<niranks; r++) {
25563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(v, "  r=%" PetscInt_FMT " iranks[r]=%d recCoordinatesPerRank[r]:\n", r, iranks[r]));
2569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(v));
2579566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(v,PETSC_COMM_SELF,&sv));
2589566063dSJacob Faibussowitsch       PetscCall(VecView(recCoordinatesPerRank[r], sv));
2599566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(v,PETSC_COMM_SELF,&sv));
2609566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(v));
2618f2c89e7SVaclav Hapla     }
2629566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(v));
2639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(v));
2648f2c89e7SVaclav Hapla   }
2658f2c89e7SVaclav Hapla 
266f84a5eb8SVaclav Hapla   /* Compare recCoordinatesPerRank with refCoordinatesPerRank */
267f84a5eb8SVaclav Hapla   for (r=0; r<niranks; r++) {
2689566063dSJacob Faibussowitsch     PetscCall(VecEqual(refCoordinatesPerRank[r], recCoordinatesPerRank[r], &same));
26928b400f6SJacob Faibussowitsch     PetscCheck(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++) {
2749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&sntCoordinatesPerRank[r]));
275f84a5eb8SVaclav Hapla   }
2769566063dSJacob Faibussowitsch   PetscCall(PetscFree(sntCoordinatesPerRank));
2779566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rmine1, rremote1));
2789566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&imsf));
279f84a5eb8SVaclav Hapla 
280f84a5eb8SVaclav Hapla   /* destroy referenced stuff */
281f84a5eb8SVaclav Hapla   for (r=0; r<niranks; r++) {
2829566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&refCoordinatesPerRank[r]));
283f84a5eb8SVaclav Hapla   }
2849566063dSJacob Faibussowitsch   PetscCall(PetscFree(refCoordinatesPerRank));
2859566063dSJacob Faibussowitsch   PetscCall(PetscFree(mine_orig_numbering));
286f84a5eb8SVaclav Hapla 
287f84a5eb8SVaclav Hapla   /* destroy received stuff */
288f84a5eb8SVaclav Hapla   for (r=0; r<niranks; r++) {
2899566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&recCoordinatesPerRank[r]));
290f84a5eb8SVaclav Hapla   }
2919566063dSJacob Faibussowitsch   PetscCall(PetscFree(recCoordinatesPerRank));
292f84a5eb8SVaclav Hapla   PetscFunctionReturn(0);
293f84a5eb8SVaclav Hapla }
294