xref: /petsc/src/dm/impls/plex/adaptors/parmmg/parmmgadapt.c (revision b2c4de1d22ff155d2c936e49c947998b8cbd29f0)
15f80ce2aSJacob Faibussowitsch #include "../mmgcommon.h" /*I      "petscdmplex.h"   I*/
22bc0b47dSJoe Wallwork #include <parmmg/libparmmg.h>
32bc0b47dSJoe Wallwork 
4*b2c4de1dSPierre Jolivet PetscBool  ParMmgCite       = PETSC_FALSE;
5*b2c4de1dSPierre Jolivet const char ParMmgCitation[] = "@techreport{cirrottola:hal-02386837,\n"
6*b2c4de1dSPierre Jolivet                               "  title       = {Parallel unstructured mesh adaptation using iterative remeshing and repartitioning},\n"
7*b2c4de1dSPierre Jolivet                               "  institution = {Inria Bordeaux},\n"
8*b2c4de1dSPierre Jolivet                               "  author      = {L. Cirrottola and A. Froehly},\n"
9*b2c4de1dSPierre Jolivet                               "  number      = {9307},\n"
10*b2c4de1dSPierre Jolivet                               "  note        = {\\url{https://hal.inria.fr/hal-02386837}},\n"
11*b2c4de1dSPierre Jolivet                               "  year        = {2019}\n}\n";
12*b2c4de1dSPierre Jolivet 
139371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew) {
140a7a67b9SPierre Jolivet   MPI_Comm           comm;
152bc0b47dSJoe Wallwork   const char        *bdName = "_boundary_";
169fe9e680SJoe Wallwork   const char        *rgName = "_regions_";
172bc0b47dSJoe Wallwork   DM                 udm, cdm;
18c1dc6da0SJoe Wallwork   DMLabel            bdLabelNew, rgLabelNew;
199fe9e680SJoe Wallwork   const char        *bdLabelName, *rgLabelName;
20c1dc6da0SJoe Wallwork   IS                 globalVertexNum;
212bc0b47dSJoe Wallwork   PetscSection       coordSection;
222bc0b47dSJoe Wallwork   Vec                coordinates;
232bc0b47dSJoe Wallwork   PetscSF            sf;
242bc0b47dSJoe Wallwork   const PetscScalar *coords, *met;
25ae8b063eSJoe Wallwork   PetscReal         *vertices, *metric, *verticesNew, *verticesNewLoc, gradationFactor, hausdorffNumber;
26c1dc6da0SJoe Wallwork   PetscInt          *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew;
27c1dc6da0SJoe Wallwork   PetscInt          *bdFaces, *faceTags, *facesNew, *faceTagsNew;
28c1dc6da0SJoe Wallwork   PetscInt          *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces;
29c1dc6da0SJoe Wallwork   PetscInt           cStart, cEnd, c, numCells, fStart, fEnd, f, numFaceTags, vStart, vEnd, v, numVertices;
30c1dc6da0SJoe Wallwork   PetscInt           dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, numIter;
312bc0b47dSJoe Wallwork   PetscInt          *numVerInterfaces, *ngbRanks, *verNgbRank, *interfaces_lv, *interfaces_gv, *intOffset;
32c1dc6da0SJoe Wallwork   PetscInt           niranks, nrranks, numNgbRanks, numVerNgbRanksTotal, count, sliceSize, p, r, n, lv, gv;
33c1dc6da0SJoe Wallwork   PetscInt          *gv_new, *owners, *verticesNewSorted, pStart, pEnd;
342bc0b47dSJoe Wallwork   PetscInt           numCellsNew, numVerticesNew, numCornersNew, numFacesNew, numVerticesNewLoc;
35c1dc6da0SJoe Wallwork   const PetscInt    *gV, *ioffset, *irootloc, *roffset, *rmine, *rremote;
3676f360caSJoe Wallwork   PetscBool          flg = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform;
372bc0b47dSJoe Wallwork   const PetscMPIInt *iranks, *rranks;
382bc0b47dSJoe Wallwork   PetscMPIInt        numProcs, rank;
392bc0b47dSJoe Wallwork   PMMG_pParMesh      parmesh = NULL;
402bc0b47dSJoe Wallwork 
412bc0b47dSJoe Wallwork   PetscFunctionBegin;
42*b2c4de1dSPierre Jolivet   PetscCall(PetscCitationsRegister(ParMmgCitation, &ParMmgCite));
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &numProcs));
459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
462bc0b47dSJoe Wallwork   if (bdLabel) {
479566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)bdLabel, &bdLabelName));
489566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(bdLabelName, bdName, &flg));
4928b400f6SJacob Faibussowitsch     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
502bc0b47dSJoe Wallwork   }
519fe9e680SJoe Wallwork   if (rgLabel) {
529566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)rgLabel, &rgLabelName));
539566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(rgLabelName, rgName, &flg));
5428b400f6SJacob Faibussowitsch     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
559fe9e680SJoe Wallwork   }
562bc0b47dSJoe Wallwork 
572bc0b47dSJoe Wallwork   /* Get mesh information */
589566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
595f80ce2aSJacob Faibussowitsch   PetscCheck(dim == 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "ParMmg only works in 3D.");
602bc0b47dSJoe Wallwork   Neq = (dim * (dim + 1)) / 2;
619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
639566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
649566063dSJacob Faibussowitsch   PetscCall(DMPlexUninterpolate(dm, &udm));
659566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
662bc0b47dSJoe Wallwork   numCells    = cEnd - cStart;
672bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
682bc0b47dSJoe Wallwork 
692bc0b47dSJoe Wallwork   /* Get cell offsets */
709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numCells * maxConeSize, &cells));
712bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
722bc0b47dSJoe Wallwork     const PetscInt *cone;
732bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
742bc0b47dSJoe Wallwork 
759566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(udm, c, &coneSize));
769566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(udm, c, &cone));
772bc0b47dSJoe Wallwork     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart + 1;
782bc0b47dSJoe Wallwork   }
792bc0b47dSJoe Wallwork 
802bc0b47dSJoe Wallwork   /* Get vertex coordinate array */
819566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
829566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(cdm, &coordSection));
839566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numVertices * Neq, &metric, dim * numVertices, &vertices));
862bc0b47dSJoe Wallwork   for (v = 0; v < vEnd - vStart; ++v) {
879566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(coordSection, v + vStart, &off));
882bc0b47dSJoe Wallwork     for (i = 0; i < dim; ++i) vertices[dim * v + i] = PetscRealPart(coords[off + i]);
892bc0b47dSJoe Wallwork   }
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
912bc0b47dSJoe Wallwork 
92c1dc6da0SJoe Wallwork   /* Get face tags */
93c1dc6da0SJoe Wallwork   if (!bdLabel) {
94c1dc6da0SJoe Wallwork     flg = PETSC_TRUE;
959566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel));
969566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
97c1dc6da0SJoe Wallwork   }
989566063dSJacob Faibussowitsch   PetscCall(DMLabelGetBounds(bdLabel, &pStart, &pEnd));
99c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
100c1dc6da0SJoe Wallwork     PetscBool hasPoint;
101c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1022bc0b47dSJoe Wallwork 
1039566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
104c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
105c1dc6da0SJoe Wallwork     numFaceTags++;
106c1dc6da0SJoe Wallwork 
1079566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1082bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize * 2; cl += 2) {
1092bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1102bc0b47dSJoe Wallwork     }
1119566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1122bc0b47dSJoe Wallwork   }
1139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags));
114c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
115c1dc6da0SJoe Wallwork     PetscBool hasPoint;
116c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1172bc0b47dSJoe Wallwork 
1189566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
119c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
120c1dc6da0SJoe Wallwork 
1219566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1222bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize * 2; cl += 2) {
1232bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart + 1;
1242bc0b47dSJoe Wallwork     }
1259566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1269566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]));
1272bc0b47dSJoe Wallwork   }
1282bc0b47dSJoe Wallwork 
1299fe9e680SJoe Wallwork   /* Get cell tags */
1309566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(numVertices, &verTags, numCells, &cellTags));
1319fe9e680SJoe Wallwork   if (rgLabel) {
1329566063dSJacob Faibussowitsch     for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelGetValue(rgLabel, c, &cellTags[c]));
1339fe9e680SJoe Wallwork   }
1349fe9e680SJoe Wallwork 
1352bc0b47dSJoe Wallwork   /* Get metric */
1369566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(vertexMetric, &met));
1389566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1399566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1402bc0b47dSJoe Wallwork   for (v = 0; v < (vEnd - vStart); ++v) {
1412bc0b47dSJoe Wallwork     for (i = 0, k = 0; i < dim; ++i) {
1420a7a67b9SPierre Jolivet       for (j = i; j < dim; ++j, ++k) {
14393520041SJoe Wallwork         if (isotropic) {
14493520041SJoe Wallwork           if (i == j) {
14593520041SJoe Wallwork             if (uniform) metric[Neq * v + k] = PetscRealPart(met[0]);
14693520041SJoe Wallwork             else metric[Neq * v + k] = PetscRealPart(met[v]);
14793520041SJoe Wallwork           } else metric[Neq * v + k] = 0.0;
1480a7a67b9SPierre Jolivet         } else metric[Neq * v + k] = PetscRealPart(met[dim * dim * v + dim * i + j]);
1492bc0b47dSJoe Wallwork       }
1502bc0b47dSJoe Wallwork     }
1512bc0b47dSJoe Wallwork   }
1529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(vertexMetric, &met));
1532bc0b47dSJoe Wallwork 
1542bc0b47dSJoe Wallwork   /* Build ParMMG communicators: the list of vertices between two partitions  */
1552bc0b47dSJoe Wallwork   niranks = nrranks = 0;
1562bc0b47dSJoe Wallwork   numNgbRanks       = 0;
1572bc0b47dSJoe Wallwork   if (numProcs > 1) {
1589566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
1599566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(sf));
1609566063dSJacob Faibussowitsch     PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
1619566063dSJacob Faibussowitsch     PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
1629566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(numProcs, &numVerInterfaces));
1632bc0b47dSJoe Wallwork 
164c51fff1bSJoe Wallwork     /* Count number of roots associated with each leaf */
1652bc0b47dSJoe Wallwork     for (r = 0; r < niranks; ++r) {
1662bc0b47dSJoe Wallwork       for (i = ioffset[r], count = 0; i < ioffset[r + 1]; ++i) {
1672bc0b47dSJoe Wallwork         if (irootloc[i] >= vStart && irootloc[i] < vEnd) count++;
1682bc0b47dSJoe Wallwork       }
1692bc0b47dSJoe Wallwork       numVerInterfaces[iranks[r]] += count;
1702bc0b47dSJoe Wallwork     }
171c51fff1bSJoe Wallwork 
172c51fff1bSJoe Wallwork     /* Count number of leaves associated with each root */
1732bc0b47dSJoe Wallwork     for (r = 0; r < nrranks; ++r) {
1742bc0b47dSJoe Wallwork       for (i = roffset[r], count = 0; i < roffset[r + 1]; ++i) {
1752bc0b47dSJoe Wallwork         if (rmine[i] >= vStart && rmine[i] < vEnd) count++;
1762bc0b47dSJoe Wallwork       }
1772bc0b47dSJoe Wallwork       numVerInterfaces[rranks[r]] += count;
1782bc0b47dSJoe Wallwork     }
179c51fff1bSJoe Wallwork 
180c51fff1bSJoe Wallwork     /* Count global number of ranks */
1812bc0b47dSJoe Wallwork     for (p = 0; p < numProcs; ++p) {
1822bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) numNgbRanks++;
1832bc0b47dSJoe Wallwork     }
184c51fff1bSJoe Wallwork 
185c51fff1bSJoe Wallwork     /* Provide numbers of vertex interfaces */
1869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(numNgbRanks, &ngbRanks, numNgbRanks, &verNgbRank));
1872bc0b47dSJoe Wallwork     for (p = 0, n = 0; p < numProcs; ++p) {
1882bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) {
1892bc0b47dSJoe Wallwork         ngbRanks[n]   = p;
1902bc0b47dSJoe Wallwork         verNgbRank[n] = numVerInterfaces[p];
1912bc0b47dSJoe Wallwork         n++;
1922bc0b47dSJoe Wallwork       }
1932bc0b47dSJoe Wallwork     }
1942bc0b47dSJoe Wallwork     numVerNgbRanksTotal = 0;
1952bc0b47dSJoe Wallwork     for (p = 0; p < numNgbRanks; ++p) numVerNgbRanksTotal += verNgbRank[p];
1962bc0b47dSJoe Wallwork 
1972bc0b47dSJoe Wallwork     /* For each neighbor, fill in interface arrays */
1989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(numVerNgbRanksTotal, &interfaces_lv, numVerNgbRanksTotal, &interfaces_gv, numNgbRanks + 1, &intOffset));
1992bc0b47dSJoe Wallwork     intOffset[0] = 0;
2002bc0b47dSJoe Wallwork     for (p = 0, r = 0, i = 0; p < numNgbRanks; ++p) {
2012bc0b47dSJoe Wallwork       intOffset[p + 1] = intOffset[p];
202c51fff1bSJoe Wallwork 
203c51fff1bSJoe Wallwork       /* Leaf case */
2042bc0b47dSJoe Wallwork       if (iranks && iranks[i] == ngbRanks[p]) {
2052bc0b47dSJoe Wallwork         /* Add the right slice of irootloc at the right place */
2062bc0b47dSJoe Wallwork         sliceSize = ioffset[i + 1] - ioffset[i];
2072bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
2085f80ce2aSJacob Faibussowitsch           PetscCheck(ioffset[i] + j < ioffset[niranks], comm, PETSC_ERR_ARG_OUTOFRANGE, "Leaf index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", ioffset[i] + j, ioffset[niranks]);
2092bc0b47dSJoe Wallwork           v = irootloc[ioffset[i] + j];
2102bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
2115f80ce2aSJacob Faibussowitsch             PetscCheck(intOffset[p + 1] + count < numVerNgbRanksTotal, comm, PETSC_ERR_ARG_OUTOFRANGE, "Leaf interface index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", intOffset[p + 1] + count, numVerNgbRanksTotal);
2122bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p + 1] + count] = v - vStart;
2132bc0b47dSJoe Wallwork             count++;
2142bc0b47dSJoe Wallwork           }
2152bc0b47dSJoe Wallwork         }
2162bc0b47dSJoe Wallwork         intOffset[p + 1] += count;
2172bc0b47dSJoe Wallwork         i++;
2182bc0b47dSJoe Wallwork       }
219c51fff1bSJoe Wallwork 
220c51fff1bSJoe Wallwork       /* Root case */
2212bc0b47dSJoe Wallwork       if (rranks && rranks[r] == ngbRanks[p]) {
2222bc0b47dSJoe Wallwork         /* Add the right slice of rmine at the right place */
2232bc0b47dSJoe Wallwork         sliceSize = roffset[r + 1] - roffset[r];
2242bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
2255f80ce2aSJacob Faibussowitsch           PetscCheck(roffset[r] + j < roffset[nrranks], comm, PETSC_ERR_ARG_OUTOFRANGE, "Root index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", roffset[r] + j, roffset[nrranks]);
2262bc0b47dSJoe Wallwork           v = rmine[roffset[r] + j];
2272bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
2285f80ce2aSJacob Faibussowitsch             PetscCheck(intOffset[p + 1] + count < numVerNgbRanksTotal, comm, PETSC_ERR_ARG_OUTOFRANGE, "Root interface index %" PetscInt_FMT " out of range (expected < %" PetscInt_FMT ")", intOffset[p + 1] + count, numVerNgbRanksTotal);
2292bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p + 1] + count] = v - vStart;
2302bc0b47dSJoe Wallwork             count++;
2312bc0b47dSJoe Wallwork           }
2322bc0b47dSJoe Wallwork         }
2332bc0b47dSJoe Wallwork         intOffset[p + 1] += count;
2342bc0b47dSJoe Wallwork         r++;
2352bc0b47dSJoe Wallwork       }
236c51fff1bSJoe Wallwork 
237c51fff1bSJoe Wallwork       /* Check validity of offsets */
2385f80ce2aSJacob Faibussowitsch       PetscCheck(intOffset[p + 1] == intOffset[p] + verNgbRank[p], comm, PETSC_ERR_ARG_OUTOFRANGE, "Missing offsets (expected %" PetscInt_FMT ", got %" PetscInt_FMT ")", intOffset[p] + verNgbRank[p], intOffset[p + 1]);
2392bc0b47dSJoe Wallwork     }
2409566063dSJacob Faibussowitsch     PetscCall(DMPlexGetVertexNumbering(udm, &globalVertexNum));
2419566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(globalVertexNum, &gV));
2422bc0b47dSJoe Wallwork     for (i = 0; i < numVerNgbRanksTotal; ++i) {
2432bc0b47dSJoe Wallwork       v                = gV[interfaces_lv[i]];
2442bc0b47dSJoe Wallwork       interfaces_gv[i] = v < 0 ? -v - 1 : v;
2452bc0b47dSJoe Wallwork       interfaces_lv[i] += 1;
2462bc0b47dSJoe Wallwork       interfaces_gv[i] += 1;
2472bc0b47dSJoe Wallwork     }
2489566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(globalVertexNum, &gV));
2499566063dSJacob Faibussowitsch     PetscCall(PetscFree(numVerInterfaces));
2502bc0b47dSJoe Wallwork   }
2519566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&udm));
2522bc0b47dSJoe Wallwork 
2532bc0b47dSJoe Wallwork   /* Send the data to ParMmg and remesh */
2549566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoInsertion(dm, &noInsert));
2559566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoSwapping(dm, &noSwap));
2569566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoMovement(dm, &noMove));
2579566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoSurf(dm, &noSurf));
2589566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetVerbosity(dm, &verbosity));
2599566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetNumIterations(dm, &numIter));
2609566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetGradationFactor(dm, &gradationFactor));
2619566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber));
2629566063dSJacob Faibussowitsch   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));
2639566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_meshSize(parmesh, numVertices, numCells, 0, numFaceTags, 0, 0));
2649566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_APImode, PMMG_APIDISTRIB_nodes));
2659566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noinsert, noInsert));
2669566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noswap, noSwap));
2679566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nomove, noMove));
2689566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nosurf, noSurf));
2699566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_verbose, verbosity));
2709566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1));
2719566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_niter, numIter));
2729566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hgrad, gradationFactor));
2739566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hausd, hausdorffNumber));
2749566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_vertices(parmesh, vertices, verTags));
2759566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_tetrahedra(parmesh, cells, cellTags));
2769566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_triangles(parmesh, bdFaces, faceTags));
2779566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_metSize(parmesh, MMG5_Vertex, numVertices, MMG5_Tensor));
2789566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_tensorMets(parmesh, metric));
2799566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_numberOfNodeCommunicators(parmesh, numNgbRanks));
2802bc0b47dSJoe Wallwork   for (c = 0; c < numNgbRanks; ++c) {
2819566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicatorSize(parmesh, c, ngbRanks[c], intOffset[c + 1] - intOffset[c]));
2829566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicator_nodes(parmesh, c, &interfaces_lv[intOffset[c]], &interfaces_gv[intOffset[c]], 1));
2832bc0b47dSJoe Wallwork   }
2849566063dSJacob Faibussowitsch   PetscCallMMG(PMMG_parmmglib_distributed(parmesh));
2859566063dSJacob Faibussowitsch   PetscCall(PetscFree(cells));
2869566063dSJacob Faibussowitsch   PetscCall(PetscFree2(metric, vertices));
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bdFaces, faceTags));
2889566063dSJacob Faibussowitsch   PetscCall(PetscFree2(verTags, cellTags));
2892bc0b47dSJoe Wallwork   if (numProcs > 1) {
2909566063dSJacob Faibussowitsch     PetscCall(PetscFree2(ngbRanks, verNgbRank));
2919566063dSJacob Faibussowitsch     PetscCall(PetscFree3(interfaces_lv, interfaces_gv, intOffset));
2922bc0b47dSJoe Wallwork   }
2932bc0b47dSJoe Wallwork 
2942f583692SJoe Wallwork   /* Retrieve mesh from Mmg */
2952bc0b47dSJoe Wallwork   numCornersNew = 4;
2969566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Get_meshSize(parmesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0));
2979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dim * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
2989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3((dim + 1) * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
2999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dim * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
3009566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Get_vertices(parmesh, verticesNew, verTagsNew, corners, requiredVer));
3019566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Get_tetrahedra(parmesh, cellsNew, cellTagsNew, requiredCells));
3029566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Get_triangles(parmesh, facesNew, faceTagsNew, requiredFaces));
3039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numVerticesNew, &owners, numVerticesNew, &gv_new));
3049566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1));
3059566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Get_verticesGloNum(parmesh, gv_new, owners));
3062bc0b47dSJoe Wallwork   for (i = 0; i < dim * numFacesNew; ++i) facesNew[i] -= 1;
3070a7a67b9SPierre Jolivet   for (i = 0; i < (dim + 1) * numCellsNew; ++i) cellsNew[i] = gv_new[cellsNew[i] - 1] - 1;
3080a7a67b9SPierre Jolivet   for (i = 0, numVerticesNewLoc = 0; i < numVerticesNew; ++i) {
3092bc0b47dSJoe Wallwork     if (owners[i] == rank) numVerticesNewLoc++;
3102bc0b47dSJoe Wallwork   }
3119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numVerticesNewLoc * dim, &verticesNewLoc, numVerticesNew, &verticesNewSorted));
3122bc0b47dSJoe Wallwork   for (i = 0, c = 0; i < numVerticesNew; i++) {
3132bc0b47dSJoe Wallwork     if (owners[i] == rank) {
3142bc0b47dSJoe Wallwork       for (j = 0; j < dim; ++j) verticesNewLoc[dim * c + j] = verticesNew[dim * i + j];
3152bc0b47dSJoe Wallwork       c++;
3162bc0b47dSJoe Wallwork     }
3172bc0b47dSJoe Wallwork   }
3182f583692SJoe Wallwork 
3192f583692SJoe Wallwork   /* Reorder for consistency with DMPlex */
3209566063dSJacob Faibussowitsch   for (i = 0; i < numCellsNew; ++i) PetscCall(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4 * i]));
3212f583692SJoe Wallwork 
3222f583692SJoe Wallwork   /* Create new plex */
3239566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNewLoc, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNewLoc, NULL, &verticesNewSorted, dmNew));
3249566063dSJacob Faibussowitsch   PetscCallMMG_NONSTANDARD(PMMG_Free_all(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_end));
3259566063dSJacob Faibussowitsch   PetscCall(PetscFree4(verticesNew, verTagsNew, corners, requiredVer));
326c1dc6da0SJoe Wallwork 
327c1dc6da0SJoe Wallwork   /* Get adapted mesh information */
3289566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
3299566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
3309566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
3312bc0b47dSJoe Wallwork 
3322bc0b47dSJoe Wallwork   /* Rebuild boundary label */
3339566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName));
3349566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew));
3352bc0b47dSJoe Wallwork   for (i = 0; i < numFacesNew; i++) {
336c1dc6da0SJoe Wallwork     PetscBool       hasTag = PETSC_FALSE;
3372bc0b47dSJoe Wallwork     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
3382bc0b47dSJoe Wallwork     const PetscInt *coveredPoints = NULL;
3392bc0b47dSJoe Wallwork 
3402bc0b47dSJoe Wallwork     for (j = 0; j < dim; ++j) {
3412bc0b47dSJoe Wallwork       lv = facesNew[i * dim + j];
3422bc0b47dSJoe Wallwork       gv = gv_new[lv] - 1;
3439566063dSJacob Faibussowitsch       PetscCall(PetscFindInt(gv, numVerticesNew, verticesNewSorted, &lv));
3442bc0b47dSJoe Wallwork       facePoints[j] = lv + vStart;
3452bc0b47dSJoe Wallwork     }
3469566063dSJacob Faibussowitsch     PetscCall(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
3472bc0b47dSJoe Wallwork     for (j = 0; j < numCoveredPoints; ++j) {
3482bc0b47dSJoe Wallwork       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
3492bc0b47dSJoe Wallwork         numFaces++;
3502bc0b47dSJoe Wallwork         f = j;
3512bc0b47dSJoe Wallwork       }
3522bc0b47dSJoe Wallwork     }
3535f80ce2aSJacob Faibussowitsch     PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces);
3549566063dSJacob Faibussowitsch     PetscCall(DMLabelHasStratum(bdLabel, faceTagsNew[i], &hasTag));
3559566063dSJacob Faibussowitsch     if (hasTag) PetscCall(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]));
3569566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
3572bc0b47dSJoe Wallwork   }
3589566063dSJacob Faibussowitsch   PetscCall(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces));
3599566063dSJacob Faibussowitsch   PetscCall(PetscFree2(owners, gv_new));
3609566063dSJacob Faibussowitsch   PetscCall(PetscFree2(verticesNewLoc, verticesNewSorted));
3619566063dSJacob Faibussowitsch   if (flg) PetscCall(DMLabelDestroy(&bdLabel));
3622bc0b47dSJoe Wallwork 
3639fe9e680SJoe Wallwork   /* Rebuild cell labels */
3649566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName));
3659566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew));
3669566063dSJacob Faibussowitsch   for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c - cStart]));
3679566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellsNew, cellTagsNew, requiredCells));
368c1dc6da0SJoe Wallwork 
3692bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
3702bc0b47dSJoe Wallwork }
371