12bc0b47dSJoe Wallwork #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 22bc0b47dSJoe Wallwork #include <mmg/libmmg.h> 32bc0b47dSJoe Wallwork 49fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew) 52bc0b47dSJoe Wallwork { 62bc0b47dSJoe Wallwork MPI_Comm comm; 72bc0b47dSJoe Wallwork const char *bdName = "_boundary_"; 89fe9e680SJoe Wallwork const char *rgName = "_regions_"; 92bc0b47dSJoe Wallwork DM udm, cdm; 108d788f11SJoe Wallwork DMLabel bdLabelNew, rgLabelNew; 119fe9e680SJoe Wallwork const char *bdLabelName, *rgLabelName; 122bc0b47dSJoe Wallwork PetscSection coordSection; 132bc0b47dSJoe Wallwork Vec coordinates; 142bc0b47dSJoe Wallwork const PetscScalar *coords, *met; 15*ae8b063eSJoe Wallwork PetscReal *vertices, *metric, *verticesNew, gradationFactor, hausdorffNumber; 168d788f11SJoe Wallwork PetscInt *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew; 178d788f11SJoe Wallwork PetscInt *bdFaces, *faceTags, *facesNew, *faceTagsNew; 182bc0b47dSJoe Wallwork PetscInt *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces; 198d788f11SJoe Wallwork PetscInt cStart, cEnd, c, numCells, fStart, fEnd, numFaceTags, f, vStart, vEnd, v, numVertices; 208d788f11SJoe Wallwork PetscInt dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, pStart, pEnd; 212bc0b47dSJoe Wallwork PetscInt numCellsNew, numVerticesNew, numCornersNew, numFacesNew; 2293520041SJoe Wallwork PetscBool flg = PETSC_FALSE, noInsert, noSwap, noMove, isotropic, uniform; 232bc0b47dSJoe Wallwork MMG5_pMesh mmg_mesh = NULL; 242bc0b47dSJoe Wallwork MMG5_pSol mmg_metric = NULL; 252bc0b47dSJoe Wallwork PetscErrorCode ierr; 262bc0b47dSJoe Wallwork 272bc0b47dSJoe Wallwork PetscFunctionBegin; 282bc0b47dSJoe Wallwork ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 292bc0b47dSJoe Wallwork if (bdLabel) { 302bc0b47dSJoe Wallwork ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr); 312bc0b47dSJoe Wallwork ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr); 322c71b3e2SJacob Faibussowitsch PetscCheckFalse(flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 332bc0b47dSJoe Wallwork } 349fe9e680SJoe Wallwork if (rgLabel) { 359fe9e680SJoe Wallwork ierr = PetscObjectGetName((PetscObject) rgLabel, &rgLabelName);CHKERRQ(ierr); 369fe9e680SJoe Wallwork ierr = PetscStrcmp(rgLabelName, rgName, &flg);CHKERRQ(ierr); 372c71b3e2SJacob Faibussowitsch PetscCheckFalse(flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName); 389fe9e680SJoe Wallwork } 392bc0b47dSJoe Wallwork 402bc0b47dSJoe Wallwork /* Get mesh information */ 412bc0b47dSJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 422bc0b47dSJoe Wallwork Neq = (dim*(dim+1))/2; 432bc0b47dSJoe Wallwork ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 448d788f11SJoe Wallwork ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 452bc0b47dSJoe Wallwork ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 462bc0b47dSJoe Wallwork ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 472bc0b47dSJoe Wallwork ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 482bc0b47dSJoe Wallwork numCells = cEnd - cStart; 492bc0b47dSJoe Wallwork numVertices = vEnd - vStart; 502bc0b47dSJoe Wallwork 512bc0b47dSJoe Wallwork /* Get cell offsets */ 522bc0b47dSJoe Wallwork ierr = PetscMalloc1(numCells*maxConeSize, &cells);CHKERRQ(ierr); 532bc0b47dSJoe Wallwork for (c = 0, coff = 0; c < numCells; ++c) { 542bc0b47dSJoe Wallwork const PetscInt *cone; 552bc0b47dSJoe Wallwork PetscInt coneSize, cl; 562bc0b47dSJoe Wallwork 572bc0b47dSJoe Wallwork ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 582bc0b47dSJoe Wallwork ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 592bc0b47dSJoe Wallwork for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart + 1; 602bc0b47dSJoe Wallwork } 612bc0b47dSJoe Wallwork 622bc0b47dSJoe Wallwork /* Get vertex coordinate array */ 632bc0b47dSJoe Wallwork ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 642bc0b47dSJoe Wallwork ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr); 652bc0b47dSJoe Wallwork ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 662bc0b47dSJoe Wallwork ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 672bc0b47dSJoe Wallwork ierr = PetscMalloc2(numVertices*Neq, &metric, dim*numVertices, &vertices);CHKERRQ(ierr); 682bc0b47dSJoe Wallwork for (v = 0; v < vEnd-vStart; ++v) { 692bc0b47dSJoe Wallwork ierr = PetscSectionGetOffset(coordSection, v+vStart, &off);CHKERRQ(ierr); 702bc0b47dSJoe Wallwork for (i = 0; i < dim; ++i) vertices[dim*v+i] = PetscRealPart(coords[off+i]); 712bc0b47dSJoe Wallwork } 722bc0b47dSJoe Wallwork ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 732bc0b47dSJoe Wallwork ierr = DMDestroy(&udm);CHKERRQ(ierr); 742bc0b47dSJoe Wallwork 758d788f11SJoe Wallwork /* Get face tags */ 768d788f11SJoe Wallwork if (!bdLabel) { 778d788f11SJoe Wallwork flg = PETSC_TRUE; 788d788f11SJoe Wallwork ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel);CHKERRQ(ierr); 798d788f11SJoe Wallwork ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabel);CHKERRQ(ierr); 808d788f11SJoe Wallwork } 818d788f11SJoe Wallwork ierr = DMLabelGetBounds(bdLabel, &pStart, &pEnd);CHKERRQ(ierr); 828d788f11SJoe Wallwork for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) { 838d788f11SJoe Wallwork PetscBool hasPoint; 848d788f11SJoe Wallwork PetscInt *closure = NULL, closureSize, cl; 852bc0b47dSJoe Wallwork 868d788f11SJoe Wallwork ierr = DMLabelHasPoint(bdLabel, f, &hasPoint);CHKERRQ(ierr); 878d788f11SJoe Wallwork if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue; 888d788f11SJoe Wallwork numFaceTags++; 898d788f11SJoe Wallwork 908d788f11SJoe Wallwork ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 912bc0b47dSJoe Wallwork for (cl = 0; cl < closureSize*2; cl += 2) { 922bc0b47dSJoe Wallwork if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 932bc0b47dSJoe Wallwork } 948d788f11SJoe Wallwork ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 952bc0b47dSJoe Wallwork } 968d788f11SJoe Wallwork ierr = PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags);CHKERRQ(ierr); 978d788f11SJoe Wallwork for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) { 988d788f11SJoe Wallwork PetscBool hasPoint; 998d788f11SJoe Wallwork PetscInt *closure = NULL, closureSize, cl; 1002bc0b47dSJoe Wallwork 1018d788f11SJoe Wallwork ierr = DMLabelHasPoint(bdLabel, f, &hasPoint);CHKERRQ(ierr); 1028d788f11SJoe Wallwork if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue; 1038d788f11SJoe Wallwork 1048d788f11SJoe Wallwork ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1052bc0b47dSJoe Wallwork for (cl = 0; cl < closureSize*2; cl += 2) { 1062bc0b47dSJoe Wallwork if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart + 1; 1072bc0b47dSJoe Wallwork } 1088d788f11SJoe Wallwork ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1098d788f11SJoe Wallwork ierr = DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++]);CHKERRQ(ierr); 1102bc0b47dSJoe Wallwork } 1118d788f11SJoe Wallwork if (flg) { ierr = DMLabelDestroy(&bdLabel);CHKERRQ(ierr); } 1122bc0b47dSJoe Wallwork 1139fe9e680SJoe Wallwork /* Get cell tags */ 1149fe9e680SJoe Wallwork ierr = PetscCalloc2(numVertices, &verTags, numCells, &cellTags);CHKERRQ(ierr); 1159fe9e680SJoe Wallwork if (rgLabel) { 1169fe9e680SJoe Wallwork for (c = cStart; c < cEnd; ++c) { ierr = DMLabelGetValue(rgLabel, c, &cellTags[c]);CHKERRQ(ierr); } 1179fe9e680SJoe Wallwork } 1189fe9e680SJoe Wallwork 1192bc0b47dSJoe Wallwork /* Get metric */ 1202bc0b47dSJoe Wallwork ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr); 1212bc0b47dSJoe Wallwork ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr); 12293520041SJoe Wallwork ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr); 12393520041SJoe Wallwork ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr); 1242bc0b47dSJoe Wallwork for (v = 0; v < (vEnd-vStart); ++v) { 1252bc0b47dSJoe Wallwork for (i = 0, k = 0; i < dim; ++i) { 1262bc0b47dSJoe Wallwork for (j = i; j < dim; ++j) { 12793520041SJoe Wallwork if (isotropic) { 12893520041SJoe Wallwork if (i == j) { 12993520041SJoe Wallwork if (uniform) metric[Neq*v+k] = PetscRealPart(met[0]); 13093520041SJoe Wallwork else metric[Neq*v+k] = PetscRealPart(met[v]); 13193520041SJoe Wallwork } else metric[Neq*v+k] = 0.0; 13293520041SJoe Wallwork } else { 1332bc0b47dSJoe Wallwork metric[Neq*v+k] = PetscRealPart(met[dim*dim*v+dim*i+j]); 13493520041SJoe Wallwork } 1352bc0b47dSJoe Wallwork k++; 1362bc0b47dSJoe Wallwork } 1372bc0b47dSJoe Wallwork } 1382bc0b47dSJoe Wallwork } 1392bc0b47dSJoe Wallwork ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr); 1402bc0b47dSJoe Wallwork 1412bc0b47dSJoe Wallwork /* Send mesh to Mmg and remesh */ 1422bc0b47dSJoe Wallwork ierr = DMPlexMetricGetVerbosity(dm, &verbosity);CHKERRQ(ierr); 1432bc0b47dSJoe Wallwork ierr = DMPlexMetricGetGradationFactor(dm, &gradationFactor);CHKERRQ(ierr); 144*ae8b063eSJoe Wallwork ierr = DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber);CHKERRQ(ierr); 1452bc0b47dSJoe Wallwork ierr = DMPlexMetricNoInsertion(dm, &noInsert);CHKERRQ(ierr); 1462bc0b47dSJoe Wallwork ierr = DMPlexMetricNoSwapping(dm, &noSwap);CHKERRQ(ierr); 1472bc0b47dSJoe Wallwork ierr = DMPlexMetricNoMovement(dm, &noMove);CHKERRQ(ierr); 1482bc0b47dSJoe Wallwork switch (dim) { 1492bc0b47dSJoe Wallwork case 2: 1502bc0b47dSJoe Wallwork ierr = MMG2D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end); 1512bc0b47dSJoe Wallwork ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noinsert, noInsert); 1522bc0b47dSJoe Wallwork ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noswap, noSwap); 1532bc0b47dSJoe Wallwork ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nomove, noMove); 1542bc0b47dSJoe Wallwork ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_verbose, verbosity); 1552bc0b47dSJoe Wallwork ierr = MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hgrad, gradationFactor); 156*ae8b063eSJoe Wallwork ierr = MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hausd, hausdorffNumber); 1578d788f11SJoe Wallwork ierr = MMG2D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags); 1582bc0b47dSJoe Wallwork ierr = MMG2D_Set_vertices(mmg_mesh, vertices, verTags); 1592bc0b47dSJoe Wallwork ierr = MMG2D_Set_triangles(mmg_mesh, cells, cellTags); 1608d788f11SJoe Wallwork ierr = MMG2D_Set_edges(mmg_mesh, bdFaces, faceTags); 1612bc0b47dSJoe Wallwork ierr = MMG2D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor); 1622bc0b47dSJoe Wallwork ierr = MMG2D_Set_tensorSols(mmg_metric, metric); 1632bc0b47dSJoe Wallwork ierr = MMG2D_mmg2dlib(mmg_mesh, mmg_metric); 1642bc0b47dSJoe Wallwork break; 1652bc0b47dSJoe Wallwork case 3: 1662bc0b47dSJoe Wallwork ierr = MMG3D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end); 1672bc0b47dSJoe Wallwork ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noinsert, noInsert); 1682bc0b47dSJoe Wallwork ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noswap, noSwap); 1692bc0b47dSJoe Wallwork ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nomove, noMove); 1702bc0b47dSJoe Wallwork ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_verbose, verbosity); 1712bc0b47dSJoe Wallwork ierr = MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG3D_DPARAM_hgrad, gradationFactor); 172*ae8b063eSJoe Wallwork ierr = MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hausd, hausdorffNumber); 1738d788f11SJoe Wallwork ierr = MMG3D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags, 0, 0); 1742bc0b47dSJoe Wallwork ierr = MMG3D_Set_vertices(mmg_mesh, vertices, verTags); 1752bc0b47dSJoe Wallwork ierr = MMG3D_Set_tetrahedra(mmg_mesh, cells, cellTags); 1768d788f11SJoe Wallwork ierr = MMG3D_Set_triangles(mmg_mesh, bdFaces, faceTags); 1772bc0b47dSJoe Wallwork ierr = MMG3D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor); 1782bc0b47dSJoe Wallwork ierr = MMG3D_Set_tensorSols(mmg_metric, metric); 1792bc0b47dSJoe Wallwork ierr = MMG3D_mmg3dlib(mmg_mesh, mmg_metric); 1802bc0b47dSJoe Wallwork break; 18198921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %D", dim); 1822bc0b47dSJoe Wallwork } 1838d788f11SJoe Wallwork ierr = PetscFree(cells);CHKERRQ(ierr); 1848d788f11SJoe Wallwork ierr = PetscFree2(metric, vertices);CHKERRQ(ierr); 1858d788f11SJoe Wallwork ierr = PetscFree2(bdFaces, faceTags);CHKERRQ(ierr); 1868d788f11SJoe Wallwork ierr = PetscFree2(verTags, cellTags);CHKERRQ(ierr); 1872bc0b47dSJoe Wallwork 18899b4fe00SJoe Wallwork /* Retrieve mesh from Mmg */ 1892bc0b47dSJoe Wallwork switch (dim) { 1902bc0b47dSJoe Wallwork case 2: 1912bc0b47dSJoe Wallwork numCornersNew = 3; 1922bc0b47dSJoe Wallwork ierr = MMG2D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew); 1932bc0b47dSJoe Wallwork ierr = PetscMalloc4(2*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer);CHKERRQ(ierr); 1942bc0b47dSJoe Wallwork ierr = PetscMalloc3(3*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells);CHKERRQ(ierr); 1952bc0b47dSJoe Wallwork ierr = PetscMalloc4(2*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces);CHKERRQ(ierr); 1962bc0b47dSJoe Wallwork ierr = MMG2D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer); 1972bc0b47dSJoe Wallwork ierr = MMG2D_Get_triangles(mmg_mesh, cellsNew, cellTagsNew, requiredCells); 1982bc0b47dSJoe Wallwork ierr = MMG2D_Get_edges(mmg_mesh, facesNew, faceTagsNew, ridges, requiredFaces); 1992bc0b47dSJoe Wallwork break; 2002bc0b47dSJoe Wallwork case 3: 2012bc0b47dSJoe Wallwork numCornersNew = 4; 2022bc0b47dSJoe Wallwork ierr = MMG3D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0); 2032bc0b47dSJoe Wallwork ierr = PetscMalloc4(3*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer);CHKERRQ(ierr); 2042bc0b47dSJoe Wallwork ierr = PetscMalloc3(4*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells);CHKERRQ(ierr); 2052bc0b47dSJoe Wallwork ierr = PetscMalloc4(3*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces);CHKERRQ(ierr); 2062bc0b47dSJoe Wallwork ierr = MMG3D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer); 2072bc0b47dSJoe Wallwork ierr = MMG3D_Get_tetrahedra(mmg_mesh, cellsNew, cellTagsNew, requiredCells); 2082bc0b47dSJoe Wallwork ierr = MMG3D_Get_triangles(mmg_mesh, facesNew, faceTagsNew, requiredFaces); 20999b4fe00SJoe Wallwork 21099b4fe00SJoe Wallwork /* Reorder for consistency with DMPlex */ 21199b4fe00SJoe Wallwork for (i = 0; i < numCellsNew; ++i) { ierr = DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4*i]);CHKERRQ(ierr); } 2122bc0b47dSJoe Wallwork break; 21399b4fe00SJoe Wallwork 21498921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %D", dim); 2152bc0b47dSJoe Wallwork } 21699b4fe00SJoe Wallwork 21799b4fe00SJoe Wallwork /* Create new Plex */ 2182bc0b47dSJoe Wallwork for (i = 0; i < (dim+1)*numCellsNew; i++) cellsNew[i] -= 1; 2192bc0b47dSJoe Wallwork for (i = 0; i < dim*numFacesNew; i++) facesNew[i] -= 1; 2202bc0b47dSJoe Wallwork ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNew, NULL, NULL, dmNew);CHKERRQ(ierr); 2212bc0b47dSJoe Wallwork switch (dim) { 2222bc0b47dSJoe Wallwork case 2: 2232bc0b47dSJoe Wallwork ierr = MMG2D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end); 2242bc0b47dSJoe Wallwork break; 2252bc0b47dSJoe Wallwork case 3: 2262bc0b47dSJoe Wallwork ierr = MMG3D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end); 2272bc0b47dSJoe Wallwork break; 22898921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %D", dim); 2292bc0b47dSJoe Wallwork } 2308d788f11SJoe Wallwork ierr = PetscFree4(verticesNew, verTagsNew, corners, requiredVer);CHKERRQ(ierr); 2318d788f11SJoe Wallwork 2328d788f11SJoe Wallwork /* Get adapted mesh information */ 2338d788f11SJoe Wallwork ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr); 2348d788f11SJoe Wallwork ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 2358d788f11SJoe Wallwork ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr); 2362bc0b47dSJoe Wallwork 2372bc0b47dSJoe Wallwork /* Rebuild boundary labels */ 2385e7de161SJoe Wallwork ierr = DMCreateLabel(*dmNew, flg ? bdName : bdLabelName);CHKERRQ(ierr); 2395e7de161SJoe Wallwork ierr = DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew);CHKERRQ(ierr); 2402bc0b47dSJoe Wallwork for (i = 0; i < numFacesNew; i++) { 2412bc0b47dSJoe Wallwork PetscInt numCoveredPoints, numFaces = 0, facePoints[3]; 2422bc0b47dSJoe Wallwork const PetscInt *coveredPoints = NULL; 2432bc0b47dSJoe Wallwork 2442bc0b47dSJoe Wallwork for (j = 0; j < dim; ++j) facePoints[j] = facesNew[i*dim+j]+vStart; 2452bc0b47dSJoe Wallwork ierr = DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr); 2462bc0b47dSJoe Wallwork for (j = 0; j < numCoveredPoints; ++j) { 2472bc0b47dSJoe Wallwork if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) { 2482bc0b47dSJoe Wallwork numFaces++; 2492bc0b47dSJoe Wallwork f = j; 2502bc0b47dSJoe Wallwork } 2512bc0b47dSJoe Wallwork } 2522c71b3e2SJacob Faibussowitsch PetscCheckFalse(numFaces != 1,comm, PETSC_ERR_ARG_OUTOFRANGE, "%d vertices cannot define more than 1 facet (%d)", dim, numFaces); 2532bc0b47dSJoe Wallwork ierr = DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]);CHKERRQ(ierr); 2542bc0b47dSJoe Wallwork ierr = DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr); 2552bc0b47dSJoe Wallwork } 2568d788f11SJoe Wallwork ierr = PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces);CHKERRQ(ierr); 2572bc0b47dSJoe Wallwork 2589fe9e680SJoe Wallwork /* Rebuild cell labels */ 2599fe9e680SJoe Wallwork ierr = DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName);CHKERRQ(ierr); 2609fe9e680SJoe Wallwork ierr = DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew);CHKERRQ(ierr); 2619fe9e680SJoe Wallwork for (c = cStart; c < cEnd; ++c) { ierr = DMLabelSetValue(rgLabelNew, c, cellTagsNew[c-cStart]);CHKERRQ(ierr); } 2622bc0b47dSJoe Wallwork ierr = PetscFree3(cellsNew, cellTagsNew, requiredCells);CHKERRQ(ierr); 2638d788f11SJoe Wallwork 2642bc0b47dSJoe Wallwork PetscFunctionReturn(0); 2652bc0b47dSJoe Wallwork } 266