1*5f80ce2aSJacob Faibussowitsch #include "../mmgcommon.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; 15ae8b063eSJoe 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; 2276f360caSJoe Wallwork PetscBool flg = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform; 232bc0b47dSJoe Wallwork MMG5_pMesh mmg_mesh = NULL; 242bc0b47dSJoe Wallwork MMG5_pSol mmg_metric = NULL; 252bc0b47dSJoe Wallwork 262bc0b47dSJoe Wallwork PetscFunctionBegin; 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 282bc0b47dSJoe Wallwork if (bdLabel) { 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) bdLabel, &bdLabelName)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(bdLabelName, bdName, &flg)); 312c71b3e2SJacob Faibussowitsch PetscCheckFalse(flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 322bc0b47dSJoe Wallwork } 339fe9e680SJoe Wallwork if (rgLabel) { 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) rgLabel, &rgLabelName)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(rgLabelName, rgName, &flg)); 362c71b3e2SJacob Faibussowitsch PetscCheckFalse(flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName); 379fe9e680SJoe Wallwork } 382bc0b47dSJoe Wallwork 392bc0b47dSJoe Wallwork /* Get mesh information */ 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 412bc0b47dSJoe Wallwork Neq = (dim*(dim+1))/2; 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexUninterpolate(dm, &udm)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetMaxSizes(udm, &maxConeSize, NULL)); 472bc0b47dSJoe Wallwork numCells = cEnd - cStart; 482bc0b47dSJoe Wallwork numVertices = vEnd - vStart; 492bc0b47dSJoe Wallwork 502bc0b47dSJoe Wallwork /* Get cell offsets */ 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numCells*maxConeSize, &cells)); 522bc0b47dSJoe Wallwork for (c = 0, coff = 0; c < numCells; ++c) { 532bc0b47dSJoe Wallwork const PetscInt *cone; 542bc0b47dSJoe Wallwork PetscInt coneSize, cl; 552bc0b47dSJoe Wallwork 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetConeSize(udm, c, &coneSize)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(udm, c, &cone)); 582bc0b47dSJoe Wallwork for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart + 1; 592bc0b47dSJoe Wallwork } 602bc0b47dSJoe Wallwork 612bc0b47dSJoe Wallwork /* Get vertex coordinate array */ 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &cdm)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(cdm, &coordSection)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(coordinates, &coords)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(numVertices*Neq, &metric, dim*numVertices, &vertices)); 672bc0b47dSJoe Wallwork for (v = 0; v < vEnd-vStart; ++v) { 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(coordSection, v+vStart, &off)); 692bc0b47dSJoe Wallwork for (i = 0; i < dim; ++i) vertices[dim*v+i] = PetscRealPart(coords[off+i]); 702bc0b47dSJoe Wallwork } 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(coordinates, &coords)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&udm)); 732bc0b47dSJoe Wallwork 748d788f11SJoe Wallwork /* Get face tags */ 758d788f11SJoe Wallwork if (!bdLabel) { 768d788f11SJoe Wallwork flg = PETSC_TRUE; 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMarkBoundaryFaces(dm, 1, bdLabel)); 798d788f11SJoe Wallwork } 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetBounds(bdLabel, &pStart, &pEnd)); 818d788f11SJoe Wallwork for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) { 828d788f11SJoe Wallwork PetscBool hasPoint; 838d788f11SJoe Wallwork PetscInt *closure = NULL, closureSize, cl; 842bc0b47dSJoe Wallwork 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelHasPoint(bdLabel, f, &hasPoint)); 868d788f11SJoe Wallwork if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue; 878d788f11SJoe Wallwork numFaceTags++; 888d788f11SJoe Wallwork 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 902bc0b47dSJoe Wallwork for (cl = 0; cl < closureSize*2; cl += 2) { 912bc0b47dSJoe Wallwork if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 922bc0b47dSJoe Wallwork } 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 942bc0b47dSJoe Wallwork } 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags)); 968d788f11SJoe Wallwork for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) { 978d788f11SJoe Wallwork PetscBool hasPoint; 988d788f11SJoe Wallwork PetscInt *closure = NULL, closureSize, cl; 992bc0b47dSJoe Wallwork 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelHasPoint(bdLabel, f, &hasPoint)); 1018d788f11SJoe Wallwork if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue; 1028d788f11SJoe Wallwork 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 1042bc0b47dSJoe Wallwork for (cl = 0; cl < closureSize*2; cl += 2) { 1052bc0b47dSJoe Wallwork if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart + 1; 1062bc0b47dSJoe Wallwork } 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++])); 1092bc0b47dSJoe Wallwork } 110*5f80ce2aSJacob Faibussowitsch if (flg) CHKERRQ(DMLabelDestroy(&bdLabel)); 1112bc0b47dSJoe Wallwork 1129fe9e680SJoe Wallwork /* Get cell tags */ 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc2(numVertices, &verTags, numCells, &cellTags)); 1149fe9e680SJoe Wallwork if (rgLabel) { 115*5f80ce2aSJacob Faibussowitsch for (c = cStart; c < cEnd; ++c) CHKERRQ(DMLabelGetValue(rgLabel, c, &cellTags[c])); 1169fe9e680SJoe Wallwork } 1179fe9e680SJoe Wallwork 1182bc0b47dSJoe Wallwork /* Get metric */ 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view")); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(vertexMetric, &met)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricIsUniform(dm, &uniform)); 1232bc0b47dSJoe Wallwork for (v = 0; v < (vEnd-vStart); ++v) { 1242bc0b47dSJoe Wallwork for (i = 0, k = 0; i < dim; ++i) { 1252bc0b47dSJoe Wallwork for (j = i; j < dim; ++j) { 12693520041SJoe Wallwork if (isotropic) { 12793520041SJoe Wallwork if (i == j) { 12893520041SJoe Wallwork if (uniform) metric[Neq*v+k] = PetscRealPart(met[0]); 12993520041SJoe Wallwork else metric[Neq*v+k] = PetscRealPart(met[v]); 13093520041SJoe Wallwork } else metric[Neq*v+k] = 0.0; 13193520041SJoe Wallwork } else { 1322bc0b47dSJoe Wallwork metric[Neq*v+k] = PetscRealPart(met[dim*dim*v+dim*i+j]); 13393520041SJoe Wallwork } 1342bc0b47dSJoe Wallwork k++; 1352bc0b47dSJoe Wallwork } 1362bc0b47dSJoe Wallwork } 1372bc0b47dSJoe Wallwork } 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(vertexMetric, &met)); 1392bc0b47dSJoe Wallwork 1402bc0b47dSJoe Wallwork /* Send mesh to Mmg and remesh */ 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricGetVerbosity(dm, &verbosity)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricGetGradationFactor(dm, &gradationFactor)); 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricNoInsertion(dm, &noInsert)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricNoSwapping(dm, &noSwap)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricNoMovement(dm, &noMove)); 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricNoSurf(dm, &noSurf)); 1482bc0b47dSJoe Wallwork switch (dim) { 1492bc0b47dSJoe Wallwork case 2: 150*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); 151*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noinsert, noInsert)); 152*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noswap, noSwap)); 153*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nomove, noMove)); 154*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nosurf, noSurf)); 155*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_verbose, verbosity)); 156*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hgrad, gradationFactor)); 157*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hausd, hausdorffNumber)); 158*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags)); 159*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_vertices(mmg_mesh, vertices, verTags)); 160*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_triangles(mmg_mesh, cells, cellTags)); 161*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_edges(mmg_mesh, bdFaces, faceTags)); 162*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor)); 163*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Set_tensorSols(mmg_metric, metric)); 164*5f80ce2aSJacob Faibussowitsch CHKERRMMG(MMG2D_mmg2dlib(mmg_mesh, mmg_metric)); 1652bc0b47dSJoe Wallwork break; 1662bc0b47dSJoe Wallwork case 3: 167*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); 168*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noinsert, noInsert)); 169*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noswap, noSwap)); 170*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nomove, noMove)); 171*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nosurf, noSurf)); 172*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_verbose, verbosity)); 173*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG3D_DPARAM_hgrad, gradationFactor)); 174*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hausd, hausdorffNumber)); 175*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags, 0, 0)); 176*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_vertices(mmg_mesh, vertices, verTags)); 177*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_tetrahedra(mmg_mesh, cells, cellTags)); 178*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_triangles(mmg_mesh, bdFaces, faceTags)); 179*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor)); 180*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Set_tensorSols(mmg_metric, metric)); 181*5f80ce2aSJacob Faibussowitsch CHKERRMMG(MMG3D_mmg3dlib(mmg_mesh, mmg_metric)); 1822bc0b47dSJoe Wallwork break; 183*5f80ce2aSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim); 1842bc0b47dSJoe Wallwork } 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cells)); 186*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(metric, vertices)); 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(bdFaces, faceTags)); 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(verTags, cellTags)); 1892bc0b47dSJoe Wallwork 19099b4fe00SJoe Wallwork /* Retrieve mesh from Mmg */ 1912bc0b47dSJoe Wallwork switch (dim) { 1922bc0b47dSJoe Wallwork case 2: 1932bc0b47dSJoe Wallwork numCornersNew = 3; 194*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Get_meshSize, mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(2*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(3*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells)); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(2*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces)); 198*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Get_vertices, mmg_mesh, verticesNew, verTagsNew, corners, requiredVer); 199*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Get_triangles, mmg_mesh, cellsNew, cellTagsNew, requiredCells); 200*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Get_edges, mmg_mesh, facesNew, faceTagsNew, ridges, requiredFaces); 2012bc0b47dSJoe Wallwork break; 2022bc0b47dSJoe Wallwork case 3: 2032bc0b47dSJoe Wallwork numCornersNew = 4; 204*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Get_meshSize, mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0); 205*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(3*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer)); 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(4*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(3*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces)); 208*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Get_vertices, mmg_mesh, verticesNew, verTagsNew, corners, requiredVer); 209*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Get_tetrahedra, mmg_mesh, cellsNew, cellTagsNew, requiredCells); 210*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Get_triangles, mmg_mesh, facesNew, faceTagsNew, requiredFaces); 21199b4fe00SJoe Wallwork 21299b4fe00SJoe Wallwork /* Reorder for consistency with DMPlex */ 213*5f80ce2aSJacob Faibussowitsch for (i = 0; i < numCellsNew; ++i) CHKERRQ(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4*i])); 2142bc0b47dSJoe Wallwork break; 21599b4fe00SJoe Wallwork 216*5f80ce2aSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim); 2172bc0b47dSJoe Wallwork } 21899b4fe00SJoe Wallwork 21999b4fe00SJoe Wallwork /* Create new Plex */ 2202bc0b47dSJoe Wallwork for (i = 0; i < (dim+1)*numCellsNew; i++) cellsNew[i] -= 1; 2212bc0b47dSJoe Wallwork for (i = 0; i < dim*numFacesNew; i++) facesNew[i] -= 1; 222*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNew, NULL, NULL, dmNew)); 2232bc0b47dSJoe Wallwork switch (dim) { 2242bc0b47dSJoe Wallwork case 2: 225*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG2D_Free_all, MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end); 2262bc0b47dSJoe Wallwork break; 2272bc0b47dSJoe Wallwork case 3: 228*5f80ce2aSJacob Faibussowitsch CHKERRMMG_NONSTANDARD(MMG3D_Free_all, MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end); 2292bc0b47dSJoe Wallwork break; 230*5f80ce2aSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim); 2312bc0b47dSJoe Wallwork } 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree4(verticesNew, verTagsNew, corners, requiredVer)); 2338d788f11SJoe Wallwork 2348d788f11SJoe Wallwork /* Get adapted mesh information */ 235*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd)); 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd)); 237*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd)); 2382bc0b47dSJoe Wallwork 2392bc0b47dSJoe Wallwork /* Rebuild boundary labels */ 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName)); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew)); 2422bc0b47dSJoe Wallwork for (i = 0; i < numFacesNew; i++) { 2432bc0b47dSJoe Wallwork PetscInt numCoveredPoints, numFaces = 0, facePoints[3]; 2442bc0b47dSJoe Wallwork const PetscInt *coveredPoints = NULL; 2452bc0b47dSJoe Wallwork 2462bc0b47dSJoe Wallwork for (j = 0; j < dim; ++j) facePoints[j] = facesNew[i*dim+j]+vStart; 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints)); 2482bc0b47dSJoe Wallwork for (j = 0; j < numCoveredPoints; ++j) { 2492bc0b47dSJoe Wallwork if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) { 2502bc0b47dSJoe Wallwork numFaces++; 2512bc0b47dSJoe Wallwork f = j; 2522bc0b47dSJoe Wallwork } 2532bc0b47dSJoe Wallwork } 254*5f80ce2aSJacob Faibussowitsch PetscCheck(numFaces == 1,comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces); 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i])); 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints)); 2572bc0b47dSJoe Wallwork } 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces)); 2592bc0b47dSJoe Wallwork 2609fe9e680SJoe Wallwork /* Rebuild cell labels */ 261*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName)); 262*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew)); 263*5f80ce2aSJacob Faibussowitsch for (c = cStart; c < cEnd; ++c) CHKERRQ(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c-cStart])); 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(cellsNew, cellTagsNew, requiredCells)); 2658d788f11SJoe Wallwork 2662bc0b47dSJoe Wallwork PetscFunctionReturn(0); 2672bc0b47dSJoe Wallwork } 268