15f80ce2aSJacob Faibussowitsch #include "../mmgcommon.h" /*I "petscdmplex.h" I*/
22bc0b47dSJoe Wallwork #include <mmg/libmmg.h>
32bc0b47dSJoe Wallwork
4b2c4de1dSPierre Jolivet PetscBool MmgCite = PETSC_FALSE;
5b2c4de1dSPierre Jolivet const char MmgCitation[] = "@article{DAPOGNY2014358,\n"
6b2c4de1dSPierre Jolivet " title = {Three-dimensional adaptive domain remeshing, implicit domain meshing, and applications to free and moving boundary problems},\n"
7b2c4de1dSPierre Jolivet " journal = {Journal of Computational Physics},\n"
8b2c4de1dSPierre Jolivet " author = {C. Dapogny and C. Dobrzynski and P. Frey},\n"
9b2c4de1dSPierre Jolivet " volume = {262},\n"
10b2c4de1dSPierre Jolivet " pages = {358--378},\n"
11b2c4de1dSPierre Jolivet " doi = {10.1016/j.jcp.2014.01.005},\n"
12b2c4de1dSPierre Jolivet " year = {2014}\n}\n";
13b2c4de1dSPierre Jolivet
DMAdaptMetric_Mmg_Plex(DM dm,Vec vertexMetric,DMLabel bdLabel,DMLabel rgLabel,DM * dmNew)14d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew)
15d71ae5a4SJacob Faibussowitsch {
162bc0b47dSJoe Wallwork MPI_Comm comm;
172bc0b47dSJoe Wallwork const char *bdName = "_boundary_";
189fe9e680SJoe Wallwork const char *rgName = "_regions_";
192bc0b47dSJoe Wallwork DM udm, cdm;
208d788f11SJoe Wallwork DMLabel bdLabelNew, rgLabelNew;
219fe9e680SJoe Wallwork const char *bdLabelName, *rgLabelName;
222bc0b47dSJoe Wallwork PetscSection coordSection;
232bc0b47dSJoe Wallwork Vec coordinates;
242bc0b47dSJoe Wallwork const PetscScalar *coords, *met;
25ae8b063eSJoe Wallwork PetscReal *vertices, *metric, *verticesNew, gradationFactor, hausdorffNumber;
268d788f11SJoe Wallwork PetscInt *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew;
278d788f11SJoe Wallwork PetscInt *bdFaces, *faceTags, *facesNew, *faceTagsNew;
28711aed9cSPierre Jolivet int *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces;
298d788f11SJoe Wallwork PetscInt cStart, cEnd, c, numCells, fStart, fEnd, numFaceTags, f, vStart, vEnd, v, numVertices;
308d788f11SJoe Wallwork PetscInt dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, pStart, pEnd;
312bc0b47dSJoe Wallwork PetscInt numCellsNew, numVerticesNew, numCornersNew, numFacesNew;
3276f360caSJoe Wallwork PetscBool flg = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform;
332bc0b47dSJoe Wallwork MMG5_pMesh mmg_mesh = NULL;
342bc0b47dSJoe Wallwork MMG5_pSol mmg_metric = NULL;
352bc0b47dSJoe Wallwork
362bc0b47dSJoe Wallwork PetscFunctionBegin;
37b2c4de1dSPierre Jolivet PetscCall(PetscCitationsRegister(MmgCitation, &MmgCite));
389566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
392bc0b47dSJoe Wallwork if (bdLabel) {
409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)bdLabel, &bdLabelName));
419566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(bdLabelName, bdName, &flg));
4228b400f6SJacob Faibussowitsch PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
432bc0b47dSJoe Wallwork }
449fe9e680SJoe Wallwork if (rgLabel) {
459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)rgLabel, &rgLabelName));
469566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(rgLabelName, rgName, &flg));
4728b400f6SJacob Faibussowitsch PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
489fe9e680SJoe Wallwork }
492bc0b47dSJoe Wallwork
502bc0b47dSJoe Wallwork /* Get mesh information */
519566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
522bc0b47dSJoe Wallwork Neq = (dim * (dim + 1)) / 2;
539566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
549566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
559566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
569566063dSJacob Faibussowitsch PetscCall(DMPlexUninterpolate(dm, &udm));
579566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
582bc0b47dSJoe Wallwork numCells = cEnd - cStart;
592bc0b47dSJoe Wallwork numVertices = vEnd - vStart;
602bc0b47dSJoe Wallwork
612bc0b47dSJoe Wallwork /* Get cell offsets */
629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells * maxConeSize, &cells));
632bc0b47dSJoe Wallwork for (c = 0, coff = 0; c < numCells; ++c) {
642bc0b47dSJoe Wallwork const PetscInt *cone;
652bc0b47dSJoe Wallwork PetscInt coneSize, cl;
662bc0b47dSJoe Wallwork
679566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(udm, c, &coneSize));
689566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(udm, c, &cone));
692bc0b47dSJoe Wallwork for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart + 1;
702bc0b47dSJoe Wallwork }
712bc0b47dSJoe Wallwork
722bc0b47dSJoe Wallwork /* Get vertex coordinate array */
739566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm));
749566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cdm, &coordSection));
759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords));
779566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numVertices * Neq, &metric, dim * numVertices, &vertices));
782bc0b47dSJoe Wallwork for (v = 0; v < vEnd - vStart; ++v) {
799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v + vStart, &off));
802bc0b47dSJoe Wallwork for (i = 0; i < dim; ++i) vertices[dim * v + i] = PetscRealPart(coords[off + i]);
812bc0b47dSJoe Wallwork }
829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords));
839566063dSJacob Faibussowitsch PetscCall(DMDestroy(&udm));
842bc0b47dSJoe Wallwork
858d788f11SJoe Wallwork /* Get face tags */
868d788f11SJoe Wallwork if (!bdLabel) {
878d788f11SJoe Wallwork flg = PETSC_TRUE;
889566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel));
899566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
908d788f11SJoe Wallwork }
919566063dSJacob Faibussowitsch PetscCall(DMLabelGetBounds(bdLabel, &pStart, &pEnd));
928d788f11SJoe Wallwork for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
938d788f11SJoe Wallwork PetscBool hasPoint;
948d788f11SJoe Wallwork PetscInt *closure = NULL, closureSize, cl;
952bc0b47dSJoe Wallwork
969566063dSJacob Faibussowitsch PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
978d788f11SJoe Wallwork if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
988d788f11SJoe Wallwork numFaceTags++;
998d788f11SJoe Wallwork
1009566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1012bc0b47dSJoe Wallwork for (cl = 0; cl < closureSize * 2; cl += 2) {
1022bc0b47dSJoe Wallwork if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1032bc0b47dSJoe Wallwork }
1049566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1052bc0b47dSJoe Wallwork }
1069566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags));
1078d788f11SJoe Wallwork for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) {
1088d788f11SJoe Wallwork PetscBool hasPoint;
1098d788f11SJoe Wallwork PetscInt *closure = NULL, closureSize, cl;
1102bc0b47dSJoe Wallwork
1119566063dSJacob Faibussowitsch PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint));
1128d788f11SJoe Wallwork if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue;
1138d788f11SJoe Wallwork
1149566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
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 }
1189566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
1199566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]));
1202bc0b47dSJoe Wallwork }
1219566063dSJacob Faibussowitsch if (flg) PetscCall(DMLabelDestroy(&bdLabel));
1222bc0b47dSJoe Wallwork
1239fe9e680SJoe Wallwork /* Get cell tags */
1249566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numVertices, &verTags, numCells, &cellTags));
1259fe9e680SJoe Wallwork if (rgLabel) {
1269566063dSJacob Faibussowitsch for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelGetValue(rgLabel, c, &cellTags[c]));
1279fe9e680SJoe Wallwork }
1289fe9e680SJoe Wallwork
1292bc0b47dSJoe Wallwork /* Get metric */
1309566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
1319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vertexMetric, &met));
1329566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1339566063dSJacob Faibussowitsch PetscCall(DMPlexMetricIsUniform(dm, &uniform));
1342bc0b47dSJoe Wallwork for (v = 0; v < (vEnd - vStart); ++v) {
1352bc0b47dSJoe Wallwork for (i = 0, k = 0; i < dim; ++i) {
1362bc0b47dSJoe Wallwork for (j = i; j < dim; ++j) {
13793520041SJoe Wallwork if (isotropic) {
13893520041SJoe Wallwork if (i == j) {
13993520041SJoe Wallwork if (uniform) metric[Neq * v + k] = PetscRealPart(met[0]);
14093520041SJoe Wallwork else metric[Neq * v + k] = PetscRealPart(met[v]);
14193520041SJoe Wallwork } else metric[Neq * v + k] = 0.0;
14293520041SJoe Wallwork } else {
1432bc0b47dSJoe Wallwork metric[Neq * v + k] = PetscRealPart(met[dim * dim * v + dim * i + j]);
14493520041SJoe Wallwork }
1452bc0b47dSJoe Wallwork k++;
1462bc0b47dSJoe Wallwork }
1472bc0b47dSJoe Wallwork }
1482bc0b47dSJoe Wallwork }
1499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vertexMetric, &met));
1502bc0b47dSJoe Wallwork
1512bc0b47dSJoe Wallwork /* Send mesh to Mmg and remesh */
1529566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetVerbosity(dm, &verbosity));
1539566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetGradationFactor(dm, &gradationFactor));
1549566063dSJacob Faibussowitsch PetscCall(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber));
1559566063dSJacob Faibussowitsch PetscCall(DMPlexMetricNoInsertion(dm, &noInsert));
1569566063dSJacob Faibussowitsch PetscCall(DMPlexMetricNoSwapping(dm, &noSwap));
1579566063dSJacob Faibussowitsch PetscCall(DMPlexMetricNoMovement(dm, &noMove));
1589566063dSJacob Faibussowitsch PetscCall(DMPlexMetricNoSurf(dm, &noSurf));
1592bc0b47dSJoe Wallwork switch (dim) {
1602bc0b47dSJoe Wallwork case 2:
161*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Init_mesh, MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
162*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter, mmg_mesh, mmg_metric, MMG2D_IPARAM_noinsert, noInsert);
163*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter, mmg_mesh, mmg_metric, MMG2D_IPARAM_noswap, noSwap);
164*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter, mmg_mesh, mmg_metric, MMG2D_IPARAM_nomove, noMove);
165*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter, mmg_mesh, mmg_metric, MMG2D_IPARAM_nosurf, noSurf);
166*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter, mmg_mesh, mmg_metric, MMG2D_IPARAM_verbose, verbosity);
167*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_dparameter, mmg_mesh, mmg_metric, MMG2D_DPARAM_hgrad, gradationFactor);
168*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_dparameter, mmg_mesh, mmg_metric, MMG2D_DPARAM_hausd, hausdorffNumber);
169*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_meshSize, mmg_mesh, numVertices, numCells, 0, numFaceTags);
170*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_vertices, mmg_mesh, vertices, verTags);
171*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_triangles, mmg_mesh, cells, cellTags);
172*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_edges, mmg_mesh, bdFaces, faceTags);
173*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_solSize, mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor);
174*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Set_tensorSols, mmg_metric, metric);
175*b3377c29SStefano Zampini PetscCallMMG(MMG2D_mmg2dlib, mmg_mesh, mmg_metric);
1762bc0b47dSJoe Wallwork break;
1772bc0b47dSJoe Wallwork case 3:
178*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Init_mesh, MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
179*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter, mmg_mesh, mmg_metric, MMG3D_IPARAM_noinsert, noInsert);
180*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter, mmg_mesh, mmg_metric, MMG3D_IPARAM_noswap, noSwap);
181*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter, mmg_mesh, mmg_metric, MMG3D_IPARAM_nomove, noMove);
182*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter, mmg_mesh, mmg_metric, MMG3D_IPARAM_nosurf, noSurf);
183*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter, mmg_mesh, mmg_metric, MMG3D_IPARAM_verbose, verbosity);
184*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_dparameter, mmg_mesh, mmg_metric, MMG3D_DPARAM_hgrad, gradationFactor);
185*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_dparameter, mmg_mesh, mmg_metric, MMG3D_DPARAM_hausd, hausdorffNumber);
186*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_meshSize, mmg_mesh, numVertices, numCells, 0, numFaceTags, 0, 0);
187*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_vertices, mmg_mesh, vertices, verTags);
188*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_tetrahedra, mmg_mesh, cells, cellTags);
189*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_triangles, mmg_mesh, bdFaces, faceTags);
190*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_solSize, mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor);
191*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Set_tensorSols, mmg_metric, metric);
192*b3377c29SStefano Zampini PetscCallMMG(MMG3D_mmg3dlib, mmg_mesh, mmg_metric);
1932bc0b47dSJoe Wallwork break;
194d71ae5a4SJacob Faibussowitsch default:
195d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim);
1962bc0b47dSJoe Wallwork }
1979566063dSJacob Faibussowitsch PetscCall(PetscFree(cells));
1989566063dSJacob Faibussowitsch PetscCall(PetscFree2(metric, vertices));
1999566063dSJacob Faibussowitsch PetscCall(PetscFree2(bdFaces, faceTags));
2009566063dSJacob Faibussowitsch PetscCall(PetscFree2(verTags, cellTags));
2012bc0b47dSJoe Wallwork
20299b4fe00SJoe Wallwork /* Retrieve mesh from Mmg */
2032bc0b47dSJoe Wallwork switch (dim) {
2042bc0b47dSJoe Wallwork case 2:
2052bc0b47dSJoe Wallwork numCornersNew = 3;
206*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Get_meshSize, mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew);
2079566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(2 * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
2089566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(3 * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
2099566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(2 * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
210*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Get_vertices, mmg_mesh, verticesNew, verTagsNew, corners, requiredVer);
211*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Get_triangles, mmg_mesh, cellsNew, cellTagsNew, requiredCells);
212*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Get_edges, mmg_mesh, facesNew, faceTagsNew, ridges, requiredFaces);
2132bc0b47dSJoe Wallwork break;
2142bc0b47dSJoe Wallwork case 3:
2152bc0b47dSJoe Wallwork numCornersNew = 4;
216*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Get_meshSize, mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0);
2179566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(3 * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer));
2189566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(4 * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells));
2199566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(3 * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces));
220*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Get_vertices, mmg_mesh, verticesNew, verTagsNew, corners, requiredVer);
221*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Get_tetrahedra, mmg_mesh, cellsNew, cellTagsNew, requiredCells);
222*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Get_triangles, mmg_mesh, facesNew, faceTagsNew, requiredFaces);
22399b4fe00SJoe Wallwork
22499b4fe00SJoe Wallwork /* Reorder for consistency with DMPlex */
2259566063dSJacob Faibussowitsch for (i = 0; i < numCellsNew; ++i) PetscCall(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4 * i]));
2262bc0b47dSJoe Wallwork break;
22799b4fe00SJoe Wallwork
228d71ae5a4SJacob Faibussowitsch default:
229d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim);
2302bc0b47dSJoe Wallwork }
23199b4fe00SJoe Wallwork
23299b4fe00SJoe Wallwork /* Create new Plex */
2332bc0b47dSJoe Wallwork for (i = 0; i < (dim + 1) * numCellsNew; i++) cellsNew[i] -= 1;
2342bc0b47dSJoe Wallwork for (i = 0; i < dim * numFacesNew; i++) facesNew[i] -= 1;
2359566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNew, NULL, NULL, dmNew));
2362bc0b47dSJoe Wallwork switch (dim) {
237d71ae5a4SJacob Faibussowitsch case 2:
238*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG2D_Free_all, MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
239d71ae5a4SJacob Faibussowitsch break;
240d71ae5a4SJacob Faibussowitsch case 3:
241*b3377c29SStefano Zampini PetscCallMMG_NONSTANDARD(MMG3D_Free_all, MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
242d71ae5a4SJacob Faibussowitsch break;
243d71ae5a4SJacob Faibussowitsch default:
244d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim);
2452bc0b47dSJoe Wallwork }
2469566063dSJacob Faibussowitsch PetscCall(PetscFree4(verticesNew, verTagsNew, corners, requiredVer));
2478d788f11SJoe Wallwork
2488d788f11SJoe Wallwork /* Get adapted mesh information */
2499566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
2509566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
2519566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
2522bc0b47dSJoe Wallwork
2532bc0b47dSJoe Wallwork /* Rebuild boundary labels */
2549566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName));
2559566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew));
2562bc0b47dSJoe Wallwork for (i = 0; i < numFacesNew; i++) {
2572bc0b47dSJoe Wallwork PetscInt numCoveredPoints, numFaces = 0, facePoints[3];
2582bc0b47dSJoe Wallwork const PetscInt *coveredPoints = NULL;
2592bc0b47dSJoe Wallwork
2602bc0b47dSJoe Wallwork for (j = 0; j < dim; ++j) facePoints[j] = facesNew[i * dim + j] + vStart;
2619566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
2622bc0b47dSJoe Wallwork for (j = 0; j < numCoveredPoints; ++j) {
2632bc0b47dSJoe Wallwork if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
2642bc0b47dSJoe Wallwork numFaces++;
2652bc0b47dSJoe Wallwork f = j;
2662bc0b47dSJoe Wallwork }
2672bc0b47dSJoe Wallwork }
2685f80ce2aSJacob Faibussowitsch PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces);
2699566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]));
2709566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints));
2712bc0b47dSJoe Wallwork }
2729566063dSJacob Faibussowitsch PetscCall(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces));
2732bc0b47dSJoe Wallwork
2749fe9e680SJoe Wallwork /* Rebuild cell labels */
2759566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName));
2769566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew));
2779566063dSJacob Faibussowitsch for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c - cStart]));
2789566063dSJacob Faibussowitsch PetscCall(PetscFree3(cellsNew, cellTagsNew, requiredCells));
2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2802bc0b47dSJoe Wallwork }
281