xref: /petsc/src/dm/impls/plex/adaptors/parmmg/parmmgadapt.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
15f80ce2aSJacob Faibussowitsch #include "../mmgcommon.h"   /*I      "petscdmplex.h"   I*/
22bc0b47dSJoe Wallwork #include <parmmg/libparmmg.h>
32bc0b47dSJoe Wallwork 
49fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew)
52bc0b47dSJoe Wallwork {
60a7a67b9SPierre Jolivet   MPI_Comm           comm;
72bc0b47dSJoe Wallwork   const char        *bdName = "_boundary_";
89fe9e680SJoe Wallwork   const char        *rgName = "_regions_";
92bc0b47dSJoe Wallwork   DM                 udm, cdm;
10c1dc6da0SJoe Wallwork   DMLabel            bdLabelNew, rgLabelNew;
119fe9e680SJoe Wallwork   const char        *bdLabelName, *rgLabelName;
12c1dc6da0SJoe Wallwork   IS                 globalVertexNum;
132bc0b47dSJoe Wallwork   PetscSection       coordSection;
142bc0b47dSJoe Wallwork   Vec                coordinates;
152bc0b47dSJoe Wallwork   PetscSF            sf;
162bc0b47dSJoe Wallwork   const PetscScalar *coords, *met;
17ae8b063eSJoe Wallwork   PetscReal         *vertices, *metric, *verticesNew, *verticesNewLoc, gradationFactor, hausdorffNumber;
18c1dc6da0SJoe Wallwork   PetscInt          *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew;
19c1dc6da0SJoe Wallwork   PetscInt          *bdFaces, *faceTags, *facesNew, *faceTagsNew;
20c1dc6da0SJoe Wallwork   PetscInt          *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces;
21c1dc6da0SJoe Wallwork   PetscInt           cStart, cEnd, c, numCells, fStart, fEnd, f, numFaceTags, vStart, vEnd, v, numVertices;
22c1dc6da0SJoe Wallwork   PetscInt           dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, numIter;
232bc0b47dSJoe Wallwork   PetscInt          *numVerInterfaces, *ngbRanks, *verNgbRank, *interfaces_lv, *interfaces_gv, *intOffset;
24c1dc6da0SJoe Wallwork   PetscInt           niranks, nrranks, numNgbRanks, numVerNgbRanksTotal, count, sliceSize, p, r, n, lv, gv;
25c1dc6da0SJoe Wallwork   PetscInt          *gv_new, *owners, *verticesNewSorted, pStart, pEnd;
262bc0b47dSJoe Wallwork   PetscInt           numCellsNew, numVerticesNew, numCornersNew, numFacesNew, numVerticesNewLoc;
27c1dc6da0SJoe Wallwork   const PetscInt    *gV, *ioffset, *irootloc, *roffset, *rmine, *rremote;
2876f360caSJoe Wallwork   PetscBool          flg = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform;
292bc0b47dSJoe Wallwork   const PetscMPIInt *iranks, *rranks;
302bc0b47dSJoe Wallwork   PetscMPIInt        numProcs, rank;
312bc0b47dSJoe Wallwork   PMMG_pParMesh      parmesh = NULL;
322bc0b47dSJoe Wallwork 
332bc0b47dSJoe Wallwork   PetscFunctionBegin;
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
355f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &numProcs));
365f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
372bc0b47dSJoe Wallwork   if (bdLabel) {
385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetName((PetscObject) bdLabel, &bdLabelName));
395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(bdLabelName, bdName, &flg));
40*28b400f6SJacob Faibussowitsch     PetscCheck(!flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
412bc0b47dSJoe Wallwork   }
429fe9e680SJoe Wallwork   if (rgLabel) {
435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetName((PetscObject) rgLabel, &rgLabelName));
445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(rgLabelName, rgName, &flg));
45*28b400f6SJacob Faibussowitsch     PetscCheck(!flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
469fe9e680SJoe Wallwork   }
472bc0b47dSJoe Wallwork 
482bc0b47dSJoe Wallwork   /* Get mesh information */
495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
505f80ce2aSJacob Faibussowitsch   PetscCheck(dim == 3,comm, PETSC_ERR_ARG_OUTOFRANGE, "ParMmg only works in 3D.");
512bc0b47dSJoe Wallwork   Neq  = (dim*(dim+1))/2;
525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexUninterpolate(dm, &udm));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
572bc0b47dSJoe Wallwork   numCells    = cEnd - cStart;
582bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
592bc0b47dSJoe Wallwork 
602bc0b47dSJoe Wallwork   /* Get cell offsets */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numCells*maxConeSize, &cells));
622bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
632bc0b47dSJoe Wallwork     const PetscInt *cone;
642bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
652bc0b47dSJoe Wallwork 
665f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetConeSize(udm, c, &coneSize));
675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCone(udm, c, &cone));
682bc0b47dSJoe Wallwork     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart+1;
692bc0b47dSJoe Wallwork   }
702bc0b47dSJoe Wallwork 
712bc0b47dSJoe Wallwork   /* Get vertex coordinate array */
725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(dm, &cdm));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(cdm, &coordSection));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(coordinates, &coords));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(numVertices*Neq, &metric, dim*numVertices, &vertices));
772bc0b47dSJoe Wallwork   for (v = 0; v < vEnd-vStart; ++v) {
785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(coordSection, v+vStart, &off));
792bc0b47dSJoe Wallwork     for (i = 0; i < dim; ++i) vertices[dim*v+i] = PetscRealPart(coords[off+i]);
802bc0b47dSJoe Wallwork   }
815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(coordinates, &coords));
822bc0b47dSJoe Wallwork 
83c1dc6da0SJoe Wallwork   /* Get face tags */
84c1dc6da0SJoe Wallwork   if (!bdLabel) {
85c1dc6da0SJoe Wallwork     flg = PETSC_TRUE;
865f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel));
875f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
88c1dc6da0SJoe Wallwork   }
895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelGetBounds(bdLabel, &pStart, &pEnd));
90c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
91c1dc6da0SJoe Wallwork     PetscBool hasPoint;
92c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
932bc0b47dSJoe Wallwork 
945f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasPoint(bdLabel, f, &hasPoint));
95c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
96c1dc6da0SJoe Wallwork     numFaceTags++;
97c1dc6da0SJoe Wallwork 
985f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
992bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1002bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1012bc0b47dSJoe Wallwork     }
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1032bc0b47dSJoe Wallwork   }
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags));
105c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
106c1dc6da0SJoe Wallwork     PetscBool hasPoint;
107c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1082bc0b47dSJoe Wallwork 
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasPoint(bdLabel, f, &hasPoint));
110c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
111c1dc6da0SJoe Wallwork 
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1132bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1142bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart+1;
1152bc0b47dSJoe Wallwork     }
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]));
1182bc0b47dSJoe Wallwork   }
1192bc0b47dSJoe Wallwork 
1209fe9e680SJoe Wallwork   /* Get cell tags */
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc2(numVertices, &verTags, numCells, &cellTags));
1229fe9e680SJoe Wallwork   if (rgLabel) {
1235f80ce2aSJacob Faibussowitsch     for (c = cStart; c < cEnd; ++c) CHKERRQ(DMLabelGetValue(rgLabel, c, &cellTags[c]));
1249fe9e680SJoe Wallwork   }
1259fe9e680SJoe Wallwork 
1262bc0b47dSJoe Wallwork   /* Get metric */
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(vertexMetric, &met));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricIsUniform(dm, &uniform));
1312bc0b47dSJoe Wallwork   for (v = 0; v < (vEnd-vStart); ++v) {
1322bc0b47dSJoe Wallwork     for (i = 0, k = 0; i < dim; ++i) {
1330a7a67b9SPierre Jolivet       for (j = i; j < dim; ++j, ++k) {
13493520041SJoe Wallwork         if (isotropic) {
13593520041SJoe Wallwork           if (i == j) {
13693520041SJoe Wallwork             if (uniform) metric[Neq*v+k] = PetscRealPart(met[0]);
13793520041SJoe Wallwork             else metric[Neq*v+k] = PetscRealPart(met[v]);
13893520041SJoe Wallwork           } else metric[Neq*v+k] = 0.0;
1390a7a67b9SPierre Jolivet         } else metric[Neq*v+k] = PetscRealPart(met[dim*dim*v+dim*i+j]);
1402bc0b47dSJoe Wallwork       }
1412bc0b47dSJoe Wallwork     }
1422bc0b47dSJoe Wallwork   }
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(vertexMetric, &met));
1442bc0b47dSJoe Wallwork 
1452bc0b47dSJoe Wallwork   /* Build ParMMG communicators: the list of vertices between two partitions  */
1462bc0b47dSJoe Wallwork   niranks = nrranks = 0;
1472bc0b47dSJoe Wallwork   numNgbRanks = 0;
1482bc0b47dSJoe Wallwork   if (numProcs > 1) {
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetPointSF(dm, &sf));
1505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFSetUp(sf));
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
1535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(numProcs, &numVerInterfaces));
1542bc0b47dSJoe Wallwork 
155c51fff1bSJoe Wallwork     /* Count number of roots associated with each leaf */
1562bc0b47dSJoe Wallwork     for (r = 0; r < niranks; ++r) {
1572bc0b47dSJoe Wallwork       for (i=ioffset[r], count=0; i<ioffset[r+1]; ++i) {
1582bc0b47dSJoe Wallwork         if (irootloc[i] >= vStart && irootloc[i] < vEnd) count++;
1592bc0b47dSJoe Wallwork       }
1602bc0b47dSJoe Wallwork       numVerInterfaces[iranks[r]] += count;
1612bc0b47dSJoe Wallwork     }
162c51fff1bSJoe Wallwork 
163c51fff1bSJoe Wallwork     /* Count number of leaves associated with each root */
1642bc0b47dSJoe Wallwork     for (r = 0; r < nrranks; ++r) {
1652bc0b47dSJoe Wallwork       for (i=roffset[r], count=0; i<roffset[r+1]; ++i) {
1662bc0b47dSJoe Wallwork         if (rmine[i] >= vStart && rmine[i] < vEnd) count++;
1672bc0b47dSJoe Wallwork       }
1682bc0b47dSJoe Wallwork       numVerInterfaces[rranks[r]] += count;
1692bc0b47dSJoe Wallwork     }
170c51fff1bSJoe Wallwork 
171c51fff1bSJoe Wallwork     /* Count global number of ranks */
1722bc0b47dSJoe Wallwork     for (p = 0; p < numProcs; ++p) {
1732bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) numNgbRanks++;
1742bc0b47dSJoe Wallwork     }
175c51fff1bSJoe Wallwork 
176c51fff1bSJoe Wallwork     /* Provide numbers of vertex interfaces */
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc2(numNgbRanks, &ngbRanks, numNgbRanks, &verNgbRank));
1782bc0b47dSJoe Wallwork     for (p = 0, n = 0; p < numProcs; ++p) {
1792bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) {
1802bc0b47dSJoe Wallwork         ngbRanks[n] = p;
1812bc0b47dSJoe Wallwork         verNgbRank[n] = numVerInterfaces[p];
1822bc0b47dSJoe Wallwork         n++;
1832bc0b47dSJoe Wallwork       }
1842bc0b47dSJoe Wallwork     }
1852bc0b47dSJoe Wallwork     numVerNgbRanksTotal = 0;
1862bc0b47dSJoe Wallwork     for (p = 0; p < numNgbRanks; ++p) numVerNgbRanksTotal += verNgbRank[p];
1872bc0b47dSJoe Wallwork 
1882bc0b47dSJoe Wallwork     /* For each neighbor, fill in interface arrays */
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc3(numVerNgbRanksTotal, &interfaces_lv, numVerNgbRanksTotal, &interfaces_gv,  numNgbRanks+1, &intOffset));
1902bc0b47dSJoe Wallwork     intOffset[0] = 0;
1912bc0b47dSJoe Wallwork     for (p = 0, r = 0, i = 0; p < numNgbRanks; ++p) {
1922bc0b47dSJoe Wallwork       intOffset[p+1] = intOffset[p];
193c51fff1bSJoe Wallwork 
194c51fff1bSJoe Wallwork       /* Leaf case */
1952bc0b47dSJoe Wallwork       if (iranks && iranks[i] == ngbRanks[p]) {
1962bc0b47dSJoe Wallwork 
1972bc0b47dSJoe Wallwork         /* Add the right slice of irootloc at the right place */
1982bc0b47dSJoe Wallwork         sliceSize = ioffset[i+1]-ioffset[i];
1992bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
2005f80ce2aSJacob 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]);
2012bc0b47dSJoe Wallwork           v = irootloc[ioffset[i]+j];
2022bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
2035f80ce2aSJacob 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);
2042bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p+1]+count] = v-vStart;
2052bc0b47dSJoe Wallwork             count++;
2062bc0b47dSJoe Wallwork           }
2072bc0b47dSJoe Wallwork         }
2082bc0b47dSJoe Wallwork         intOffset[p+1] += count;
2092bc0b47dSJoe Wallwork         i++;
2102bc0b47dSJoe Wallwork       }
211c51fff1bSJoe Wallwork 
212c51fff1bSJoe Wallwork       /* Root case */
2132bc0b47dSJoe Wallwork       if (rranks && rranks[r] == ngbRanks[p]) {
2142bc0b47dSJoe Wallwork 
2152bc0b47dSJoe Wallwork         /* Add the right slice of rmine at the right place */
2162bc0b47dSJoe Wallwork         sliceSize = roffset[r+1]-roffset[r];
2172bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
2185f80ce2aSJacob 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]);
2192bc0b47dSJoe Wallwork           v = rmine[roffset[r]+j];
2202bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
2215f80ce2aSJacob 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);
2222bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p+1]+count] = v-vStart;
2232bc0b47dSJoe Wallwork             count++;
2242bc0b47dSJoe Wallwork           }
2252bc0b47dSJoe Wallwork         }
2262bc0b47dSJoe Wallwork         intOffset[p+1] += count;
2272bc0b47dSJoe Wallwork         r++;
2282bc0b47dSJoe Wallwork       }
229c51fff1bSJoe Wallwork 
230c51fff1bSJoe Wallwork       /* Check validity of offsets */
2315f80ce2aSJacob 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]);
2322bc0b47dSJoe Wallwork     }
2335f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetVertexNumbering(udm, &globalVertexNum));
2345f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(globalVertexNum, &gV));
2352bc0b47dSJoe Wallwork     for (i = 0; i < numVerNgbRanksTotal; ++i) {
2362bc0b47dSJoe Wallwork       v = gV[interfaces_lv[i]];
2372bc0b47dSJoe Wallwork       interfaces_gv[i] = v < 0 ? -v-1 : v;
2382bc0b47dSJoe Wallwork       interfaces_lv[i] += 1;
2392bc0b47dSJoe Wallwork       interfaces_gv[i] += 1;
2402bc0b47dSJoe Wallwork     }
2415f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(globalVertexNum, &gV));
2425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(numVerInterfaces));
2432bc0b47dSJoe Wallwork   }
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&udm));
2452bc0b47dSJoe Wallwork 
2462bc0b47dSJoe Wallwork   /* Send the data to ParMmg and remesh */
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricNoInsertion(dm, &noInsert));
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricNoSwapping(dm, &noSwap));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricNoMovement(dm, &noMove));
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricNoSurf(dm, &noSurf));
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricGetVerbosity(dm, &verbosity));
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricGetNumIterations(dm, &numIter));
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricGetGradationFactor(dm, &gradationFactor));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber));
2555f80ce2aSJacob Faibussowitsch   CHKERRMMG_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));
2565f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_meshSize(parmesh, numVertices, numCells, 0, numFaceTags, 0, 0));
2575f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_APImode, PMMG_APIDISTRIB_nodes));
2585f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noinsert, noInsert));
2595f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noswap, noSwap));
2605f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nomove, noMove));
2615f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nosurf, noSurf));
2625f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_verbose, verbosity));
2635f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1));
2645f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter(parmesh, PMMG_IPARAM_niter, numIter));
2655f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hgrad, gradationFactor));
2665f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hausd, hausdorffNumber));
2675f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_vertices(parmesh, vertices, verTags));
2685f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_tetrahedra(parmesh, cells, cellTags));
2695f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_triangles(parmesh, bdFaces, faceTags));
2705f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_metSize(parmesh, MMG5_Vertex, numVertices, MMG5_Tensor));
2715f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_tensorMets(parmesh, metric));
2725f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_numberOfNodeCommunicators(parmesh, numNgbRanks));
2732bc0b47dSJoe Wallwork   for (c = 0; c < numNgbRanks; ++c) {
2745f80ce2aSJacob Faibussowitsch     CHKERRMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicatorSize, parmesh, c, ngbRanks[c], intOffset[c+1]-intOffset[c]);
2755f80ce2aSJacob Faibussowitsch     CHKERRMMG_NONSTANDARD(PMMG_Set_ithNodeCommunicator_nodes, parmesh, c, &interfaces_lv[intOffset[c]], &interfaces_gv[intOffset[c]], 1);
2762bc0b47dSJoe Wallwork   }
2775f80ce2aSJacob Faibussowitsch   CHKERRMMG(PMMG_parmmglib_distributed, parmesh);
2785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cells));
2795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(metric, vertices));
2805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(bdFaces, faceTags));
2815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(verTags, cellTags));
2822bc0b47dSJoe Wallwork   if (numProcs > 1) {
2835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(ngbRanks, verNgbRank));
2845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree3(interfaces_lv, interfaces_gv, intOffset));
2852bc0b47dSJoe Wallwork   }
2862bc0b47dSJoe Wallwork 
2872f583692SJoe Wallwork   /* Retrieve mesh from Mmg */
2882bc0b47dSJoe Wallwork   numCornersNew = 4;
2895f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Get_meshSize, parmesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0);
2905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(dim*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3((dim+1)*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
2925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(dim*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
2935f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Get_vertices, parmesh, verticesNew, verTagsNew, corners, requiredVer);
2945f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Get_tetrahedra, parmesh, cellsNew, cellTagsNew, requiredCells);
2955f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Get_triangles, parmesh, facesNew, faceTagsNew, requiredFaces);
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(numVerticesNew, &owners, numVerticesNew, &gv_new));
2975f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Set_iparameter, parmesh, PMMG_IPARAM_globalNum, 1);
2985f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Get_verticesGloNum, parmesh, gv_new, owners);
2992bc0b47dSJoe Wallwork   for (i = 0; i < dim*numFacesNew; ++i) facesNew[i] -= 1;
3000a7a67b9SPierre Jolivet   for (i = 0; i < (dim+1)*numCellsNew; ++i) cellsNew[i] = gv_new[cellsNew[i]-1]-1;
3010a7a67b9SPierre Jolivet   for (i = 0, numVerticesNewLoc = 0; i < numVerticesNew; ++i) {
3022bc0b47dSJoe Wallwork     if (owners[i] == rank) numVerticesNewLoc++;
3032bc0b47dSJoe Wallwork   }
3045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(numVerticesNewLoc*dim, &verticesNewLoc, numVerticesNew, &verticesNewSorted));
3052bc0b47dSJoe Wallwork   for (i = 0, c = 0; i < numVerticesNew; i++) {
3062bc0b47dSJoe Wallwork     if (owners[i] == rank) {
3072bc0b47dSJoe Wallwork       for (j=0; j<dim; ++j) verticesNewLoc[dim*c+j] = verticesNew[dim*i+j];
3082bc0b47dSJoe Wallwork       c++;
3092bc0b47dSJoe Wallwork     }
3102bc0b47dSJoe Wallwork   }
3112f583692SJoe Wallwork 
3122f583692SJoe Wallwork   /* Reorder for consistency with DMPlex */
3135f80ce2aSJacob Faibussowitsch   for (i = 0; i < numCellsNew; ++i) CHKERRQ(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4*i]));
3142f583692SJoe Wallwork 
3152f583692SJoe Wallwork   /* Create new plex */
3165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNewLoc, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNewLoc, NULL, &verticesNewSorted, dmNew));
3175f80ce2aSJacob Faibussowitsch   CHKERRMMG_NONSTANDARD(PMMG_Free_all, PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_end);
3185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(verticesNew, verTagsNew, corners, requiredVer));
319c1dc6da0SJoe Wallwork 
320c1dc6da0SJoe Wallwork   /* Get adapted mesh information */
3215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
3225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
3235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
3242bc0b47dSJoe Wallwork 
3252bc0b47dSJoe Wallwork   /* Rebuild boundary label */
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName));
3275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew));
3282bc0b47dSJoe Wallwork   for (i = 0; i < numFacesNew; i++) {
329c1dc6da0SJoe Wallwork     PetscBool       hasTag = PETSC_FALSE;
3302bc0b47dSJoe Wallwork     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
3312bc0b47dSJoe Wallwork     const PetscInt *coveredPoints = NULL;
3322bc0b47dSJoe Wallwork 
3332bc0b47dSJoe Wallwork     for (j = 0; j < dim; ++j) {
3342bc0b47dSJoe Wallwork       lv = facesNew[i*dim+j];
3352bc0b47dSJoe Wallwork       gv = gv_new[lv]-1;
3365f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFindInt(gv, numVerticesNew, verticesNewSorted, &lv));
3372bc0b47dSJoe Wallwork       facePoints[j] = lv+vStart;
3382bc0b47dSJoe Wallwork     }
3395f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
3402bc0b47dSJoe Wallwork     for (j = 0; j < numCoveredPoints; ++j) {
3412bc0b47dSJoe Wallwork       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
3422bc0b47dSJoe Wallwork         numFaces++;
3432bc0b47dSJoe Wallwork         f = j;
3442bc0b47dSJoe Wallwork       }
3452bc0b47dSJoe Wallwork     }
3465f80ce2aSJacob Faibussowitsch     PetscCheck(numFaces == 1,comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces);
3475f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelHasStratum(bdLabel, faceTagsNew[i], &hasTag));
3485f80ce2aSJacob Faibussowitsch     if (hasTag) CHKERRQ(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]));
3495f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
3502bc0b47dSJoe Wallwork   }
3515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces));
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(owners, gv_new));
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(verticesNewLoc, verticesNewSorted));
3545f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(DMLabelDestroy(&bdLabel));
3552bc0b47dSJoe Wallwork 
3569fe9e680SJoe Wallwork   /* Rebuild cell labels */
3575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName));
3585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew));
3595f80ce2aSJacob Faibussowitsch   for (c = cStart; c < cEnd; ++c) CHKERRQ(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c-cStart]));
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(cellsNew, cellTagsNew, requiredCells));
361c1dc6da0SJoe Wallwork 
3622bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
3632bc0b47dSJoe Wallwork }
364