xref: /petsc/src/dm/impls/plex/adaptors/mmg/mmgadapt.c (revision b2c4de1d22ff155d2c936e49c947998b8cbd29f0)
15f80ce2aSJacob Faibussowitsch #include "../mmgcommon.h" /*I      "petscdmplex.h"   I*/
22bc0b47dSJoe Wallwork #include <mmg/libmmg.h>
32bc0b47dSJoe Wallwork 
4*b2c4de1dSPierre Jolivet PetscBool  MmgCite       = PETSC_FALSE;
5*b2c4de1dSPierre Jolivet const char MmgCitation[] = "@article{DAPOGNY2014358,\n"
6*b2c4de1dSPierre Jolivet                            "  title   = {Three-dimensional adaptive domain remeshing, implicit domain meshing, and applications to free and moving boundary problems},\n"
7*b2c4de1dSPierre Jolivet                            "  journal = {Journal of Computational Physics},\n"
8*b2c4de1dSPierre Jolivet                            "  author  = {C. Dapogny and C. Dobrzynski and P. Frey},\n"
9*b2c4de1dSPierre Jolivet                            "  volume  = {262},\n"
10*b2c4de1dSPierre Jolivet                            "  pages   = {358--378},\n"
11*b2c4de1dSPierre Jolivet                            "  doi     = {10.1016/j.jcp.2014.01.005},\n"
12*b2c4de1dSPierre Jolivet                            "  year    = {2014}\n}\n";
13*b2c4de1dSPierre Jolivet 
149371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew) {
152bc0b47dSJoe Wallwork   MPI_Comm           comm;
162bc0b47dSJoe Wallwork   const char        *bdName = "_boundary_";
179fe9e680SJoe Wallwork   const char        *rgName = "_regions_";
182bc0b47dSJoe Wallwork   DM                 udm, cdm;
198d788f11SJoe Wallwork   DMLabel            bdLabelNew, rgLabelNew;
209fe9e680SJoe Wallwork   const char        *bdLabelName, *rgLabelName;
212bc0b47dSJoe Wallwork   PetscSection       coordSection;
222bc0b47dSJoe Wallwork   Vec                coordinates;
232bc0b47dSJoe Wallwork   const PetscScalar *coords, *met;
24ae8b063eSJoe Wallwork   PetscReal         *vertices, *metric, *verticesNew, gradationFactor, hausdorffNumber;
258d788f11SJoe Wallwork   PetscInt          *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew;
268d788f11SJoe Wallwork   PetscInt          *bdFaces, *faceTags, *facesNew, *faceTagsNew;
272bc0b47dSJoe Wallwork   PetscInt          *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces;
288d788f11SJoe Wallwork   PetscInt           cStart, cEnd, c, numCells, fStart, fEnd, numFaceTags, f, vStart, vEnd, v, numVertices;
298d788f11SJoe Wallwork   PetscInt           dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, pStart, pEnd;
302bc0b47dSJoe Wallwork   PetscInt           numCellsNew, numVerticesNew, numCornersNew, numFacesNew;
3176f360caSJoe Wallwork   PetscBool          flg        = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform;
322bc0b47dSJoe Wallwork   MMG5_pMesh         mmg_mesh   = NULL;
332bc0b47dSJoe Wallwork   MMG5_pSol          mmg_metric = NULL;
342bc0b47dSJoe Wallwork 
352bc0b47dSJoe Wallwork   PetscFunctionBegin;
36*b2c4de1dSPierre Jolivet   PetscCall(PetscCitationsRegister(MmgCitation, &MmgCite));
379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
382bc0b47dSJoe Wallwork   if (bdLabel) {
399566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)bdLabel, &bdLabelName));
409566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(bdLabelName, bdName, &flg));
4128b400f6SJacob Faibussowitsch     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
422bc0b47dSJoe Wallwork   }
439fe9e680SJoe Wallwork   if (rgLabel) {
449566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)rgLabel, &rgLabelName));
459566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(rgLabelName, rgName, &flg));
4628b400f6SJacob Faibussowitsch     PetscCheck(!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 */
509566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
512bc0b47dSJoe Wallwork   Neq = (dim * (dim + 1)) / 2;
529566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
539566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
549566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
559566063dSJacob Faibussowitsch   PetscCall(DMPlexUninterpolate(dm, &udm));
569566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
572bc0b47dSJoe Wallwork   numCells    = cEnd - cStart;
582bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
592bc0b47dSJoe Wallwork 
602bc0b47dSJoe Wallwork   /* Get cell offsets */
619566063dSJacob Faibussowitsch   PetscCall(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 
669566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(udm, c, &coneSize));
679566063dSJacob Faibussowitsch     PetscCall(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 */
729566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
739566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(cdm, &coordSection));
749566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numVertices * Neq, &metric, dim * numVertices, &vertices));
772bc0b47dSJoe Wallwork   for (v = 0; v < vEnd - vStart; ++v) {
789566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(coordSection, v + vStart, &off));
792bc0b47dSJoe Wallwork     for (i = 0; i < dim; ++i) vertices[dim * v + i] = PetscRealPart(coords[off + i]);
802bc0b47dSJoe Wallwork   }
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
829566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&udm));
832bc0b47dSJoe Wallwork 
848d788f11SJoe Wallwork   /* Get face tags */
858d788f11SJoe Wallwork   if (!bdLabel) {
868d788f11SJoe Wallwork     flg = PETSC_TRUE;
879566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel));
889566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
898d788f11SJoe Wallwork   }
909566063dSJacob Faibussowitsch   PetscCall(DMLabelGetBounds(bdLabel, &pStart, &pEnd));
918d788f11SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
928d788f11SJoe Wallwork     PetscBool hasPoint;
938d788f11SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
942bc0b47dSJoe Wallwork 
959566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
968d788f11SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
978d788f11SJoe Wallwork     numFaceTags++;
988d788f11SJoe Wallwork 
999566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1002bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize * 2; cl += 2) {
1012bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1022bc0b47dSJoe Wallwork     }
1039566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1042bc0b47dSJoe Wallwork   }
1059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags));
1068d788f11SJoe Wallwork   for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
1078d788f11SJoe Wallwork     PetscBool hasPoint;
1088d788f11SJoe Wallwork     PetscInt *closure = NULL, closureSize, cl;
1092bc0b47dSJoe Wallwork 
1109566063dSJacob Faibussowitsch     PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
1118d788f11SJoe Wallwork     if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
1128d788f11SJoe Wallwork 
1139566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
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     }
1179566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1189566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]));
1192bc0b47dSJoe Wallwork   }
1209566063dSJacob Faibussowitsch   if (flg) PetscCall(DMLabelDestroy(&bdLabel));
1212bc0b47dSJoe Wallwork 
1229fe9e680SJoe Wallwork   /* Get cell tags */
1239566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(numVertices, &verTags, numCells, &cellTags));
1249fe9e680SJoe Wallwork   if (rgLabel) {
1259566063dSJacob Faibussowitsch     for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelGetValue(rgLabel, c, &cellTags[c]));
1269fe9e680SJoe Wallwork   }
1279fe9e680SJoe Wallwork 
1282bc0b47dSJoe Wallwork   /* Get metric */
1299566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(vertexMetric, &met));
1319566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1329566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1332bc0b47dSJoe Wallwork   for (v = 0; v < (vEnd - vStart); ++v) {
1342bc0b47dSJoe Wallwork     for (i = 0, k = 0; i < dim; ++i) {
1352bc0b47dSJoe Wallwork       for (j = i; j < dim; ++j) {
13693520041SJoe Wallwork         if (isotropic) {
13793520041SJoe Wallwork           if (i == j) {
13893520041SJoe Wallwork             if (uniform) metric[Neq * v + k] = PetscRealPart(met[0]);
13993520041SJoe Wallwork             else metric[Neq * v + k] = PetscRealPart(met[v]);
14093520041SJoe Wallwork           } else metric[Neq * v + k] = 0.0;
14193520041SJoe Wallwork         } else {
1422bc0b47dSJoe Wallwork           metric[Neq * v + k] = PetscRealPart(met[dim * dim * v + dim * i + j]);
14393520041SJoe Wallwork         }
1442bc0b47dSJoe Wallwork         k++;
1452bc0b47dSJoe Wallwork       }
1462bc0b47dSJoe Wallwork     }
1472bc0b47dSJoe Wallwork   }
1489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(vertexMetric, &met));
1492bc0b47dSJoe Wallwork 
1502bc0b47dSJoe Wallwork   /* Send mesh to Mmg and remesh */
1519566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetVerbosity(dm, &verbosity));
1529566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetGradationFactor(dm, &gradationFactor));
1539566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber));
1549566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoInsertion(dm, &noInsert));
1559566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoSwapping(dm, &noSwap));
1569566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoMovement(dm, &noMove));
1579566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNoSurf(dm, &noSurf));
1582bc0b47dSJoe Wallwork   switch (dim) {
1592bc0b47dSJoe Wallwork   case 2:
1609566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end));
1619566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noinsert, noInsert));
1629566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noswap, noSwap));
1639566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nomove, noMove));
1649566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nosurf, noSurf));
1659566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_verbose, verbosity));
1669566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hgrad, gradationFactor));
1679566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hausd, hausdorffNumber));
1689566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags));
1699566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_vertices(mmg_mesh, vertices, verTags));
1709566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_triangles(mmg_mesh, cells, cellTags));
1719566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_edges(mmg_mesh, bdFaces, faceTags));
1729566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor));
1739566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Set_tensorSols(mmg_metric, metric));
1749566063dSJacob Faibussowitsch     PetscCallMMG(MMG2D_mmg2dlib(mmg_mesh, mmg_metric));
1752bc0b47dSJoe Wallwork     break;
1762bc0b47dSJoe Wallwork   case 3:
1779566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end));
1789566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noinsert, noInsert));
1799566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noswap, noSwap));
1809566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nomove, noMove));
1819566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nosurf, noSurf));
1829566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_verbose, verbosity));
1839566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG3D_DPARAM_hgrad, gradationFactor));
184c387a134SPierre Jolivet     PetscCallMMG_NONSTANDARD(MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG3D_DPARAM_hausd, hausdorffNumber));
1859566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags, 0, 0));
1869566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_vertices(mmg_mesh, vertices, verTags));
1879566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_tetrahedra(mmg_mesh, cells, cellTags));
1889566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_triangles(mmg_mesh, bdFaces, faceTags));
1899566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor));
1909566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Set_tensorSols(mmg_metric, metric));
1919566063dSJacob Faibussowitsch     PetscCallMMG(MMG3D_mmg3dlib(mmg_mesh, mmg_metric));
1922bc0b47dSJoe Wallwork     break;
1935f80ce2aSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim);
1942bc0b47dSJoe Wallwork   }
1959566063dSJacob Faibussowitsch   PetscCall(PetscFree(cells));
1969566063dSJacob Faibussowitsch   PetscCall(PetscFree2(metric, vertices));
1979566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bdFaces, faceTags));
1989566063dSJacob Faibussowitsch   PetscCall(PetscFree2(verTags, cellTags));
1992bc0b47dSJoe Wallwork 
20099b4fe00SJoe Wallwork   /* Retrieve mesh from Mmg */
2012bc0b47dSJoe Wallwork   switch (dim) {
2022bc0b47dSJoe Wallwork   case 2:
2032bc0b47dSJoe Wallwork     numCornersNew = 3;
2049566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew));
2059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(2 * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
2069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(3 * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
2079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(2 * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
2089566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer));
2099566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Get_triangles(mmg_mesh, cellsNew, cellTagsNew, requiredCells));
2109566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG2D_Get_edges(mmg_mesh, facesNew, faceTagsNew, ridges, requiredFaces));
2112bc0b47dSJoe Wallwork     break;
2122bc0b47dSJoe Wallwork   case 3:
2132bc0b47dSJoe Wallwork     numCornersNew = 4;
2149566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0));
2159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(3 * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
2169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(4 * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
2179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(3 * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
2189566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer));
2199566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Get_tetrahedra(mmg_mesh, cellsNew, cellTagsNew, requiredCells));
2209566063dSJacob Faibussowitsch     PetscCallMMG_NONSTANDARD(MMG3D_Get_triangles(mmg_mesh, facesNew, faceTagsNew, requiredFaces));
22199b4fe00SJoe Wallwork 
22299b4fe00SJoe Wallwork     /* Reorder for consistency with DMPlex */
2239566063dSJacob Faibussowitsch     for (i = 0; i < numCellsNew; ++i) PetscCall(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4 * i]));
2242bc0b47dSJoe Wallwork     break;
22599b4fe00SJoe Wallwork 
2265f80ce2aSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim);
2272bc0b47dSJoe Wallwork   }
22899b4fe00SJoe Wallwork 
22999b4fe00SJoe Wallwork   /* Create new Plex */
2302bc0b47dSJoe Wallwork   for (i = 0; i < (dim + 1) * numCellsNew; i++) cellsNew[i] -= 1;
2312bc0b47dSJoe Wallwork   for (i = 0; i < dim * numFacesNew; i++) facesNew[i] -= 1;
2329566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNew, NULL, NULL, dmNew));
2332bc0b47dSJoe Wallwork   switch (dim) {
2349371c9d4SSatish Balay   case 2: PetscCallMMG_NONSTANDARD(MMG2D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); break;
2359371c9d4SSatish Balay   case 3: PetscCallMMG_NONSTANDARD(MMG3D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); break;
2365f80ce2aSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim);
2372bc0b47dSJoe Wallwork   }
2389566063dSJacob Faibussowitsch   PetscCall(PetscFree4(verticesNew, verTagsNew, corners, requiredVer));
2398d788f11SJoe Wallwork 
2408d788f11SJoe Wallwork   /* Get adapted mesh information */
2419566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
2429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
2439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
2442bc0b47dSJoe Wallwork 
2452bc0b47dSJoe Wallwork   /* Rebuild boundary labels */
2469566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName));
2479566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew));
2482bc0b47dSJoe Wallwork   for (i = 0; i < numFacesNew; i++) {
2492bc0b47dSJoe Wallwork     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
2502bc0b47dSJoe Wallwork     const PetscInt *coveredPoints = NULL;
2512bc0b47dSJoe Wallwork 
2522bc0b47dSJoe Wallwork     for (j = 0; j < dim; ++j) facePoints[j] = facesNew[i * dim + j] + vStart;
2539566063dSJacob Faibussowitsch     PetscCall(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
2542bc0b47dSJoe Wallwork     for (j = 0; j < numCoveredPoints; ++j) {
2552bc0b47dSJoe Wallwork       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
2562bc0b47dSJoe Wallwork         numFaces++;
2572bc0b47dSJoe Wallwork         f = j;
2582bc0b47dSJoe Wallwork       }
2592bc0b47dSJoe Wallwork     }
2605f80ce2aSJacob Faibussowitsch     PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces);
2619566063dSJacob Faibussowitsch     PetscCall(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]));
2629566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
2632bc0b47dSJoe Wallwork   }
2649566063dSJacob Faibussowitsch   PetscCall(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces));
2652bc0b47dSJoe Wallwork 
2669fe9e680SJoe Wallwork   /* Rebuild cell labels */
2679566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName));
2689566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew));
2699566063dSJacob Faibussowitsch   for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c - cStart]));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cellsNew, cellTagsNew, requiredCells));
2718d788f11SJoe Wallwork 
2722bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
2732bc0b47dSJoe Wallwork }
274