xref: /petsc/src/dm/impls/plex/adaptors/parmmg/parmmgadapt.c (revision 546078ac55d0f4e36dc638bf7d35f3433941479a)
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 {
62bc0b47dSJoe Wallwork   MPI_Comm           comm, tmpComm;
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;
17c1dc6da0SJoe Wallwork   PetscReal         *vertices, *metric, *verticesNew, *verticesNewLoc, gradationFactor;
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;
28c1dc6da0SJoe Wallwork   PetscBool          flg = PETSC_FALSE, noInsert, noSwap, noMove;
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   ierr = MPI_Comm_dup(comm, &tmpComm);CHKERRMPI(ierr);
392bc0b47dSJoe Wallwork   if (bdLabel) {
402bc0b47dSJoe Wallwork     ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr);
412bc0b47dSJoe Wallwork     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
422bc0b47dSJoe Wallwork     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
432bc0b47dSJoe Wallwork   }
449fe9e680SJoe Wallwork   if (rgLabel) {
459fe9e680SJoe Wallwork     ierr = PetscObjectGetName((PetscObject) rgLabel, &rgLabelName);CHKERRQ(ierr);
469fe9e680SJoe Wallwork     ierr = PetscStrcmp(rgLabelName, rgName, &flg);CHKERRQ(ierr);
479fe9e680SJoe Wallwork     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
489fe9e680SJoe Wallwork   }
492bc0b47dSJoe Wallwork 
502bc0b47dSJoe Wallwork   /* Get mesh information */
512bc0b47dSJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
52*546078acSJacob Faibussowitsch   if (dim != 3) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "ParMmg only works in 3D.");
532bc0b47dSJoe Wallwork   Neq  = (dim*(dim+1))/2;
542bc0b47dSJoe Wallwork   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
55c1dc6da0SJoe Wallwork   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
562bc0b47dSJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
572bc0b47dSJoe Wallwork   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
582bc0b47dSJoe Wallwork   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
592bc0b47dSJoe Wallwork   numCells    = cEnd - cStart;
602bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
612bc0b47dSJoe Wallwork 
622bc0b47dSJoe Wallwork   /* Get cell offsets */
632bc0b47dSJoe Wallwork   ierr = PetscMalloc1(numCells*maxConeSize, &cells);CHKERRQ(ierr);
642bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
652bc0b47dSJoe Wallwork     const PetscInt *cone;
662bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
672bc0b47dSJoe Wallwork 
682bc0b47dSJoe Wallwork     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
692bc0b47dSJoe Wallwork     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
702bc0b47dSJoe Wallwork     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart+1;
712bc0b47dSJoe Wallwork   }
722bc0b47dSJoe Wallwork 
732bc0b47dSJoe Wallwork   /* Get vertex coordinate array */
742bc0b47dSJoe Wallwork   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
752bc0b47dSJoe Wallwork   ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
762bc0b47dSJoe Wallwork   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
772bc0b47dSJoe Wallwork   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
782bc0b47dSJoe Wallwork   ierr = PetscMalloc2(numVertices*Neq, &metric, dim*numVertices, &vertices);CHKERRQ(ierr);
792bc0b47dSJoe Wallwork   for (v = 0; v < vEnd-vStart; ++v) {
802bc0b47dSJoe Wallwork     ierr = PetscSectionGetOffset(coordSection, v+vStart, &off);CHKERRQ(ierr);
812bc0b47dSJoe Wallwork     for (i = 0; i < dim; ++i) vertices[dim*v+i] = PetscRealPart(coords[off+i]);
822bc0b47dSJoe Wallwork   }
832bc0b47dSJoe Wallwork   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
842bc0b47dSJoe Wallwork 
85c1dc6da0SJoe Wallwork   /* Get face tags */
86c1dc6da0SJoe Wallwork   if (!bdLabel) {
87c1dc6da0SJoe Wallwork     flg = PETSC_TRUE;
88c1dc6da0SJoe Wallwork     ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel);CHKERRQ(ierr);
89c1dc6da0SJoe Wallwork     ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabel);CHKERRQ(ierr);
90c1dc6da0SJoe Wallwork   }
91c1dc6da0SJoe Wallwork   ierr = DMLabelGetBounds(bdLabel, &pStart, &pEnd);CHKERRQ(ierr);
92c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
93c1dc6da0SJoe Wallwork     PetscBool hasPoint;
94c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
952bc0b47dSJoe Wallwork 
96c1dc6da0SJoe Wallwork     ierr = DMLabelHasPoint(bdLabel, f, &hasPoint);CHKERRQ(ierr);
97c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
98c1dc6da0SJoe Wallwork     numFaceTags++;
99c1dc6da0SJoe Wallwork 
100c1dc6da0SJoe Wallwork     ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1012bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1022bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1032bc0b47dSJoe Wallwork     }
104c1dc6da0SJoe Wallwork     ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1052bc0b47dSJoe Wallwork   }
106c1dc6da0SJoe Wallwork   ierr = PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags);CHKERRQ(ierr);
107c1dc6da0SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
108c1dc6da0SJoe Wallwork     PetscBool hasPoint;
109c1dc6da0SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1102bc0b47dSJoe Wallwork 
111c1dc6da0SJoe Wallwork     ierr = DMLabelHasPoint(bdLabel, f, &hasPoint);CHKERRQ(ierr);
112c1dc6da0SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
113c1dc6da0SJoe Wallwork 
114c1dc6da0SJoe Wallwork     ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1152bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1162bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart+1;
1172bc0b47dSJoe Wallwork     }
118c1dc6da0SJoe Wallwork     ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
119c1dc6da0SJoe Wallwork     ierr = DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]);CHKERRQ(ierr);
1202bc0b47dSJoe Wallwork   }
1212bc0b47dSJoe Wallwork 
1229fe9e680SJoe Wallwork   /* Get cell tags */
1239fe9e680SJoe Wallwork   ierr = PetscCalloc2(numVertices, &verTags, numCells, &cellTags);CHKERRQ(ierr);
1249fe9e680SJoe Wallwork   if (rgLabel) {
1259fe9e680SJoe Wallwork     for (c = cStart; c < cEnd; ++c) { ierr = DMLabelGetValue(rgLabel, c, &cellTags[c]);CHKERRQ(ierr); }
1269fe9e680SJoe Wallwork   }
1279fe9e680SJoe Wallwork 
1282bc0b47dSJoe Wallwork   /* Get metric */
1292bc0b47dSJoe Wallwork   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
1302bc0b47dSJoe Wallwork   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
1312bc0b47dSJoe Wallwork   for (v = 0; v < (vEnd-vStart); ++v) {
1322bc0b47dSJoe Wallwork     for (i = 0, k = 0; i < dim; ++i) {
1332bc0b47dSJoe Wallwork       for (j = i; j < dim; ++j) {
1342bc0b47dSJoe Wallwork         metric[Neq*v+k] = PetscRealPart(met[dim*dim*v+dim*i+j]);
1352bc0b47dSJoe Wallwork         k++;
1362bc0b47dSJoe Wallwork       }
1372bc0b47dSJoe Wallwork     }
1382bc0b47dSJoe Wallwork   }
1392bc0b47dSJoe Wallwork   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
1402bc0b47dSJoe Wallwork 
1412bc0b47dSJoe Wallwork   /* Build ParMMG communicators: the list of vertices between two partitions  */
1422bc0b47dSJoe Wallwork   niranks = nrranks = 0;
1432bc0b47dSJoe Wallwork   numNgbRanks = 0;
1442bc0b47dSJoe Wallwork   if (numProcs > 1) {
1452bc0b47dSJoe Wallwork     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
1462bc0b47dSJoe Wallwork     ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1472bc0b47dSJoe Wallwork     ierr = PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc);CHKERRQ(ierr);
1482bc0b47dSJoe Wallwork     ierr = PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote);CHKERRQ(ierr);
1492bc0b47dSJoe Wallwork     ierr = PetscCalloc1(numProcs, &numVerInterfaces);CHKERRQ(ierr);
1502bc0b47dSJoe Wallwork 
1512bc0b47dSJoe Wallwork     /* Counting */
1522bc0b47dSJoe Wallwork     for (r = 0; r < niranks; ++r) {
1532bc0b47dSJoe Wallwork       for (i=ioffset[r], count=0; i<ioffset[r+1]; ++i) {
1542bc0b47dSJoe Wallwork         if (irootloc[i] >= vStart && irootloc[i] < vEnd) count++;
1552bc0b47dSJoe Wallwork       }
1562bc0b47dSJoe Wallwork       numVerInterfaces[iranks[r]] += count;
1572bc0b47dSJoe Wallwork     }
1582bc0b47dSJoe Wallwork     for (r = 0; r < nrranks; ++r) {
1592bc0b47dSJoe Wallwork       for (i=roffset[r], count=0; i<roffset[r+1]; ++i) {
1602bc0b47dSJoe Wallwork         if (rmine[i] >= vStart && rmine[i] < vEnd) count++;
1612bc0b47dSJoe Wallwork       }
1622bc0b47dSJoe Wallwork       numVerInterfaces[rranks[r]] += count;
1632bc0b47dSJoe Wallwork     }
1642bc0b47dSJoe Wallwork     for (p = 0; p < numProcs; ++p) {
1652bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) numNgbRanks++;
1662bc0b47dSJoe Wallwork     }
1672bc0b47dSJoe Wallwork     ierr = PetscMalloc2(numNgbRanks, &ngbRanks, numNgbRanks, &verNgbRank);CHKERRQ(ierr);
1682bc0b47dSJoe Wallwork     for (p = 0, n = 0; p < numProcs; ++p) {
1692bc0b47dSJoe Wallwork       if (numVerInterfaces[p]) {
1702bc0b47dSJoe Wallwork         ngbRanks[n] = p;
1712bc0b47dSJoe Wallwork         verNgbRank[n] = numVerInterfaces[p];
1722bc0b47dSJoe Wallwork         n++;
1732bc0b47dSJoe Wallwork       }
1742bc0b47dSJoe Wallwork     }
1752bc0b47dSJoe Wallwork     numVerNgbRanksTotal = 0;
1762bc0b47dSJoe Wallwork     for (p = 0; p < numNgbRanks; ++p) numVerNgbRanksTotal += verNgbRank[p];
1772bc0b47dSJoe Wallwork 
1782bc0b47dSJoe Wallwork     /* For each neighbor, fill in interface arrays */
1792bc0b47dSJoe Wallwork     ierr = PetscMalloc3(numVerNgbRanksTotal, &interfaces_lv, numVerNgbRanksTotal, &interfaces_gv,  numNgbRanks+1, &intOffset);CHKERRQ(ierr);
1802bc0b47dSJoe Wallwork     intOffset[0] = 0;
1812bc0b47dSJoe Wallwork     for (p = 0, r = 0, i = 0; p < numNgbRanks; ++p) {
1822bc0b47dSJoe Wallwork       intOffset[p+1] = intOffset[p];
1832bc0b47dSJoe Wallwork       if (iranks && iranks[i] == ngbRanks[p]) {
1842bc0b47dSJoe Wallwork 
1852bc0b47dSJoe Wallwork         /* Add the right slice of irootloc at the right place */
1862bc0b47dSJoe Wallwork         sliceSize = ioffset[i+1]-ioffset[i];
1872bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
1882bc0b47dSJoe Wallwork           if (ioffset[i]+j >= ioffset[niranks]) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
1892bc0b47dSJoe Wallwork           v = irootloc[ioffset[i]+j];
1902bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
1912bc0b47dSJoe Wallwork             if (intOffset[p+1]+count >= numVerNgbRanksTotal) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
1922bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p+1]+count] = v-vStart;
1932bc0b47dSJoe Wallwork             count++;
1942bc0b47dSJoe Wallwork           }
1952bc0b47dSJoe Wallwork         }
1962bc0b47dSJoe Wallwork         intOffset[p+1] += count;
1972bc0b47dSJoe Wallwork         i++;
1982bc0b47dSJoe Wallwork       }
1992bc0b47dSJoe Wallwork       if (rranks && rranks[r] == ngbRanks[p]) {
2002bc0b47dSJoe Wallwork 
2012bc0b47dSJoe Wallwork         /* Add the right slice of rmine at the right place */
2022bc0b47dSJoe Wallwork         sliceSize = roffset[r+1]-roffset[r];
2032bc0b47dSJoe Wallwork         for (j = 0, count = 0; j < sliceSize; ++j) {
2042bc0b47dSJoe Wallwork           if (roffset[r]+j >= roffset[nrranks]) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
2052bc0b47dSJoe Wallwork           v = rmine[roffset[r]+j];
2062bc0b47dSJoe Wallwork           if (v >= vStart && v < vEnd) {
2072bc0b47dSJoe Wallwork             if (intOffset[p+1]+count >= numVerNgbRanksTotal) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
2082bc0b47dSJoe Wallwork             interfaces_lv[intOffset[p+1]+count] = v-vStart;
2092bc0b47dSJoe Wallwork             count++;
2102bc0b47dSJoe Wallwork           }
2112bc0b47dSJoe Wallwork         }
2122bc0b47dSJoe Wallwork         intOffset[p+1] += count;
2132bc0b47dSJoe Wallwork         r++;
2142bc0b47dSJoe Wallwork       }
2152bc0b47dSJoe Wallwork       if (intOffset[p+1] != intOffset[p] + verNgbRank[p]) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unequal offsets");
2162bc0b47dSJoe Wallwork     }
2172bc0b47dSJoe Wallwork     ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
2182bc0b47dSJoe Wallwork     ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
2192bc0b47dSJoe Wallwork     for (i = 0; i < numVerNgbRanksTotal; ++i) {
2202bc0b47dSJoe Wallwork       v = gV[interfaces_lv[i]];
2212bc0b47dSJoe Wallwork       interfaces_gv[i] = v < 0 ? -v-1 : v;
2222bc0b47dSJoe Wallwork       interfaces_lv[i] += 1;
2232bc0b47dSJoe Wallwork       interfaces_gv[i] += 1;
2242bc0b47dSJoe Wallwork     }
2252bc0b47dSJoe Wallwork     ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
2262bc0b47dSJoe Wallwork     ierr = PetscFree(numVerInterfaces);CHKERRQ(ierr);
2272bc0b47dSJoe Wallwork   }
2282bc0b47dSJoe Wallwork   ierr = DMDestroy(&udm);CHKERRQ(ierr);
2292bc0b47dSJoe Wallwork 
2302bc0b47dSJoe Wallwork   /* Send the data to ParMmg and remesh */
2312bc0b47dSJoe Wallwork   ierr = DMPlexMetricNoInsertion(dm, &noInsert);CHKERRQ(ierr);
2322bc0b47dSJoe Wallwork   ierr = DMPlexMetricNoSwapping(dm, &noSwap);CHKERRQ(ierr);
2332bc0b47dSJoe Wallwork   ierr = DMPlexMetricNoMovement(dm, &noMove);CHKERRQ(ierr);
2342bc0b47dSJoe Wallwork   ierr = DMPlexMetricGetVerbosity(dm, &verbosity);CHKERRQ(ierr);
2352bc0b47dSJoe Wallwork   ierr = DMPlexMetricGetNumIterations(dm, &numIter);CHKERRQ(ierr);
2362bc0b47dSJoe Wallwork   ierr = DMPlexMetricGetGradationFactor(dm, &gradationFactor);CHKERRQ(ierr);
2372bc0b47dSJoe Wallwork   ierr = PMMG_Init_parMesh(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_pMesh, PMMG_ARG_pMet, PMMG_ARG_dim, 3, PMMG_ARG_MPIComm, tmpComm, PMMG_ARG_end);
238c1dc6da0SJoe Wallwork   ierr = PMMG_Set_meshSize(parmesh, numVertices, numCells, 0, numFaceTags, 0, 0);
2392bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_APImode, PMMG_APIDISTRIB_nodes);
2402bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noinsert, noInsert);
2412bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noswap, noSwap);
2422bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nomove, noMove);
2432bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_verbose, verbosity);
2442bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1);
2452bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_niter, numIter);
2462bc0b47dSJoe Wallwork   ierr = PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hgrad, gradationFactor);
2472bc0b47dSJoe Wallwork   ierr = PMMG_Set_vertices(parmesh, vertices, verTags);
248cf6859d0SJoe Wallwork   ierr = PMMG_Set_tetrahedra(parmesh, cells, cellTags);
249c1dc6da0SJoe Wallwork   ierr = PMMG_Set_triangles(parmesh, bdFaces, faceTags);
2502bc0b47dSJoe Wallwork   ierr = PMMG_Set_metSize(parmesh, MMG5_Vertex, numVertices, MMG5_Tensor);
251cf6859d0SJoe Wallwork   ierr = PMMG_Set_tensorMets(parmesh, metric);
2522bc0b47dSJoe Wallwork   ierr = PMMG_Set_numberOfNodeCommunicators(parmesh, numNgbRanks);
2532bc0b47dSJoe Wallwork   for (c = 0; c < numNgbRanks; ++c) {
2542bc0b47dSJoe Wallwork     ierr = PMMG_Set_ithNodeCommunicatorSize(parmesh, c, ngbRanks[c], intOffset[c+1]-intOffset[c]);
2552bc0b47dSJoe Wallwork     ierr = PMMG_Set_ithNodeCommunicator_nodes(parmesh, c, &interfaces_lv[intOffset[c]], &interfaces_gv[intOffset[c]], 1);
2562bc0b47dSJoe Wallwork   }
2572bc0b47dSJoe Wallwork   ierr = PMMG_parmmglib_distributed(parmesh);
2582bc0b47dSJoe Wallwork   ierr = PetscFree(cells);CHKERRQ(ierr);
2592bc0b47dSJoe Wallwork   ierr = PetscFree2(metric, vertices);CHKERRQ(ierr);
260c1dc6da0SJoe Wallwork   ierr = PetscFree2(bdFaces, faceTags);CHKERRQ(ierr);
2612bc0b47dSJoe Wallwork   ierr = PetscFree2(verTags, cellTags);CHKERRQ(ierr);
2622bc0b47dSJoe Wallwork   if (numProcs > 1) {
2632bc0b47dSJoe Wallwork     ierr = PetscFree2(ngbRanks, verNgbRank);CHKERRQ(ierr);
2642bc0b47dSJoe Wallwork     ierr = PetscFree3(interfaces_lv, interfaces_gv, intOffset);CHKERRQ(ierr);
2652bc0b47dSJoe Wallwork   }
2662bc0b47dSJoe Wallwork 
2672f583692SJoe Wallwork   /* Retrieve mesh from Mmg */
2682bc0b47dSJoe Wallwork   numCornersNew = 4;
2692bc0b47dSJoe Wallwork   ierr = PMMG_Get_meshSize(parmesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0);
2702bc0b47dSJoe Wallwork   ierr = PetscMalloc4(dim*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer);CHKERRQ(ierr);
2712bc0b47dSJoe Wallwork   ierr = PetscMalloc3((dim+1)*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells);CHKERRQ(ierr);
2722bc0b47dSJoe Wallwork   ierr = PetscMalloc4(dim*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces);CHKERRQ(ierr);
2732bc0b47dSJoe Wallwork   ierr = PMMG_Get_vertices(parmesh, verticesNew, verTagsNew, corners, requiredVer);
2742bc0b47dSJoe Wallwork   ierr = PMMG_Get_tetrahedra(parmesh, cellsNew, cellTagsNew, requiredCells);
2752bc0b47dSJoe Wallwork   ierr = PMMG_Get_triangles(parmesh, facesNew, faceTagsNew, requiredFaces);
2762bc0b47dSJoe Wallwork   ierr = PetscMalloc2(numVerticesNew, &owners, numVerticesNew, &gv_new);
2772bc0b47dSJoe Wallwork   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1);
2782bc0b47dSJoe Wallwork   ierr = PMMG_Get_verticesGloNum(parmesh, gv_new, owners);
2792bc0b47dSJoe Wallwork   for (i = 0; i < dim*numFacesNew; ++i) facesNew[i] -= 1;
2802bc0b47dSJoe Wallwork   for (i = 0; i < (dim+1)*numCellsNew; ++i) {
2812bc0b47dSJoe Wallwork     cellsNew[i] = gv_new[cellsNew[i]-1]-1;
2822bc0b47dSJoe Wallwork   }
2832bc0b47dSJoe Wallwork   numVerticesNewLoc = 0;
2842bc0b47dSJoe Wallwork   for (i = 0; i < numVerticesNew; ++i) {
2852bc0b47dSJoe Wallwork     if (owners[i] == rank) numVerticesNewLoc++;
2862bc0b47dSJoe Wallwork   }
2872bc0b47dSJoe Wallwork   ierr = PetscMalloc2(numVerticesNewLoc*dim, &verticesNewLoc, numVerticesNew, &verticesNewSorted);CHKERRQ(ierr);
2882bc0b47dSJoe Wallwork   for (i = 0, c = 0; i < numVerticesNew; i++) {
2892bc0b47dSJoe Wallwork     if (owners[i] == rank) {
2902bc0b47dSJoe Wallwork       for (j=0; j<dim; ++j) verticesNewLoc[dim*c+j] = verticesNew[dim*i+j];
2912bc0b47dSJoe Wallwork         c++;
2922bc0b47dSJoe Wallwork     }
2932bc0b47dSJoe Wallwork   }
2942f583692SJoe Wallwork 
2952f583692SJoe Wallwork   /* Reorder for consistency with DMPlex */
2962f583692SJoe Wallwork   for (i = 0; i < numCellsNew; ++i) { ierr = DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4*i]);CHKERRQ(ierr); }
2972f583692SJoe Wallwork 
2982f583692SJoe Wallwork   /* Create new plex */
2992bc0b47dSJoe Wallwork   ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNewLoc, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNewLoc, NULL, &verticesNewSorted, dmNew);CHKERRQ(ierr);
3002bc0b47dSJoe Wallwork   ierr = PMMG_Free_all(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_end);
301c1dc6da0SJoe Wallwork   ierr = PetscFree4(verticesNew, verTagsNew, corners, requiredVer);CHKERRQ(ierr);
302c1dc6da0SJoe Wallwork 
303c1dc6da0SJoe Wallwork   /* Get adapted mesh information */
304c1dc6da0SJoe Wallwork   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
305c1dc6da0SJoe Wallwork   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
306c1dc6da0SJoe Wallwork   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
3072bc0b47dSJoe Wallwork 
3082bc0b47dSJoe Wallwork   /* Rebuild boundary label */
3097881f211SJoe Wallwork   ierr = DMCreateLabel(*dmNew, flg ? bdName : bdLabelName);CHKERRQ(ierr);
3107881f211SJoe Wallwork   ierr = DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew);CHKERRQ(ierr);
3112bc0b47dSJoe Wallwork   for (i = 0; i < numFacesNew; i++) {
312c1dc6da0SJoe Wallwork     PetscBool       hasTag = PETSC_FALSE;
3132bc0b47dSJoe Wallwork     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
3142bc0b47dSJoe Wallwork     const PetscInt *coveredPoints = NULL;
3152bc0b47dSJoe Wallwork 
3162bc0b47dSJoe Wallwork     for (j = 0; j < dim; ++j) {
3172bc0b47dSJoe Wallwork       lv = facesNew[i*dim+j];
3182bc0b47dSJoe Wallwork       gv = gv_new[lv]-1;
3192bc0b47dSJoe Wallwork       ierr = PetscFindInt(gv, numVerticesNew, verticesNewSorted, &lv);CHKERRQ(ierr);
3202bc0b47dSJoe Wallwork       facePoints[j] = lv+vStart;
3212bc0b47dSJoe Wallwork     }
3222bc0b47dSJoe Wallwork     ierr = DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
3232bc0b47dSJoe Wallwork     for (j = 0; j < numCoveredPoints; ++j) {
3242bc0b47dSJoe Wallwork       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
3252bc0b47dSJoe Wallwork         numFaces++;
3262bc0b47dSJoe Wallwork         f = j;
3272bc0b47dSJoe Wallwork       }
3282bc0b47dSJoe Wallwork     }
3292bc0b47dSJoe Wallwork     if (numFaces != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "%d vertices cannot define more than 1 facet (%d)", dim, numFaces);
330c1dc6da0SJoe Wallwork     ierr = DMLabelHasStratum(bdLabel, faceTagsNew[i], &hasTag);CHKERRQ(ierr);
331c1dc6da0SJoe Wallwork     if (hasTag) { ierr = DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]);CHKERRQ(ierr); }
3322bc0b47dSJoe Wallwork     ierr = DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
3332bc0b47dSJoe Wallwork   }
334c1dc6da0SJoe Wallwork   ierr = PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces);CHKERRQ(ierr);
335c1dc6da0SJoe Wallwork   ierr = PetscFree2(owners, gv_new);CHKERRQ(ierr);
336c1dc6da0SJoe Wallwork   ierr = PetscFree2(verticesNewLoc, verticesNewSorted);CHKERRQ(ierr);
337c1dc6da0SJoe Wallwork   if (flg) { ierr = DMLabelDestroy(&bdLabel);CHKERRQ(ierr); }
3382bc0b47dSJoe Wallwork 
3399fe9e680SJoe Wallwork   /* Rebuild cell labels */
3409fe9e680SJoe Wallwork   ierr = DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName);CHKERRQ(ierr);
3419fe9e680SJoe Wallwork   ierr = DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew);CHKERRQ(ierr);
3429fe9e680SJoe Wallwork   for (c = cStart; c < cEnd; ++c) { ierr = DMLabelSetValue(rgLabelNew, c, cellTagsNew[c-cStart]);CHKERRQ(ierr); }
3432bc0b47dSJoe Wallwork   ierr = PetscFree3(cellsNew, cellTagsNew, requiredCells);CHKERRQ(ierr);
344c1dc6da0SJoe Wallwork 
3452bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
3462bc0b47dSJoe Wallwork }
347