1 #include "../mmgcommon.h" /*I "petscdmplex.h" I*/ 2 #include <mmg/libmmg.h> 3 4 PetscBool MmgCite = PETSC_FALSE; 5 const char MmgCitation[] = "@article{DAPOGNY2014358,\n" 6 " title = {Three-dimensional adaptive domain remeshing, implicit domain meshing, and applications to free and moving boundary problems},\n" 7 " journal = {Journal of Computational Physics},\n" 8 " author = {C. Dapogny and C. Dobrzynski and P. Frey},\n" 9 " volume = {262},\n" 10 " pages = {358--378},\n" 11 " doi = {10.1016/j.jcp.2014.01.005},\n" 12 " year = {2014}\n}\n"; 13 14 PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew) { 15 MPI_Comm comm; 16 const char *bdName = "_boundary_"; 17 const char *rgName = "_regions_"; 18 DM udm, cdm; 19 DMLabel bdLabelNew, rgLabelNew; 20 const char *bdLabelName, *rgLabelName; 21 PetscSection coordSection; 22 Vec coordinates; 23 const PetscScalar *coords, *met; 24 PetscReal *vertices, *metric, *verticesNew, gradationFactor, hausdorffNumber; 25 PetscInt *cells, *cellsNew, *cellTags, *cellTagsNew, *verTags, *verTagsNew; 26 PetscInt *bdFaces, *faceTags, *facesNew, *faceTagsNew; 27 PetscInt *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces; 28 PetscInt cStart, cEnd, c, numCells, fStart, fEnd, numFaceTags, f, vStart, vEnd, v, numVertices; 29 PetscInt dim, off, coff, maxConeSize, bdSize, i, j, k, Neq, verbosity, pStart, pEnd; 30 PetscInt numCellsNew, numVerticesNew, numCornersNew, numFacesNew; 31 PetscBool flg = PETSC_FALSE, noInsert, noSwap, noMove, noSurf, isotropic, uniform; 32 MMG5_pMesh mmg_mesh = NULL; 33 MMG5_pSol mmg_metric = NULL; 34 35 PetscFunctionBegin; 36 PetscCall(PetscCitationsRegister(MmgCitation, &MmgCite)); 37 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 38 if (bdLabel) { 39 PetscCall(PetscObjectGetName((PetscObject)bdLabel, &bdLabelName)); 40 PetscCall(PetscStrcmp(bdLabelName, bdName, &flg)); 41 PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 42 } 43 if (rgLabel) { 44 PetscCall(PetscObjectGetName((PetscObject)rgLabel, &rgLabelName)); 45 PetscCall(PetscStrcmp(rgLabelName, rgName, &flg)); 46 PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName); 47 } 48 49 /* Get mesh information */ 50 PetscCall(DMGetDimension(dm, &dim)); 51 Neq = (dim * (dim + 1)) / 2; 52 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 53 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 54 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 55 PetscCall(DMPlexUninterpolate(dm, &udm)); 56 PetscCall(DMPlexGetMaxSizes(udm, &maxConeSize, NULL)); 57 numCells = cEnd - cStart; 58 numVertices = vEnd - vStart; 59 60 /* Get cell offsets */ 61 PetscCall(PetscMalloc1(numCells * maxConeSize, &cells)); 62 for (c = 0, coff = 0; c < numCells; ++c) { 63 const PetscInt *cone; 64 PetscInt coneSize, cl; 65 66 PetscCall(DMPlexGetConeSize(udm, c, &coneSize)); 67 PetscCall(DMPlexGetCone(udm, c, &cone)); 68 for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart + 1; 69 } 70 71 /* Get vertex coordinate array */ 72 PetscCall(DMGetCoordinateDM(dm, &cdm)); 73 PetscCall(DMGetLocalSection(cdm, &coordSection)); 74 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 75 PetscCall(VecGetArrayRead(coordinates, &coords)); 76 PetscCall(PetscMalloc2(numVertices * Neq, &metric, dim * numVertices, &vertices)); 77 for (v = 0; v < vEnd - vStart; ++v) { 78 PetscCall(PetscSectionGetOffset(coordSection, v + vStart, &off)); 79 for (i = 0; i < dim; ++i) vertices[dim * v + i] = PetscRealPart(coords[off + i]); 80 } 81 PetscCall(VecRestoreArrayRead(coordinates, &coords)); 82 PetscCall(DMDestroy(&udm)); 83 84 /* Get face tags */ 85 if (!bdLabel) { 86 flg = PETSC_TRUE; 87 PetscCall(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabel)); 88 PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel)); 89 } 90 PetscCall(DMLabelGetBounds(bdLabel, &pStart, &pEnd)); 91 for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) { 92 PetscBool hasPoint; 93 PetscInt *closure = NULL, closureSize, cl; 94 95 PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint)); 96 if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue; 97 numFaceTags++; 98 99 PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 100 for (cl = 0; cl < closureSize * 2; cl += 2) { 101 if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 102 } 103 PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 104 } 105 PetscCall(PetscMalloc2(bdSize, &bdFaces, numFaceTags, &faceTags)); 106 for (f = pStart, bdSize = 0, numFaceTags = 0; f < pEnd; ++f) { 107 PetscBool hasPoint; 108 PetscInt *closure = NULL, closureSize, cl; 109 110 PetscCall(DMLabelHasPoint(bdLabel, f, &hasPoint)); 111 if ((!hasPoint) || (f < fStart) || (f >= fEnd)) continue; 112 113 PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 114 for (cl = 0; cl < closureSize * 2; cl += 2) { 115 if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart + 1; 116 } 117 PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure)); 118 PetscCall(DMLabelGetValue(bdLabel, f, &faceTags[numFaceTags++])); 119 } 120 if (flg) PetscCall(DMLabelDestroy(&bdLabel)); 121 122 /* Get cell tags */ 123 PetscCall(PetscCalloc2(numVertices, &verTags, numCells, &cellTags)); 124 if (rgLabel) { 125 for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelGetValue(rgLabel, c, &cellTags[c])); 126 } 127 128 /* Get metric */ 129 PetscCall(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view")); 130 PetscCall(VecGetArrayRead(vertexMetric, &met)); 131 PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic)); 132 PetscCall(DMPlexMetricIsUniform(dm, &uniform)); 133 for (v = 0; v < (vEnd - vStart); ++v) { 134 for (i = 0, k = 0; i < dim; ++i) { 135 for (j = i; j < dim; ++j) { 136 if (isotropic) { 137 if (i == j) { 138 if (uniform) metric[Neq * v + k] = PetscRealPart(met[0]); 139 else metric[Neq * v + k] = PetscRealPart(met[v]); 140 } else metric[Neq * v + k] = 0.0; 141 } else { 142 metric[Neq * v + k] = PetscRealPart(met[dim * dim * v + dim * i + j]); 143 } 144 k++; 145 } 146 } 147 } 148 PetscCall(VecRestoreArrayRead(vertexMetric, &met)); 149 150 /* Send mesh to Mmg and remesh */ 151 PetscCall(DMPlexMetricGetVerbosity(dm, &verbosity)); 152 PetscCall(DMPlexMetricGetGradationFactor(dm, &gradationFactor)); 153 PetscCall(DMPlexMetricGetHausdorffNumber(dm, &hausdorffNumber)); 154 PetscCall(DMPlexMetricNoInsertion(dm, &noInsert)); 155 PetscCall(DMPlexMetricNoSwapping(dm, &noSwap)); 156 PetscCall(DMPlexMetricNoMovement(dm, &noMove)); 157 PetscCall(DMPlexMetricNoSurf(dm, &noSurf)); 158 switch (dim) { 159 case 2: 160 PetscCallMMG_NONSTANDARD(MMG2D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); 161 PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noinsert, noInsert)); 162 PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noswap, noSwap)); 163 PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nomove, noMove)); 164 PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nosurf, noSurf)); 165 PetscCallMMG_NONSTANDARD(MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_verbose, verbosity)); 166 PetscCallMMG_NONSTANDARD(MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hgrad, gradationFactor)); 167 PetscCallMMG_NONSTANDARD(MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hausd, hausdorffNumber)); 168 PetscCallMMG_NONSTANDARD(MMG2D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags)); 169 PetscCallMMG_NONSTANDARD(MMG2D_Set_vertices(mmg_mesh, vertices, verTags)); 170 PetscCallMMG_NONSTANDARD(MMG2D_Set_triangles(mmg_mesh, cells, cellTags)); 171 PetscCallMMG_NONSTANDARD(MMG2D_Set_edges(mmg_mesh, bdFaces, faceTags)); 172 PetscCallMMG_NONSTANDARD(MMG2D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor)); 173 PetscCallMMG_NONSTANDARD(MMG2D_Set_tensorSols(mmg_metric, metric)); 174 PetscCallMMG(MMG2D_mmg2dlib(mmg_mesh, mmg_metric)); 175 break; 176 case 3: 177 PetscCallMMG_NONSTANDARD(MMG3D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); 178 PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noinsert, noInsert)); 179 PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noswap, noSwap)); 180 PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nomove, noMove)); 181 PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nosurf, noSurf)); 182 PetscCallMMG_NONSTANDARD(MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_verbose, verbosity)); 183 PetscCallMMG_NONSTANDARD(MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG3D_DPARAM_hgrad, gradationFactor)); 184 PetscCallMMG_NONSTANDARD(MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG3D_DPARAM_hausd, hausdorffNumber)); 185 PetscCallMMG_NONSTANDARD(MMG3D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numFaceTags, 0, 0)); 186 PetscCallMMG_NONSTANDARD(MMG3D_Set_vertices(mmg_mesh, vertices, verTags)); 187 PetscCallMMG_NONSTANDARD(MMG3D_Set_tetrahedra(mmg_mesh, cells, cellTags)); 188 PetscCallMMG_NONSTANDARD(MMG3D_Set_triangles(mmg_mesh, bdFaces, faceTags)); 189 PetscCallMMG_NONSTANDARD(MMG3D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor)); 190 PetscCallMMG_NONSTANDARD(MMG3D_Set_tensorSols(mmg_metric, metric)); 191 PetscCallMMG(MMG3D_mmg3dlib(mmg_mesh, mmg_metric)); 192 break; 193 default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim); 194 } 195 PetscCall(PetscFree(cells)); 196 PetscCall(PetscFree2(metric, vertices)); 197 PetscCall(PetscFree2(bdFaces, faceTags)); 198 PetscCall(PetscFree2(verTags, cellTags)); 199 200 /* Retrieve mesh from Mmg */ 201 switch (dim) { 202 case 2: 203 numCornersNew = 3; 204 PetscCallMMG_NONSTANDARD(MMG2D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew)); 205 PetscCall(PetscMalloc4(2 * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer)); 206 PetscCall(PetscMalloc3(3 * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells)); 207 PetscCall(PetscMalloc4(2 * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces)); 208 PetscCallMMG_NONSTANDARD(MMG2D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer)); 209 PetscCallMMG_NONSTANDARD(MMG2D_Get_triangles(mmg_mesh, cellsNew, cellTagsNew, requiredCells)); 210 PetscCallMMG_NONSTANDARD(MMG2D_Get_edges(mmg_mesh, facesNew, faceTagsNew, ridges, requiredFaces)); 211 break; 212 case 3: 213 numCornersNew = 4; 214 PetscCallMMG_NONSTANDARD(MMG3D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0)); 215 PetscCall(PetscMalloc4(3 * numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer)); 216 PetscCall(PetscMalloc3(4 * numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells)); 217 PetscCall(PetscMalloc4(3 * numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces)); 218 PetscCallMMG_NONSTANDARD(MMG3D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer)); 219 PetscCallMMG_NONSTANDARD(MMG3D_Get_tetrahedra(mmg_mesh, cellsNew, cellTagsNew, requiredCells)); 220 PetscCallMMG_NONSTANDARD(MMG3D_Get_triangles(mmg_mesh, facesNew, faceTagsNew, requiredFaces)); 221 222 /* Reorder for consistency with DMPlex */ 223 for (i = 0; i < numCellsNew; ++i) PetscCall(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cellsNew[4 * i])); 224 break; 225 226 default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim); 227 } 228 229 /* Create new Plex */ 230 for (i = 0; i < (dim + 1) * numCellsNew; i++) cellsNew[i] -= 1; 231 for (i = 0; i < dim * numFacesNew; i++) facesNew[i] -= 1; 232 PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNew, NULL, NULL, dmNew)); 233 switch (dim) { 234 case 2: PetscCallMMG_NONSTANDARD(MMG2D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); break; 235 case 3: PetscCallMMG_NONSTANDARD(MMG3D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end)); break; 236 default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %" PetscInt_FMT, dim); 237 } 238 PetscCall(PetscFree4(verticesNew, verTagsNew, corners, requiredVer)); 239 240 /* Get adapted mesh information */ 241 PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd)); 242 PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd)); 243 PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd)); 244 245 /* Rebuild boundary labels */ 246 PetscCall(DMCreateLabel(*dmNew, flg ? bdName : bdLabelName)); 247 PetscCall(DMGetLabel(*dmNew, flg ? bdName : bdLabelName, &bdLabelNew)); 248 for (i = 0; i < numFacesNew; i++) { 249 PetscInt numCoveredPoints, numFaces = 0, facePoints[3]; 250 const PetscInt *coveredPoints = NULL; 251 252 for (j = 0; j < dim; ++j) facePoints[j] = facesNew[i * dim + j] + vStart; 253 PetscCall(DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints)); 254 for (j = 0; j < numCoveredPoints; ++j) { 255 if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) { 256 numFaces++; 257 f = j; 258 } 259 } 260 PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT " vertices cannot define more than 1 facet (%" PetscInt_FMT ")", dim, numFaces); 261 PetscCall(DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i])); 262 PetscCall(DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints)); 263 } 264 PetscCall(PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces)); 265 266 /* Rebuild cell labels */ 267 PetscCall(DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName)); 268 PetscCall(DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew)); 269 for (c = cStart; c < cEnd; ++c) PetscCall(DMLabelSetValue(rgLabelNew, c, cellTagsNew[c - cStart])); 270 PetscCall(PetscFree3(cellsNew, cellTagsNew, requiredCells)); 271 272 PetscFunctionReturn(0); 273 } 274