xref: /petsc/src/dm/impls/plex/adaptors/parmmg/parmmgadapt.c (revision d9126033840fa4b9e125cadad06c85976e6bf099)
15f80ce2aSJacob Faibussowitsch #include "../mmgcommon.h" /*I      "petscdmplex.h"   I*/
22bc0b47dSJoe Wallwork #include <parmmg/libparmmg.h>
32bc0b47dSJoe Wallwork 
4b2c4de1dSPierre Jolivet PetscBool  ParMmgCite       = PETSC_FALSE;
5b2c4de1dSPierre Jolivet const char ParMmgCitation[] = "@techreport{cirrottola:hal-02386837,\n"
6b2c4de1dSPierre Jolivet                               "  title       = {Parallel unstructured mesh adaptation using iterative remeshing and repartitioning},\n"
7b2c4de1dSPierre Jolivet                               "  institution = {Inria Bordeaux},\n"
8b2c4de1dSPierre Jolivet                               "  author      = {L. Cirrottola and A. Froehly},\n"
9b2c4de1dSPierre Jolivet                               "  number      = {9307},\n"
10b2c4de1dSPierre Jolivet                               "  note        = {\\url{https://hal.inria.fr/hal-02386837}},\n"
11b2c4de1dSPierre Jolivet                               "  year        = {2019}\n}\n";
12b2c4de1dSPierre Jolivet 
138a7f912aSStephan Kramer /*
148a7f912aSStephan Kramer  Coupling code for the ParMmg metric-based mesh adaptation package.
158a7f912aSStephan Kramer */
DMAdaptMetric_ParMmg_Plex(DM dm,Vec vertexMetric,DMLabel bdLabel,DMLabel rgLabel,DM * dmNew)16d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew)
17d71ae5a4SJacob Faibussowitsch {
180a7a67b9SPierre Jolivet   MPI_Comm           comm;
192bc0b47dSJoe Wallwork   const char        *bdName = "_boundary_";
209fe9e680SJoe Wallwork   const char        *rgName = "_regions_";
212bc0b47dSJoe Wallwork   DM                 udm, cdm;
22c1dc6da0SJoe Wallwork   DMLabel            bdLabelNew, rgLabelNew;
239fe9e680SJoe Wallwork   const char        *bdLabelName, *rgLabelName;
24c1dc6da0SJoe Wallwork   IS                 globalVertexNum;
252bc0b47dSJoe Wallwork   PetscSection       coordSection;
262bc0b47dSJoe Wallwork   Vec                coordinates;
272bc0b47dSJoe Wallwork   PetscSF            sf;
282bc0b47dSJoe Wallwork   const PetscScalar *coords, *met;
29ae8b063eSJoe Wallwork   PetscReal         *vertices, *metric, *verticesNew, *verticesNewLoc, gradationFactor, hausdorffNumber;
30c1dc6da0SJoe Wallwork   PetscInt          *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew;
31c1dc6da0SJoe Wallwork   PetscInt          *bdFaces, *faceTags, *facesNew, *faceTagsNew;
32c1dc6da0SJoe Wallwork   PetscInt          *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces;
33c1dc6da0SJoe Wallwork   PetscInt           cStart, cEnd, c, numCells, fStart, fEnd, f, numFaceTags, vStart, vEnd, v, numVertices;
348a7f912aSStephan Kramer   PetscInt           numCellsNotShared, *cIsLeaf, numUsedVertices, *vertexNumber, *fIsIncluded;
35c1dc6da0SJoe Wallwork   PetscInt           dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, numIter;
368a7f912aSStephan Kramer   PetscInt          *interfaces_lv, *interfaces_gv, *interfacesOffset;
376497c311SBarry Smith   PetscMPIInt        niranks, nrranks, numNgbRanks;
386497c311SBarry Smith   PetscInt           r, lv, gv;
39c1dc6da0SJoe Wallwork   PetscInt          *gv_new, *owners, *verticesNewSorted, pStart, pEnd;
402bc0b47dSJoe Wallwork   PetscInt           numCellsNew, numVerticesNew, numCornersNew, numFacesNew, numVerticesNewLoc;
41c1dc6da0SJoe Wallwork   const PetscInt    *gV, *ioffset, *irootloc, *roffset, *rmine, *rremote;
4276f360caSJoe Wallwork   PetscBool          flg = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform;
432bc0b47dSJoe Wallwork   const PetscMPIInt *iranks, *rranks;
442bc0b47dSJoe Wallwork   PetscMPIInt        numProcs, rank;
452bc0b47dSJoe Wallwork   PMMG_pParMesh      parmesh = NULL;
462bc0b47dSJoe Wallwork 
478a7f912aSStephan Kramer   // DEVELOPER NOTE: ParMmg wants to know the rank of every process which is sharing a given point and
488a7f912aSStephan Kramer   //                 for this information to be conveyed to every process that is sharing that point.
498a7f912aSStephan Kramer 
502bc0b47dSJoe Wallwork   PetscFunctionBegin;
51b2c4de1dSPierre Jolivet   PetscCall(PetscCitationsRegister(ParMmgCitation, &ParMmgCite));
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &numProcs));
549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
552bc0b47dSJoe Wallwork   if (bdLabel) {
569566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)bdLabel, &bdLabelName));
579566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(bdLabelName, bdName, &flg));
5828b400f6SJacob Faibussowitsch     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
592bc0b47dSJoe Wallwork   }
609fe9e680SJoe Wallwork   if (rgLabel) {
619566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)rgLabel, &rgLabelName));
629566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(rgLabelName, rgName, &flg));
6328b400f6SJacob Faibussowitsch     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
649fe9e680SJoe Wallwork   }
652bc0b47dSJoe Wallwork 
662bc0b47dSJoe Wallwork   /* Get mesh information */
679566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
685f80ce2aSJacob Faibussowitsch   PetscCheck(dim == 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "ParMmg only works in 3D.");
692bc0b47dSJoe Wallwork   Neq = (dim * (dim + 1)) / 2;
709566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
719566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
729566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
739566063dSJacob Faibussowitsch   PetscCall(DMPlexUninterpolate(dm, &udm));
749566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
752bc0b47dSJoe Wallwork   numCells    = cEnd - cStart;
762bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
772bc0b47dSJoe Wallwork 
788a7f912aSStephan Kramer   /* Get parallel data; work out which cells are owned and which are leaves */
798a7f912aSStephan Kramer   PetscCall(PetscCalloc1(numCells, &cIsLeaf));
808a7f912aSStephan Kramer   numCellsNotShared = numCells;
818a7f912aSStephan Kramer   niranks = nrranks = 0;
828a7f912aSStephan Kramer   if (numProcs > 1) {
838a7f912aSStephan Kramer     PetscCall(DMGetPointSF(dm, &sf));
848a7f912aSStephan Kramer     PetscCall(PetscSFSetUp(sf));
858a7f912aSStephan Kramer     PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
868a7f912aSStephan Kramer     PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
878a7f912aSStephan Kramer     for (r = 0; r < nrranks; ++r) {
888a7f912aSStephan Kramer       for (i = roffset[r]; i < roffset[r + 1]; ++i) {
898a7f912aSStephan Kramer         if (rmine[i] >= cStart && rmine[i] < cEnd) {
908a7f912aSStephan Kramer           cIsLeaf[rmine[i] - cStart] = 1;
918a7f912aSStephan Kramer           numCellsNotShared--;
928a7f912aSStephan Kramer         }
938a7f912aSStephan Kramer       }
948a7f912aSStephan Kramer     }
958a7f912aSStephan Kramer   }
968a7f912aSStephan Kramer 
978a7f912aSStephan Kramer   /*
988a7f912aSStephan Kramer     Create a vertex numbering for ParMmg starting at 1. Vertices not included in any
998a7f912aSStephan Kramer     owned cell remain 0 and will be removed. Using this numbering, create cells.
1008a7f912aSStephan Kramer   */
1018a7f912aSStephan Kramer   numUsedVertices = 0;
1028a7f912aSStephan Kramer   PetscCall(PetscCalloc1(numVertices, &vertexNumber));
1038a7f912aSStephan Kramer   PetscCall(PetscMalloc1(numCellsNotShared * maxConeSize, &cells));
1042bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
1052bc0b47dSJoe Wallwork     const PetscInt *cone;
1062bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
1072bc0b47dSJoe Wallwork 
1088a7f912aSStephan Kramer     if (!cIsLeaf[c]) {
1098a7f912aSStephan Kramer       PetscCall(DMPlexGetConeSize(udm, cStart + c, &coneSize));
1108a7f912aSStephan Kramer       PetscCall(DMPlexGetCone(udm, cStart + c, &cone));
1118a7f912aSStephan Kramer       for (cl = 0; cl < coneSize; ++cl) {
1128a7f912aSStephan Kramer         if (!vertexNumber[cone[cl] - vStart]) vertexNumber[cone[cl] - vStart] = ++numUsedVertices;
1138a7f912aSStephan Kramer         cells[coff++] = vertexNumber[cone[cl] - vStart];
1148a7f912aSStephan Kramer       }
1158a7f912aSStephan Kramer     }
1162bc0b47dSJoe Wallwork   }
1172bc0b47dSJoe Wallwork 
1188a7f912aSStephan Kramer   /* Get array of vertex coordinates */
1199566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
1209566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(cdm, &coordSection));
1219566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
1238a7f912aSStephan Kramer   PetscCall(PetscMalloc2(numUsedVertices * Neq, &metric, dim * numUsedVertices, &vertices));
1242bc0b47dSJoe Wallwork   for (v = 0; v < vEnd - vStart; ++v) {
1259566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(coordSection, v + vStart, &off));
1268a7f912aSStephan Kramer     if (vertexNumber[v]) {
1278a7f912aSStephan Kramer       for (i = 0; i < dim; ++i) vertices[dim * (vertexNumber[v] - 1) + i] = PetscRealPart(coords[off + i]);
1288a7f912aSStephan Kramer     }
1292bc0b47dSJoe Wallwork   }
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
1312bc0b47dSJoe Wallwork 
132c1dc6da0SJoe Wallwork   /* Get face tags */
133c1dc6da0SJoe Wallwork   if (!bdLabel) {
134c1dc6da0SJoe Wallwork     flg = PETSC_TRUE;
1359566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel));
1369566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
137c1dc6da0SJoe Wallwork   }
1389566063dSJacob Faibussowitsch   PetscCall(DMLabelGetBounds(bdLabel, &pStart, &pEnd));
1398a7f912aSStephan Kramer   PetscCall(PetscCalloc1(pEnd - pStart, &fIsIncluded));
140c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
141c1dc6da0SJoe Wallwork     PetscBool hasPoint;
142c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1432bc0b47dSJoe Wallwork 
1449566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
145c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
1468a7f912aSStephan Kramer 
1478a7f912aSStephan Kramer     /* Only faces adjacent to an owned (non-leaf) cell are included */
1488a7f912aSStephan Kramer     PetscInt        nnbrs;
1498a7f912aSStephan Kramer     const PetscInt *nbrs;
1508a7f912aSStephan Kramer     PetscCall(DMPlexGetSupportSize(dm, f, &nnbrs));
1518a7f912aSStephan Kramer     PetscCall(DMPlexGetSupport(dm, f, &nbrs));
1528a7f912aSStephan Kramer     for (c = 0; c < nnbrs; ++c) fIsIncluded[f - pStart] = fIsIncluded[f - pStart] || !cIsLeaf[nbrs[c]];
1538a7f912aSStephan Kramer     if (!fIsIncluded[f - pStart]) continue;
1548a7f912aSStephan Kramer 
155c1dc6da0SJoe Wallwork     numFaceTags++;
156c1dc6da0SJoe Wallwork 
1579566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1582bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize * 2; cl += 2) {
1592bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1602bc0b47dSJoe Wallwork     }
1619566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1622bc0b47dSJoe Wallwork   }
1639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags));
164c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
165c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1662bc0b47dSJoe Wallwork 
1678a7f912aSStephan Kramer     if (!fIsIncluded[f - pStart]) continue;
168c1dc6da0SJoe Wallwork 
1699566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1702bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize * 2; cl += 2) {
1718a7f912aSStephan Kramer       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = vertexNumber[closure[cl] - vStart];
1722bc0b47dSJoe Wallwork     }
1739566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1749566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]));
1752bc0b47dSJoe Wallwork   }
1768a7f912aSStephan Kramer   PetscCall(PetscFree(fIsIncluded));
1772bc0b47dSJoe Wallwork 
1789fe9e680SJoe Wallwork   /* Get cell tags */
1798a7f912aSStephan Kramer   PetscCall(PetscCalloc2(numUsedVertices, &verTags, numCellsNotShared, &cellTags));
1809fe9e680SJoe Wallwork   if (rgLabel) {
1818a7f912aSStephan Kramer     for (c = cStart, coff = 0; c < cEnd; ++c) {
1823a7d0413SPierre Jolivet       if (!cIsLeaf[c - cStart]) PetscCall(DMLabelGetValue(rgLabel, c, &cellTags[coff++]));
1839fe9e680SJoe Wallwork     }
1848a7f912aSStephan Kramer   }
1858a7f912aSStephan Kramer   PetscCall(PetscFree(cIsLeaf));
1869fe9e680SJoe Wallwork 
1878a7f912aSStephan Kramer   /* Get metric, using only the upper triangular part */
1889566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
1899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(vertexMetric, &met));
1909566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1919566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1922bc0b47dSJoe Wallwork   for (v = 0; v < (vEnd - vStart); ++v) {
1938a7f912aSStephan Kramer     PetscInt vv = vertexNumber[v];
1948a7f912aSStephan Kramer     if (!vv--) continue;
1952bc0b47dSJoe Wallwork     for (i = 0, k = 0; i < dim; ++i) {
1960a7a67b9SPierre Jolivet       for (j = i; j < dim; ++j, ++k) {
19793520041SJoe Wallwork         if (isotropic) {
19893520041SJoe Wallwork           if (i == j) {
1998a7f912aSStephan Kramer             if (uniform) metric[Neq * vv + k] = PetscRealPart(met[0]);
2008a7f912aSStephan Kramer             else metric[Neq * vv + k] = PetscRealPart(met[v]);
2018a7f912aSStephan Kramer           } else metric[Neq * vv + k] = 0.0;
2028a7f912aSStephan Kramer         } else metric[Neq * vv + k] = PetscRealPart(met[dim * dim * v + dim * i + j]);
2032bc0b47dSJoe Wallwork       }
2042bc0b47dSJoe Wallwork     }
2052bc0b47dSJoe Wallwork   }
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(vertexMetric, &met));
2072bc0b47dSJoe Wallwork 
2088a7f912aSStephan Kramer   /* Build ParMmg communicators: the list of vertices between two partitions  */
2092bc0b47dSJoe Wallwork   numNgbRanks = 0;
2102bc0b47dSJoe Wallwork   if (numProcs > 1) {
2118a7f912aSStephan Kramer     DM              rankdm;
2128a7f912aSStephan Kramer     PetscSection    rankSection, rankGlobalSection;
2138a7f912aSStephan Kramer     PetscSF         rankSF;
2148a7f912aSStephan Kramer     const PetscInt *degree;
2158a7f912aSStephan Kramer     PetscInt       *rankOfUsedVertices, *rankOfUsedMultiRootLeaves, *usedCopies;
2168a7f912aSStephan Kramer     PetscInt       *rankArray, *rankGlobalArray, *interfacesPerRank;
2178a7f912aSStephan Kramer     PetscInt        offset, mrl, rootDegreeCnt, s, shareCnt, gv;
2182bc0b47dSJoe Wallwork 
2198a7f912aSStephan Kramer     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
2208a7f912aSStephan Kramer     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
2218a7f912aSStephan Kramer     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
2228a7f912aSStephan Kramer     for (i = 0, rootDegreeCnt = 0; i < pEnd - pStart; ++i) rootDegreeCnt += degree[i];
2238a7f912aSStephan Kramer 
2248a7f912aSStephan Kramer     /* rankOfUsedVertices, point-array: rank+1 if vertex and in use */
2258a7f912aSStephan Kramer     PetscCall(PetscCalloc1(pEnd - pStart, &rankOfUsedVertices));
2268a7f912aSStephan Kramer     for (i = 0; i < pEnd - pStart; ++i) rankOfUsedVertices[i] = -1;
2278a7f912aSStephan Kramer     for (i = vStart; i < vEnd; ++i) {
2288a7f912aSStephan Kramer       if (vertexNumber[i - vStart]) rankOfUsedVertices[i] = rank;
2292bc0b47dSJoe Wallwork     }
230c51fff1bSJoe Wallwork 
2318a7f912aSStephan Kramer     /* rankOfUsedMultiRootLeaves, multiroot-array: rank if vertex and in use, else -1 */
2328a7f912aSStephan Kramer     PetscCall(PetscMalloc1(rootDegreeCnt, &rankOfUsedMultiRootLeaves));
2338a7f912aSStephan Kramer     PetscCall(PetscSFGatherBegin(sf, MPIU_INT, rankOfUsedVertices, rankOfUsedMultiRootLeaves));
2348a7f912aSStephan Kramer     PetscCall(PetscSFGatherEnd(sf, MPIU_INT, rankOfUsedVertices, rankOfUsedMultiRootLeaves));
2358a7f912aSStephan Kramer     PetscCall(PetscFree(rankOfUsedVertices));
2368a7f912aSStephan Kramer 
2378a7f912aSStephan Kramer     /* usedCopies, point-array: if vertex, shared by how many processes */
2388a7f912aSStephan Kramer     PetscCall(PetscCalloc1(pEnd - pStart, &usedCopies));
2398a7f912aSStephan Kramer     for (i = 0, mrl = 0; i < vStart - pStart; i++) mrl += degree[i];
2408a7f912aSStephan Kramer     for (i = vStart - pStart; i < vEnd - pStart; ++i) {
2418a7f912aSStephan Kramer       for (j = 0; j < degree[i]; j++, mrl++) {
2428a7f912aSStephan Kramer         if (rankOfUsedMultiRootLeaves[mrl] != -1) usedCopies[i]++;
2432bc0b47dSJoe Wallwork       }
2448a7f912aSStephan Kramer       if (vertexNumber[i - vStart + pStart]) usedCopies[i]++;
2458a7f912aSStephan Kramer     }
2468a7f912aSStephan Kramer     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, usedCopies, usedCopies, MPI_REPLACE));
2478a7f912aSStephan Kramer     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, usedCopies, usedCopies, MPI_REPLACE));
2488a7f912aSStephan Kramer 
2498a7f912aSStephan Kramer     /* Create a section to store ranks of vertices shared by more than one process */
2508a7f912aSStephan Kramer     PetscCall(PetscSectionCreate(comm, &rankSection));
2518a7f912aSStephan Kramer     PetscCall(PetscSectionSetNumFields(rankSection, 1));
2528a7f912aSStephan Kramer     PetscCall(PetscSectionSetChart(rankSection, pStart, pEnd));
2538a7f912aSStephan Kramer     for (i = vStart - pStart; i < vEnd - pStart; ++i) {
2543a7d0413SPierre Jolivet       if (usedCopies[i] > 1) PetscCall(PetscSectionSetDof(rankSection, i + pStart, usedCopies[i]));
2558a7f912aSStephan Kramer     }
2568a7f912aSStephan Kramer     PetscCall(PetscSectionSetUp(rankSection));
2578a7f912aSStephan Kramer     PetscCall(PetscSectionCreateGlobalSection(rankSection, sf, PETSC_FALSE, PETSC_FALSE, PETSC_TRUE, &rankGlobalSection));
2588a7f912aSStephan Kramer 
2598a7f912aSStephan Kramer     PetscCall(PetscSectionGetStorageSize(rankGlobalSection, &s));
2608a7f912aSStephan Kramer     PetscCall(PetscMalloc1(s, &rankGlobalArray));
2618a7f912aSStephan Kramer     for (i = 0, mrl = 0; i < vStart - pStart; i++) mrl += degree[i];
2628a7f912aSStephan Kramer     for (i = vStart - pStart, k = 0; i < vEnd - pStart; ++i) {
2638a7f912aSStephan Kramer       if (usedCopies[i] > 1 && degree[i]) {
2648a7f912aSStephan Kramer         PetscCall(PetscSectionGetOffset(rankSection, k, &offset));
2658a7f912aSStephan Kramer         if (vertexNumber[i - vStart + pStart]) rankGlobalArray[k++] = rank;
2668a7f912aSStephan Kramer         for (j = 0; j < degree[i]; j++, mrl++) {
267ac530a7eSPierre Jolivet           if (rankOfUsedMultiRootLeaves[mrl] != -1) rankGlobalArray[k++] = rankOfUsedMultiRootLeaves[mrl];
2688a7f912aSStephan Kramer         }
2698a7f912aSStephan Kramer       } else mrl += degree[i];
2708a7f912aSStephan Kramer     }
2718a7f912aSStephan Kramer     PetscCall(PetscFree(rankOfUsedMultiRootLeaves));
2728a7f912aSStephan Kramer     PetscCall(PetscFree(usedCopies));
2738a7f912aSStephan Kramer     PetscCall(PetscSectionDestroy(&rankGlobalSection));
2748a7f912aSStephan Kramer 
2758a7f912aSStephan Kramer     /*
2768a7f912aSStephan Kramer       Broadcast the array of ranks.
2778a7f912aSStephan Kramer         (We want all processes to know all the ranks that are looking at each point.
2788a7f912aSStephan Kramer         Above, we tell the roots. Here, the roots tell the leaves.)
2798a7f912aSStephan Kramer     */
2808a7f912aSStephan Kramer     PetscCall(DMClone(dm, &rankdm));
2818a7f912aSStephan Kramer     PetscCall(DMSetLocalSection(rankdm, rankSection));
2828a7f912aSStephan Kramer     PetscCall(DMGetSectionSF(rankdm, &rankSF));
2838a7f912aSStephan Kramer     PetscCall(PetscSectionGetStorageSize(rankSection, &s));
2848a7f912aSStephan Kramer     PetscCall(PetscMalloc1(s, &rankArray));
2858a7f912aSStephan Kramer     PetscCall(PetscSFBcastBegin(rankSF, MPI_INT, rankGlobalArray, rankArray, MPI_REPLACE));
2868a7f912aSStephan Kramer     PetscCall(PetscSFBcastEnd(rankSF, MPI_INT, rankGlobalArray, rankArray, MPI_REPLACE));
2878a7f912aSStephan Kramer     PetscCall(PetscFree(rankGlobalArray));
2888a7f912aSStephan Kramer     PetscCall(DMDestroy(&rankdm));
2898a7f912aSStephan Kramer 
2908a7f912aSStephan Kramer     /* Count the number of interfaces per rank, not including those on the root */
2918a7f912aSStephan Kramer     PetscCall(PetscCalloc1(numProcs, &interfacesPerRank));
2928a7f912aSStephan Kramer     for (v = vStart; v < vEnd; v++) {
2938a7f912aSStephan Kramer       if (vertexNumber[v - vStart]) {
2948a7f912aSStephan Kramer         PetscCall(PetscSectionGetDof(rankSection, v, &shareCnt));
2958a7f912aSStephan Kramer         if (shareCnt) {
2968a7f912aSStephan Kramer           PetscCall(PetscSectionGetOffset(rankSection, v, &offset));
297ac530a7eSPierre Jolivet           for (j = 0; j < shareCnt; j++) interfacesPerRank[rankArray[offset + j]]++;
2988a7f912aSStephan Kramer         }
2998a7f912aSStephan Kramer       }
3008a7f912aSStephan Kramer     }
3018a7f912aSStephan Kramer     for (r = 0, k = 0, interfacesPerRank[rank] = 0; r < numProcs; r++) k += interfacesPerRank[r];
3028a7f912aSStephan Kramer 
3038a7f912aSStephan Kramer     /* Get the degree of the vertex */
3048a7f912aSStephan Kramer     PetscCall(PetscMalloc3(k, &interfaces_lv, k, &interfaces_gv, numProcs + 1, &interfacesOffset));
3058a7f912aSStephan Kramer     interfacesOffset[0] = 0;
3068a7f912aSStephan Kramer     for (r = 0; r < numProcs; r++) {
3078a7f912aSStephan Kramer       interfacesOffset[r + 1] = interfacesOffset[r] + interfacesPerRank[r];
3088a7f912aSStephan Kramer       if (interfacesPerRank[r]) numNgbRanks++;
3098a7f912aSStephan Kramer       interfacesPerRank[r] = 0;
3102bc0b47dSJoe Wallwork     }
311c51fff1bSJoe Wallwork 
3128a7f912aSStephan Kramer     /* Get the local and global vertex numbers at interfaces */
3138a7f912aSStephan Kramer     PetscCall(DMPlexGetVertexNumbering(dm, &globalVertexNum));
3149566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(globalVertexNum, &gV));
3158a7f912aSStephan Kramer     for (v = vStart; v < vEnd; v++) {
3168a7f912aSStephan Kramer       if (vertexNumber[v - vStart]) {
3178a7f912aSStephan Kramer         PetscCall(PetscSectionGetDof(rankSection, v, &shareCnt));
3188a7f912aSStephan Kramer         if (shareCnt) {
3198a7f912aSStephan Kramer           PetscCall(PetscSectionGetOffset(rankSection, v, &offset));
3208a7f912aSStephan Kramer           for (j = 0; j < shareCnt; j++) {
3218a7f912aSStephan Kramer             r = rankArray[offset + j];
3228a7f912aSStephan Kramer             if (r == rank) continue;
3238a7f912aSStephan Kramer             k                = interfacesOffset[r] + interfacesPerRank[r]++;
3248a7f912aSStephan Kramer             interfaces_lv[k] = vertexNumber[v - vStart];
3258a7f912aSStephan Kramer             gv               = gV[v - vStart];
3268a7f912aSStephan Kramer             interfaces_gv[k] = gv < 0 ? -gv : gv + 1;
3272bc0b47dSJoe Wallwork           }
3288a7f912aSStephan Kramer         }
3298a7f912aSStephan Kramer       }
3308a7f912aSStephan Kramer     }
3318a7f912aSStephan Kramer     PetscCall(PetscFree(interfacesPerRank));
3328a7f912aSStephan Kramer     PetscCall(PetscFree(rankArray));
3339566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(globalVertexNum, &gV));
3348a7f912aSStephan Kramer     PetscCall(PetscSectionDestroy(&rankSection));
3352bc0b47dSJoe Wallwork   }
3369566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&udm));
3378a7f912aSStephan Kramer   PetscCall(PetscFree(vertexNumber));
3382bc0b47dSJoe Wallwork 
3392bc0b47dSJoe Wallwork   /* Send the data to ParMmg and remesh */
3409566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoInsertion(dm, &noInsert));
3419566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoSwapping(dm, &noSwap));
3429566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoMovement(dm, &noMove));
3439566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoSurf(dm, &noSurf));
3449566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetVerbosity(dm, &verbosity));
3459566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNumIterations(dm, &numIter));
3469566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetGradationFactor(dm, &gradationFactor));
3479566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber));
348b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Init_parMesh, PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_pMesh, PMMG_ARG_pMet, PMMG_ARG_dim, 3, PMMG_ARG_MPIComm, comm, PMMG_ARG_end);
3498a7f912aSStephan Kramer   PetscCallMMG_NONSTANDARD(PMMG_Set_meshSize, parmesh, numUsedVertices, numCellsNotShared, 0, numFaceTags, 0, 0);
350b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_APImode, PMMG_APIDISTRIB_nodes);
351b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_noinsert, noInsert);
352b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_noswap, noSwap);
353b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_nomove, noMove);
354b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_nosurf, noSurf);
355b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_verbose, verbosity);
356b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_globalNum, 1);
357b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_niter, numIter);
358b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_dparameter, parmesh, PMMG_DPARAM_hgrad, gradationFactor);
359b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_dparameter, parmesh, PMMG_DPARAM_hausd, hausdorffNumber);
360b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_vertices, parmesh, vertices, verTags);
361b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_tetrahedra, parmesh, cells, cellTags);
362b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_triangles, parmesh, bdFaces, faceTags);
3638a7f912aSStephan Kramer   PetscCallMMG_NONSTANDARD(PMMG_Set_metSize, parmesh, MMG5_Vertex, numUsedVertices, MMG5_Tensor);
364b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_tensorMets, parmesh, metric);
365b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_numberOfNodeCommunicators, parmesh, numNgbRanks);
3668a7f912aSStephan Kramer   for (r = 0, c = 0; r < numProcs; ++r) {
3678a7f912aSStephan Kramer     if (interfacesOffset[r + 1] > interfacesOffset[r]) {
3688a7f912aSStephan Kramer       PetscCallMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicatorSize, parmesh, c, r, interfacesOffset[r + 1] - interfacesOffset[r]);
3698a7f912aSStephan Kramer       PetscCallMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicator_nodes, parmesh, c++, &interfaces_lv[interfacesOffset[r]], &interfaces_gv[interfacesOffset[r]], 1);
3708a7f912aSStephan Kramer     }
3712bc0b47dSJoe Wallwork   }
372b3377c29SStefano Zampini   PetscCallMMG(PMMG_parmmglib_distributed, parmesh);
3739566063dSJacob Faibussowitsch   PetscCall(PetscFree(cells));
3749566063dSJacob Faibussowitsch   PetscCall(PetscFree2(metric, vertices));
3759566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bdFaces, faceTags));
3769566063dSJacob Faibussowitsch   PetscCall(PetscFree2(verTags, cellTags));
3773a7d0413SPierre Jolivet   if (numProcs > 1) PetscCall(PetscFree3(interfaces_lv, interfaces_gv, interfacesOffset));
3782bc0b47dSJoe Wallwork 
3792f583692SJoe Wallwork   /* Retrieve mesh from Mmg */
3802bc0b47dSJoe Wallwork   numCornersNew = 4;
381b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Get_meshSize, parmesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0);
3829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dim * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
3839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3((dim + 1) * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
3849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dim * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
385b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Get_vertices, parmesh, verticesNew, verTagsNew, corners, requiredVer);
386b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Get_tetrahedra, parmesh, cellsNew, cellTagsNew, requiredCells);
387b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Get_triangles, parmesh, facesNew, faceTagsNew, requiredFaces);
3889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numVerticesNew, &owners, numVerticesNew, &gv_new));
389b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_globalNum, 1);
390b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Get_verticesGloNum, parmesh, gv_new, owners);
3912bc0b47dSJoe Wallwork   for (i = 0; i < dim * numFacesNew; ++i) facesNew[i] -= 1;
3920a7a67b9SPierre Jolivet   for (i = 0; i < (dim + 1) * numCellsNew; ++i) cellsNew[i] = gv_new[cellsNew[i] - 1] - 1;
3930a7a67b9SPierre Jolivet   for (i = 0, numVerticesNewLoc = 0; i < numVerticesNew; ++i) {
3942bc0b47dSJoe Wallwork     if (owners[i] == rank) numVerticesNewLoc++;
3952bc0b47dSJoe Wallwork   }
396*ec93882dSJames Wright   PetscCall(PetscMalloc1(numVerticesNewLoc * dim, &verticesNewLoc));
3972bc0b47dSJoe Wallwork   for (i = 0, c = 0; i < numVerticesNew; i++) {
3982bc0b47dSJoe Wallwork     if (owners[i] == rank) {
3992bc0b47dSJoe Wallwork       for (j = 0; j < dim; ++j) verticesNewLoc[dim * c + j] = verticesNew[dim * i + j];
4002bc0b47dSJoe Wallwork       c++;
4012bc0b47dSJoe Wallwork     }
4022bc0b47dSJoe Wallwork   }
4032f583692SJoe Wallwork 
4042f583692SJoe Wallwork   /* Reorder for consistency with DMPlex */
4059566063dSJacob Faibussowitsch   for (i = 0; i < numCellsNew; ++i) PetscCall(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4 * i]));
4062f583692SJoe Wallwork 
4072f583692SJoe Wallwork   /* Create new plex */
4089566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNewLoc, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNewLoc, NULL, &verticesNewSorted, dmNew));
409b3377c29SStefano Zampini   PetscCallMMG_NONSTANDARD(PMMG_Free_all, PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_end);
4109566063dSJacob Faibussowitsch   PetscCall(PetscFree4(verticesNew, verTagsNew, corners, requiredVer));
411c1dc6da0SJoe Wallwork 
412c1dc6da0SJoe Wallwork   /* Get adapted mesh information */
4139566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
4149566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
4159566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
4162bc0b47dSJoe Wallwork 
4172bc0b47dSJoe Wallwork   /* Rebuild boundary label */
4189566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName));
4199566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew));
4202bc0b47dSJoe Wallwork   for (i = 0; i < numFacesNew; i++) {
421c1dc6da0SJoe Wallwork     PetscBool       hasTag = PETSC_FALSE;
4222bc0b47dSJoe Wallwork     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
4232bc0b47dSJoe Wallwork     const PetscInt *coveredPoints = NULL;
4242bc0b47dSJoe Wallwork 
4252bc0b47dSJoe Wallwork     for (j = 0; j < dim; ++j) {
4262bc0b47dSJoe Wallwork       lv = facesNew[i * dim + j];
4272bc0b47dSJoe Wallwork       gv = gv_new[lv] - 1;
4289566063dSJacob Faibussowitsch       PetscCall(PetscFindInt(gv, numVerticesNew, verticesNewSorted, &lv));
4292bc0b47dSJoe Wallwork       facePoints[j] = lv + vStart;
4302bc0b47dSJoe Wallwork     }
4319566063dSJacob Faibussowitsch     PetscCall(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
4322bc0b47dSJoe Wallwork     for (j = 0; j < numCoveredPoints; ++j) {
4332bc0b47dSJoe Wallwork       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
4342bc0b47dSJoe Wallwork         numFaces++;
4352bc0b47dSJoe Wallwork         f = j;
4362bc0b47dSJoe Wallwork       }
4372bc0b47dSJoe Wallwork     }
4385f80ce2aSJacob Faibussowitsch     PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces);
4399566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(bdLabel, faceTagsNew[i], &hasTag));
4409566063dSJacob Faibussowitsch     if (hasTag) PetscCall(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]));
4419566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
4422bc0b47dSJoe Wallwork   }
4439566063dSJacob Faibussowitsch   PetscCall(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces));
4449566063dSJacob Faibussowitsch   PetscCall(PetscFree2(owners, gv_new));
445*ec93882dSJames Wright   PetscCall(PetscFree(verticesNewLoc));
4469566063dSJacob Faibussowitsch   if (flg) PetscCall(DMLabelDestroy(&bdLabel));
4472bc0b47dSJoe Wallwork 
4489fe9e680SJoe Wallwork   /* Rebuild cell labels */
4499566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName));
4509566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew));
4519566063dSJacob Faibussowitsch   for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c - cStart]));
4529566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellsNew, cellTagsNew, requiredCells));
4533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4542bc0b47dSJoe Wallwork }
455