xref: /petsc/src/dm/impls/plex/adaptors/parmmg/parmmgadapt.c (revision ae8b063e4f7e8ac6d2de6909aa6d539b126b38d9)
12bc0b47dSJoe Wallwork #include <petsc/private/dmpleximpl.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;
17*ae8b063eSJoe 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;
2893520041SJoe Wallwork   PetscBool          flg = PETSC_FALSE, noInsert, noSwap, noMove, isotropic, uniform;
292bc0b47dSJoe Wallwork   const PetscMPIInt *iranks, *rranks;
302bc0b47dSJoe Wallwork   PetscMPIInt        numProcs, rank;
312bc0b47dSJoe Wallwork   PMMG_pParMesh      parmesh = NULL;
322bc0b47dSJoe Wallwork   PetscErrorCode     ierr;
332bc0b47dSJoe Wallwork 
342bc0b47dSJoe Wallwork   PetscFunctionBegin;
352bc0b47dSJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
362bc0b47dSJoe Wallwork   ierr = MPI_Comm_size(comm, &numProcs);CHKERRMPI(ierr);
372bc0b47dSJoe Wallwork   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
382bc0b47dSJoe Wallwork   if (bdLabel) {
392bc0b47dSJoe Wallwork     ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr);
402bc0b47dSJoe Wallwork     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
412c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
422bc0b47dSJoe Wallwork   }
439fe9e680SJoe Wallwork   if (rgLabel) {
449fe9e680SJoe Wallwork     ierr = PetscObjectGetName((PetscObject) rgLabel, &rgLabelName);CHKERRQ(ierr);
459fe9e680SJoe Wallwork     ierr = PetscStrcmp(rgLabelName, rgName, &flg);CHKERRQ(ierr);
462c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
479fe9e680SJoe Wallwork   }
482bc0b47dSJoe Wallwork 
492bc0b47dSJoe Wallwork   /* Get mesh information */
502bc0b47dSJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
512c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dim != 3,comm, PETSC_ERR_ARG_OUTOFRANGE, "ParMmg only works in 3D.");
522bc0b47dSJoe Wallwork   Neq  = (dim*(dim+1))/2;
532bc0b47dSJoe Wallwork   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
54c1dc6da0SJoe Wallwork   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
552bc0b47dSJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
562bc0b47dSJoe Wallwork   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
572bc0b47dSJoe Wallwork   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
582bc0b47dSJoe Wallwork   numCells    = cEnd - cStart;
592bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
602bc0b47dSJoe Wallwork 
612bc0b47dSJoe Wallwork   /* Get cell offsets */
622bc0b47dSJoe Wallwork   ierr = PetscMalloc1(numCells*maxConeSize, &cells);CHKERRQ(ierr);
632bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
642bc0b47dSJoe Wallwork     const PetscInt *cone;
652bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
662bc0b47dSJoe Wallwork 
672bc0b47dSJoe Wallwork     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
682bc0b47dSJoe Wallwork     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
692bc0b47dSJoe Wallwork     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart+1;
702bc0b47dSJoe Wallwork   }
712bc0b47dSJoe Wallwork 
722bc0b47dSJoe Wallwork   /* Get vertex coordinate array */
732bc0b47dSJoe Wallwork   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
742bc0b47dSJoe Wallwork   ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
752bc0b47dSJoe Wallwork   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
762bc0b47dSJoe Wallwork   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
772bc0b47dSJoe Wallwork   ierr = PetscMalloc2(numVertices*Neq, &metric, dim*numVertices, &vertices);CHKERRQ(ierr);
782bc0b47dSJoe Wallwork   for (v = 0; v < vEnd-vStart; ++v) {
792bc0b47dSJoe Wallwork     ierr = PetscSectionGetOffset(coordSection, v+vStart, &off);CHKERRQ(ierr);
802bc0b47dSJoe Wallwork     for (i = 0; i < dim; ++i) vertices[dim*v+i] = PetscRealPart(coords[off+i]);
812bc0b47dSJoe Wallwork   }
822bc0b47dSJoe Wallwork   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
832bc0b47dSJoe Wallwork 
84c1dc6da0SJoe Wallwork   /* Get face tags */
85c1dc6da0SJoe Wallwork   if (!bdLabel) {
86c1dc6da0SJoe Wallwork     flg = PETSC_TRUE;
87c1dc6da0SJoe Wallwork     ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel);CHKERRQ(ierr);
88c1dc6da0SJoe Wallwork     ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabel);CHKERRQ(ierr);
89c1dc6da0SJoe Wallwork   }
90c1dc6da0SJoe Wallwork   ierr = DMLabelGetBounds(bdLabel, &pStart, &pEnd);CHKERRQ(ierr);
91c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
92c1dc6da0SJoe Wallwork     PetscBool hasPoint;
93c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
942bc0b47dSJoe Wallwork 
95c1dc6da0SJoe Wallwork     ierr = DMLabelHasPoint(bdLabel, f, &hasPoint);CHKERRQ(ierr);
96c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
97c1dc6da0SJoe Wallwork     numFaceTags++;
98c1dc6da0SJoe Wallwork 
99c1dc6da0SJoe Wallwork     ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1002bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1012bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1022bc0b47dSJoe Wallwork     }
103c1dc6da0SJoe Wallwork     ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1042bc0b47dSJoe Wallwork   }
105c1dc6da0SJoe Wallwork   ierr = PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags);CHKERRQ(ierr);
106c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
107c1dc6da0SJoe Wallwork     PetscBool hasPoint;
108c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1092bc0b47dSJoe Wallwork 
110c1dc6da0SJoe Wallwork     ierr = DMLabelHasPoint(bdLabel, f, &hasPoint);CHKERRQ(ierr);
111c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
112c1dc6da0SJoe Wallwork 
113c1dc6da0SJoe Wallwork     ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1142bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1152bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart+1;
1162bc0b47dSJoe Wallwork     }
117c1dc6da0SJoe Wallwork     ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
118c1dc6da0SJoe Wallwork     ierr = DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]);CHKERRQ(ierr);
1192bc0b47dSJoe Wallwork   }
1202bc0b47dSJoe Wallwork 
1219fe9e680SJoe Wallwork   /* Get cell tags */
1229fe9e680SJoe Wallwork   ierr = PetscCalloc2(numVertices, &verTags, numCells, &cellTags);CHKERRQ(ierr);
1239fe9e680SJoe Wallwork   if (rgLabel) {
1249fe9e680SJoe Wallwork     for (c = cStart; c < cEnd; ++c) { ierr = DMLabelGetValue(rgLabel, c, &cellTags[c]);CHKERRQ(ierr); }
1259fe9e680SJoe Wallwork   }
1269fe9e680SJoe Wallwork 
1272bc0b47dSJoe Wallwork   /* Get metric */
1282bc0b47dSJoe Wallwork   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
1292bc0b47dSJoe Wallwork   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
13093520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
13193520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
1322bc0b47dSJoe Wallwork   for (v = 0; v < (vEnd-vStart); ++v) {
1332bc0b47dSJoe Wallwork     for (i = 0, k = 0; i < dim; ++i) {
1340a7a67b9SPierre Jolivet       for (j = i; j < dim; ++j, ++k) {
13593520041SJoe Wallwork         if (isotropic) {
13693520041SJoe Wallwork           if (i == j) {
13793520041SJoe Wallwork             if (uniform) metric[Neq*v+k] = PetscRealPart(met[0]);
13893520041SJoe Wallwork             else metric[Neq*v+k] = PetscRealPart(met[v]);
13993520041SJoe Wallwork           } else metric[Neq*v+k] = 0.0;
1400a7a67b9SPierre Jolivet         } else metric[Neq*v+k] = PetscRealPart(met[dim*dim*v+dim*i+j]);
1412bc0b47dSJoe Wallwork       }
1422bc0b47dSJoe Wallwork     }
1432bc0b47dSJoe Wallwork   }
1442bc0b47dSJoe Wallwork   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
1452bc0b47dSJoe Wallwork 
1462bc0b47dSJoe Wallwork   /* Build ParMMG communicators: the list of vertices between two partitions  */
1472bc0b47dSJoe Wallwork   niranks = nrranks = 0;
1482bc0b47dSJoe Wallwork   numNgbRanks = 0;
1492bc0b47dSJoe Wallwork   if (numProcs > 1) {
1502bc0b47dSJoe Wallwork     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
1512bc0b47dSJoe Wallwork     ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1522bc0b47dSJoe Wallwork     ierr = PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc);CHKERRQ(ierr);
1532bc0b47dSJoe Wallwork     ierr = PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote);CHKERRQ(ierr);
1542bc0b47dSJoe Wallwork     ierr = PetscCalloc1(numProcs, &numVerInterfaces);CHKERRQ(ierr);
1552bc0b47dSJoe Wallwork 
156c51fff1bSJoe Wallwork     /* Count number of roots associated with each leaf */
1572bc0b47dSJoe Wallwork     for (r = 0; r < niranks; ++r) {
1582bc0b47dSJoe Wallwork       for (i=ioffset[r], count=0; i<ioffset[r+1]; ++i) {
1592bc0b47dSJoe Wallwork         if (irootloc[i] >= vStart && irootloc[i] < vEnd) count++;
1602bc0b47dSJoe Wallwork       }
1612bc0b47dSJoe Wallwork       numVerInterfaces[iranks[r]] += count;
1622bc0b47dSJoe Wallwork     }
163c51fff1bSJoe Wallwork 
164c51fff1bSJoe Wallwork     /* Count number of leaves associated with each root */
1652bc0b47dSJoe Wallwork     for (r = 0; r < nrranks; ++r) {
1662bc0b47dSJoe Wallwork       for (i=roffset[r], count=0; i<roffset[r+1]; ++i) {
1672bc0b47dSJoe Wallwork         if (rmine[i] >= vStart && rmine[i] < vEnd) count++;
1682bc0b47dSJoe Wallwork       }
1692bc0b47dSJoe Wallwork       numVerInterfaces[rranks[r]] += count;
1702bc0b47dSJoe Wallwork     }
171c51fff1bSJoe Wallwork 
172c51fff1bSJoe Wallwork     /* Count global number of ranks */
1732bc0b47dSJoe Wallwork     for (p = 0; p < numProcs; ++p) {
1742bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) numNgbRanks++;
1752bc0b47dSJoe Wallwork     }
176c51fff1bSJoe Wallwork 
177c51fff1bSJoe Wallwork     /* Provide numbers of vertex interfaces */
1782bc0b47dSJoe Wallwork     ierr = PetscMalloc2(numNgbRanks, &ngbRanks, numNgbRanks, &verNgbRank);CHKERRQ(ierr);
1792bc0b47dSJoe Wallwork     for (p = 0, n = 0; p < numProcs; ++p) {
1802bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) {
1812bc0b47dSJoe Wallwork         ngbRanks[n] = p;
1822bc0b47dSJoe Wallwork         verNgbRank[n] = numVerInterfaces[p];
1832bc0b47dSJoe Wallwork         n++;
1842bc0b47dSJoe Wallwork       }
1852bc0b47dSJoe Wallwork     }
1862bc0b47dSJoe Wallwork     numVerNgbRanksTotal = 0;
1872bc0b47dSJoe Wallwork     for (p = 0; p < numNgbRanks; ++p) numVerNgbRanksTotal += verNgbRank[p];
1882bc0b47dSJoe Wallwork 
1892bc0b47dSJoe Wallwork     /* For each neighbor, fill in interface arrays */
1902bc0b47dSJoe Wallwork     ierr = PetscMalloc3(numVerNgbRanksTotal, &interfaces_lv, numVerNgbRanksTotal, &interfaces_gv,  numNgbRanks+1, &intOffset);CHKERRQ(ierr);
1912bc0b47dSJoe Wallwork     intOffset[0] = 0;
1922bc0b47dSJoe Wallwork     for (p = 0, r = 0, i = 0; p < numNgbRanks; ++p) {
1932bc0b47dSJoe Wallwork       intOffset[p+1] = intOffset[p];
194c51fff1bSJoe Wallwork 
195c51fff1bSJoe Wallwork       /* Leaf case */
1962bc0b47dSJoe Wallwork       if (iranks && iranks[i] == ngbRanks[p]) {
1972bc0b47dSJoe Wallwork 
1982bc0b47dSJoe Wallwork         /* Add the right slice of irootloc at the right place */
1992bc0b47dSJoe Wallwork         sliceSize = ioffset[i+1]-ioffset[i];
2002bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
2012c71b3e2SJacob Faibussowitsch           PetscCheckFalse(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]);
2022bc0b47dSJoe Wallwork           v = irootloc[ioffset[i]+j];
2032bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
2042c71b3e2SJacob Faibussowitsch             PetscCheckFalse(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);
2052bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p+1]+count] = v-vStart;
2062bc0b47dSJoe Wallwork             count++;
2072bc0b47dSJoe Wallwork           }
2082bc0b47dSJoe Wallwork         }
2092bc0b47dSJoe Wallwork         intOffset[p+1] += count;
2102bc0b47dSJoe Wallwork         i++;
2112bc0b47dSJoe Wallwork       }
212c51fff1bSJoe Wallwork 
213c51fff1bSJoe Wallwork       /* Root case */
2142bc0b47dSJoe Wallwork       if (rranks && rranks[r] == ngbRanks[p]) {
2152bc0b47dSJoe Wallwork 
2162bc0b47dSJoe Wallwork         /* Add the right slice of rmine at the right place */
2172bc0b47dSJoe Wallwork         sliceSize = roffset[r+1]-roffset[r];
2182bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
2192c71b3e2SJacob Faibussowitsch           PetscCheckFalse(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]);
2202bc0b47dSJoe Wallwork           v = rmine[roffset[r]+j];
2212bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
2222c71b3e2SJacob Faibussowitsch             PetscCheckFalse(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);
2232bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p+1]+count] = v-vStart;
2242bc0b47dSJoe Wallwork             count++;
2252bc0b47dSJoe Wallwork           }
2262bc0b47dSJoe Wallwork         }
2272bc0b47dSJoe Wallwork         intOffset[p+1] += count;
2282bc0b47dSJoe Wallwork         r++;
2292bc0b47dSJoe Wallwork       }
230c51fff1bSJoe Wallwork 
231c51fff1bSJoe Wallwork       /* Check validity of offsets */
2322c71b3e2SJacob Faibussowitsch       PetscCheckFalse(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]);
2332bc0b47dSJoe Wallwork     }
2342bc0b47dSJoe Wallwork     ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
2352bc0b47dSJoe Wallwork     ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
2362bc0b47dSJoe Wallwork     for (i = 0; i < numVerNgbRanksTotal; ++i) {
2372bc0b47dSJoe Wallwork       v = gV[interfaces_lv[i]];
2382bc0b47dSJoe Wallwork       interfaces_gv[i] = v < 0 ? -v-1 : v;
2392bc0b47dSJoe Wallwork       interfaces_lv[i] += 1;
2402bc0b47dSJoe Wallwork       interfaces_gv[i] += 1;
2412bc0b47dSJoe Wallwork     }
2422bc0b47dSJoe Wallwork     ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
2432bc0b47dSJoe Wallwork     ierr = PetscFree(numVerInterfaces);CHKERRQ(ierr);
2442bc0b47dSJoe Wallwork   }
2452bc0b47dSJoe Wallwork   ierr = DMDestroy(&udm);CHKERRQ(ierr);
2462bc0b47dSJoe Wallwork 
2472bc0b47dSJoe Wallwork   /* Send the data to ParMmg and remesh */
2482bc0b47dSJoe Wallwork   ierr = DMPlexMetricNoInsertion(dm, &noInsert);CHKERRQ(ierr);
2492bc0b47dSJoe Wallwork   ierr = DMPlexMetricNoSwapping(dm, &noSwap);CHKERRQ(ierr);
2502bc0b47dSJoe Wallwork   ierr = DMPlexMetricNoMovement(dm, &noMove);CHKERRQ(ierr);
2512bc0b47dSJoe Wallwork   ierr = DMPlexMetricGetVerbosity(dm, &verbosity);CHKERRQ(ierr);
2522bc0b47dSJoe Wallwork   ierr = DMPlexMetricGetNumIterations(dm, &numIter);CHKERRQ(ierr);
2532bc0b47dSJoe Wallwork   ierr = DMPlexMetricGetGradationFactor(dm, &gradationFactor);CHKERRQ(ierr);
254*ae8b063eSJoe Wallwork   ierr = DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber);CHKERRQ(ierr);
2550a7a67b9SPierre Jolivet   ierr = 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);
256c1dc6da0SJoe Wallwork   ierr = PMMG_Set_meshSize(parmesh, numVertices, numCells, 0, numFaceTags, 0, 0);
2572bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_APImode, PMMG_APIDISTRIB_nodes);
2582bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noinsert, noInsert);
2592bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noswap, noSwap);
2602bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nomove, noMove);
2612bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_verbose, verbosity);
2622bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1);
2632bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_niter, numIter);
2642bc0b47dSJoe Wallwork   ierr = PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hgrad, gradationFactor);
265*ae8b063eSJoe Wallwork   ierr = PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hausd, hausdorffNumber);
2662bc0b47dSJoe Wallwork   ierr = PMMG_Set_vertices(parmesh, vertices, verTags);
267cf6859d0SJoe Wallwork   ierr = PMMG_Set_tetrahedra(parmesh, cells, cellTags);
268c1dc6da0SJoe Wallwork   ierr = PMMG_Set_triangles(parmesh, bdFaces, faceTags);
2692bc0b47dSJoe Wallwork   ierr = PMMG_Set_metSize(parmesh, MMG5_Vertex, numVertices, MMG5_Tensor);
270cf6859d0SJoe Wallwork   ierr = PMMG_Set_tensorMets(parmesh, metric);
2712bc0b47dSJoe Wallwork   ierr = PMMG_Set_numberOfNodeCommunicators(parmesh, numNgbRanks);
2722bc0b47dSJoe Wallwork   for (c = 0; c < numNgbRanks; ++c) {
2732bc0b47dSJoe Wallwork     ierr = PMMG_Set_ithNodeCommunicatorSize(parmesh, c, ngbRanks[c], intOffset[c+1]-intOffset[c]);
2742bc0b47dSJoe Wallwork     ierr = PMMG_Set_ithNodeCommunicator_nodes(parmesh, c, &interfaces_lv[intOffset[c]], &interfaces_gv[intOffset[c]], 1);
2752bc0b47dSJoe Wallwork   }
2762bc0b47dSJoe Wallwork   ierr = PMMG_parmmglib_distributed(parmesh);
2772bc0b47dSJoe Wallwork   ierr = PetscFree(cells);CHKERRQ(ierr);
2782bc0b47dSJoe Wallwork   ierr = PetscFree2(metric, vertices);CHKERRQ(ierr);
279c1dc6da0SJoe Wallwork   ierr = PetscFree2(bdFaces, faceTags);CHKERRQ(ierr);
2802bc0b47dSJoe Wallwork   ierr = PetscFree2(verTags, cellTags);CHKERRQ(ierr);
2812bc0b47dSJoe Wallwork   if (numProcs > 1) {
2822bc0b47dSJoe Wallwork     ierr = PetscFree2(ngbRanks, verNgbRank);CHKERRQ(ierr);
2832bc0b47dSJoe Wallwork     ierr = PetscFree3(interfaces_lv, interfaces_gv, intOffset);CHKERRQ(ierr);
2842bc0b47dSJoe Wallwork   }
2852bc0b47dSJoe Wallwork 
2862f583692SJoe Wallwork   /* Retrieve mesh from Mmg */
2872bc0b47dSJoe Wallwork   numCornersNew = 4;
2882bc0b47dSJoe Wallwork   ierr = PMMG_Get_meshSize(parmesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0);
2892bc0b47dSJoe Wallwork   ierr = PetscMalloc4(dim*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer);CHKERRQ(ierr);
2902bc0b47dSJoe Wallwork   ierr = PetscMalloc3((dim+1)*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells);CHKERRQ(ierr);
2912bc0b47dSJoe Wallwork   ierr = PetscMalloc4(dim*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces);CHKERRQ(ierr);
2922bc0b47dSJoe Wallwork   ierr = PMMG_Get_vertices(parmesh, verticesNew, verTagsNew, corners, requiredVer);
2932bc0b47dSJoe Wallwork   ierr = PMMG_Get_tetrahedra(parmesh, cellsNew, cellTagsNew, requiredCells);
2942bc0b47dSJoe Wallwork   ierr = PMMG_Get_triangles(parmesh, facesNew, faceTagsNew, requiredFaces);
2952bc0b47dSJoe Wallwork   ierr = PetscMalloc2(numVerticesNew, &owners, numVerticesNew, &gv_new);
2962bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1);
2972bc0b47dSJoe Wallwork   ierr = PMMG_Get_verticesGloNum(parmesh, gv_new, owners);
2982bc0b47dSJoe Wallwork   for (i = 0; i < dim*numFacesNew; ++i) facesNew[i] -= 1;
2990a7a67b9SPierre Jolivet   for (i = 0; i < (dim+1)*numCellsNew; ++i) cellsNew[i] = gv_new[cellsNew[i]-1]-1;
3000a7a67b9SPierre Jolivet   for (i = 0, numVerticesNewLoc = 0; i < numVerticesNew; ++i) {
3012bc0b47dSJoe Wallwork     if (owners[i] == rank) numVerticesNewLoc++;
3022bc0b47dSJoe Wallwork   }
3032bc0b47dSJoe Wallwork   ierr = PetscMalloc2(numVerticesNewLoc*dim, &verticesNewLoc, numVerticesNew, &verticesNewSorted);CHKERRQ(ierr);
3042bc0b47dSJoe Wallwork   for (i = 0, c = 0; i < numVerticesNew; i++) {
3052bc0b47dSJoe Wallwork     if (owners[i] == rank) {
3062bc0b47dSJoe Wallwork       for (j=0; j<dim; ++j) verticesNewLoc[dim*c+j] = verticesNew[dim*i+j];
3072bc0b47dSJoe Wallwork       c++;
3082bc0b47dSJoe Wallwork     }
3092bc0b47dSJoe Wallwork   }
3102f583692SJoe Wallwork 
3112f583692SJoe Wallwork   /* Reorder for consistency with DMPlex */
3122f583692SJoe Wallwork   for (i = 0; i < numCellsNew; ++i) { ierr = DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4*i]);CHKERRQ(ierr); }
3132f583692SJoe Wallwork 
3142f583692SJoe Wallwork   /* Create new plex */
3152bc0b47dSJoe Wallwork   ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNewLoc, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNewLoc, NULL, &verticesNewSorted, dmNew);CHKERRQ(ierr);
3162bc0b47dSJoe Wallwork   ierr = PMMG_Free_all(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_end);
317c1dc6da0SJoe Wallwork   ierr = PetscFree4(verticesNew, verTagsNew, corners, requiredVer);CHKERRQ(ierr);
318c1dc6da0SJoe Wallwork 
319c1dc6da0SJoe Wallwork   /* Get adapted mesh information */
320c1dc6da0SJoe Wallwork   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
321c1dc6da0SJoe Wallwork   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
322c1dc6da0SJoe Wallwork   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
3232bc0b47dSJoe Wallwork 
3242bc0b47dSJoe Wallwork   /* Rebuild boundary label */
3257881f211SJoe Wallwork   ierr = DMCreateLabel(*dmNew, flg ? bdName : bdLabelName);CHKERRQ(ierr);
3267881f211SJoe Wallwork   ierr = DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew);CHKERRQ(ierr);
3272bc0b47dSJoe Wallwork   for (i = 0; i < numFacesNew; i++) {
328c1dc6da0SJoe Wallwork     PetscBool       hasTag = PETSC_FALSE;
3292bc0b47dSJoe Wallwork     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
3302bc0b47dSJoe Wallwork     const PetscInt *coveredPoints = NULL;
3312bc0b47dSJoe Wallwork 
3322bc0b47dSJoe Wallwork     for (j = 0; j < dim; ++j) {
3332bc0b47dSJoe Wallwork       lv = facesNew[i*dim+j];
3342bc0b47dSJoe Wallwork       gv = gv_new[lv]-1;
3352bc0b47dSJoe Wallwork       ierr = PetscFindInt(gv, numVerticesNew, verticesNewSorted, &lv);CHKERRQ(ierr);
3362bc0b47dSJoe Wallwork       facePoints[j] = lv+vStart;
3372bc0b47dSJoe Wallwork     }
3382bc0b47dSJoe Wallwork     ierr = DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
3392bc0b47dSJoe Wallwork     for (j = 0; j < numCoveredPoints; ++j) {
3402bc0b47dSJoe Wallwork       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
3412bc0b47dSJoe Wallwork         numFaces++;
3422bc0b47dSJoe Wallwork         f = j;
3432bc0b47dSJoe Wallwork       }
3442bc0b47dSJoe Wallwork     }
3452c71b3e2SJacob Faibussowitsch     PetscCheckFalse(numFaces != 1,comm, PETSC_ERR_ARG_OUTOFRANGE, "%d vertices cannot define more than 1 facet (%d)", dim, numFaces);
346c1dc6da0SJoe Wallwork     ierr = DMLabelHasStratum(bdLabel, faceTagsNew[i], &hasTag);CHKERRQ(ierr);
347c1dc6da0SJoe Wallwork     if (hasTag) { ierr = DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]);CHKERRQ(ierr); }
3482bc0b47dSJoe Wallwork     ierr = DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
3492bc0b47dSJoe Wallwork   }
350c1dc6da0SJoe Wallwork   ierr = PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces);CHKERRQ(ierr);
351c1dc6da0SJoe Wallwork   ierr = PetscFree2(owners, gv_new);CHKERRQ(ierr);
352c1dc6da0SJoe Wallwork   ierr = PetscFree2(verticesNewLoc, verticesNewSorted);CHKERRQ(ierr);
353c1dc6da0SJoe Wallwork   if (flg) { ierr = DMLabelDestroy(&bdLabel);CHKERRQ(ierr); }
3542bc0b47dSJoe Wallwork 
3559fe9e680SJoe Wallwork   /* Rebuild cell labels */
3569fe9e680SJoe Wallwork   ierr = DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName);CHKERRQ(ierr);
3579fe9e680SJoe Wallwork   ierr = DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew);CHKERRQ(ierr);
3589fe9e680SJoe Wallwork   for (c = cStart; c < cEnd; ++c) { ierr = DMLabelSetValue(rgLabelNew, c, cellTagsNew[c-cStart]);CHKERRQ(ierr); }
3592bc0b47dSJoe Wallwork   ierr = PetscFree3(cellsNew, cellTagsNew, requiredCells);CHKERRQ(ierr);
360c1dc6da0SJoe Wallwork 
3612bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
3622bc0b47dSJoe Wallwork }
363