xref: /petsc/src/dm/impls/plex/adaptors/mmg/mmgadapt.c (revision 8a2871f65f293ab823050fa896a854ffce5d8a35)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <mmg/libmmg.h>
3 
4 PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew)
5 {
6   MPI_Comm           comm;
7   const char        *bdName = "_boundary_";
8   const char        *rgName = "_regions_";
9   DM                 udm, cdm;
10   DMLabel            bdLabelFull;
11   const char        *bdLabelName, *rgLabelName;
12   IS                 bdIS;
13   PetscSection       coordSection;
14   Vec                coordinates;
15   const PetscScalar *coords, *met;
16   const PetscInt    *bdFacesFull;
17   PetscInt          *bdFaces, *bdFaceIds;
18   PetscReal         *vertices, *metric, *verticesNew, gradationFactor;
19   PetscInt          *cells, *cellsNew;
20   PetscInt          *facesNew, *faceTagsNew;
21   PetscInt          *verTags, *cellTags, *verTagsNew, *cellTagsNew;
22   PetscInt          *corners, *requiredCells, *requiredVer, *ridges, *requiredFaces;
23   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, v;
24   PetscInt           off, maxConeSize, numBdFaces, f, bdSize, fStart, fEnd, i, j, k, Neq;
25   PetscInt           numCellsNew, numVerticesNew, numCornersNew, numFacesNew;
26   PetscInt           verbosity;
27   PetscBool          flg, noInsert, noSwap, noMove;
28   DMLabel            bdLabelNew, rgLabelNew;
29   MMG5_pMesh         mmg_mesh = NULL;
30   MMG5_pSol          mmg_metric = NULL;
31   PetscErrorCode     ierr;
32 
33   PetscFunctionBegin;
34   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
35   if (bdLabel) {
36     ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr);
37     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
38     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
39   }
40   if (rgLabel) {
41     ierr = PetscObjectGetName((PetscObject) rgLabel, &rgLabelName);CHKERRQ(ierr);
42     ierr = PetscStrcmp(rgLabelName, rgName, &flg);CHKERRQ(ierr);
43     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for element tags", rgLabelName);
44   }
45 
46   /* Get mesh information */
47   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
48   Neq  = (dim*(dim+1))/2;
49   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
50   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
51   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
52   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
53   numCells    = cEnd - cStart;
54   numVertices = vEnd - vStart;
55 
56   /* Get cell offsets */
57   ierr = PetscMalloc1(numCells*maxConeSize, &cells);CHKERRQ(ierr);
58   for (c = 0, coff = 0; c < numCells; ++c) {
59     const PetscInt *cone;
60     PetscInt        coneSize, cl;
61 
62     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
63     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
64     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart + 1;
65   }
66 
67   /* Get vertex coordinate array */
68   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
69   ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
70   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
71   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
72   ierr = PetscMalloc2(numVertices*Neq, &metric, dim*numVertices, &vertices);CHKERRQ(ierr);
73   for (v = 0; v < vEnd-vStart; ++v) {
74     ierr = PetscSectionGetOffset(coordSection, v+vStart, &off);CHKERRQ(ierr);
75     for (i = 0; i < dim; ++i) vertices[dim*v+i] = PetscRealPart(coords[off+i]);
76   }
77   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
78   ierr = DMDestroy(&udm);CHKERRQ(ierr);
79 
80   /* Get boundary mesh */
81   ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull);CHKERRQ(ierr);
82   ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull);CHKERRQ(ierr);
83   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
84   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
85   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
86   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
87     PetscInt *closure = NULL;
88     PetscInt  closureSize, cl;
89 
90     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
91     for (cl = 0; cl < closureSize*2; cl += 2) {
92       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
93     }
94     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
95   }
96   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
97   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
98     PetscInt *closure = NULL;
99     PetscInt  closureSize, cl;
100 
101     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
102     for (cl = 0; cl < closureSize*2; cl += 2) {
103       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart + 1;
104     }
105     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
106     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
107     else         {bdFaceIds[f] = 1;}
108   }
109   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
110   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
111 
112   /* Get cell tags */
113   ierr = PetscCalloc2(numVertices, &verTags, numCells, &cellTags);CHKERRQ(ierr);
114   if (rgLabel) {
115     for (c = cStart; c < cEnd; ++c) { ierr = DMLabelGetValue(rgLabel, c, &cellTags[c]);CHKERRQ(ierr); }
116   }
117 
118   /* Get metric */
119   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
120   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
121   for (v = 0; v < (vEnd-vStart); ++v) {
122     for (i = 0, k = 0; i < dim; ++i) {
123       for (j = i; j < dim; ++j) {
124         metric[Neq*v+k] = PetscRealPart(met[dim*dim*v+dim*i+j]);
125         k++;
126       }
127     }
128   }
129   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
130 
131   /* Send mesh to Mmg and remesh */
132   ierr = DMPlexMetricGetVerbosity(dm, &verbosity);CHKERRQ(ierr);
133   ierr = DMPlexMetricGetGradationFactor(dm, &gradationFactor);CHKERRQ(ierr);
134   ierr = DMPlexMetricNoInsertion(dm, &noInsert);CHKERRQ(ierr);
135   ierr = DMPlexMetricNoSwapping(dm, &noSwap);CHKERRQ(ierr);
136   ierr = DMPlexMetricNoMovement(dm, &noMove);CHKERRQ(ierr);
137   switch (dim) {
138   case 2:
139     ierr = MMG2D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
140     ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noinsert, noInsert);
141     ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_noswap, noSwap);
142     ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_nomove, noMove);
143     ierr = MMG2D_Set_iparameter(mmg_mesh, mmg_metric, MMG2D_IPARAM_verbose, verbosity);
144     ierr = MMG2D_Set_dparameter(mmg_mesh, mmg_metric, MMG2D_DPARAM_hgrad, gradationFactor);
145     ierr = MMG2D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numBdFaces);
146     ierr = MMG2D_Set_vertices(mmg_mesh, vertices, verTags);
147     ierr = MMG2D_Set_triangles(mmg_mesh, cells, cellTags);
148     ierr = MMG2D_Set_edges(mmg_mesh, bdFaces, bdFaceIds);
149     ierr = MMG2D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor);
150     ierr = MMG2D_Set_tensorSols(mmg_metric, metric);
151     ierr = MMG2D_mmg2dlib(mmg_mesh, mmg_metric);
152     break;
153   case 3:
154     ierr = MMG3D_Init_mesh(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
155     ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noinsert, noInsert);
156     ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_noswap, noSwap);
157     ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_nomove, noMove);
158     ierr = MMG3D_Set_iparameter(mmg_mesh, mmg_metric, MMG3D_IPARAM_verbose, verbosity);
159     ierr = MMG3D_Set_dparameter(mmg_mesh, mmg_metric, MMG3D_DPARAM_hgrad, gradationFactor);
160     ierr = MMG3D_Set_meshSize(mmg_mesh, numVertices, numCells, 0, numBdFaces, 0, 0);
161     ierr = MMG3D_Set_vertices(mmg_mesh, vertices, verTags);
162     ierr = MMG3D_Set_tetrahedra(mmg_mesh, cells, cellTags);
163     ierr = MMG3D_Set_triangles(mmg_mesh, bdFaces, bdFaceIds);
164     ierr = MMG3D_Set_solSize(mmg_mesh, mmg_metric, MMG5_Vertex, numVertices, MMG5_Tensor);
165     ierr = MMG3D_Set_tensorSols(mmg_metric, metric);
166     ierr = MMG3D_mmg3dlib(mmg_mesh, mmg_metric);
167     break;
168   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %D", dim);
169   }
170 
171   /* Retrieve mesh from Mmg and create new Plex*/
172   switch (dim) {
173   case 2:
174     numCornersNew = 3;
175     ierr = MMG2D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew);
176     ierr = PetscMalloc4(2*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer);CHKERRQ(ierr);
177     ierr = PetscMalloc3(3*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells);CHKERRQ(ierr);
178     ierr = PetscMalloc4(2*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces);CHKERRQ(ierr);
179     ierr = MMG2D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer);
180     ierr = MMG2D_Get_triangles(mmg_mesh, cellsNew, cellTagsNew, requiredCells);
181     ierr = MMG2D_Get_edges(mmg_mesh, facesNew, faceTagsNew, ridges, requiredFaces);
182     break;
183   case 3:
184     numCornersNew = 4;
185     ierr = MMG3D_Get_meshSize(mmg_mesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0);
186     ierr = PetscMalloc4(3*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer);CHKERRQ(ierr);
187     ierr = PetscMalloc3(4*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells);CHKERRQ(ierr);
188     ierr = PetscMalloc4(3*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces);CHKERRQ(ierr);
189     ierr = MMG3D_Get_vertices(mmg_mesh, verticesNew, verTagsNew, corners, requiredVer);
190     ierr = MMG3D_Get_tetrahedra(mmg_mesh, cellsNew, cellTagsNew, requiredCells);
191     ierr = MMG3D_Get_triangles(mmg_mesh, facesNew, faceTagsNew, requiredFaces);
192     break;
193   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %D", dim);
194   }
195   for (i = 0; i < (dim+1)*numCellsNew; i++) cellsNew[i] -= 1;
196   for (i = 0; i < dim*numFacesNew; i++) facesNew[i] -= 1;
197   ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNew, NULL, NULL, dmNew);CHKERRQ(ierr);
198   switch (dim) {
199   case 2:
200     ierr = MMG2D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
201     break;
202   case 3:
203     ierr = MMG3D_Free_all(MMG5_ARG_start, MMG5_ARG_ppMesh, &mmg_mesh, MMG5_ARG_ppMet, &mmg_metric, MMG5_ARG_end);
204     break;
205   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Mmg adaptation defined for dimension %D", dim);
206   }
207 
208   /* Rebuild boundary labels */
209   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
210   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
211   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
212   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
213   for (i = 0; i < numFacesNew; i++) {
214     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
215     const PetscInt *coveredPoints = NULL;
216 
217     for (j = 0; j < dim; ++j) facePoints[j] = facesNew[i*dim+j]+vStart;
218     ierr = DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
219     for (j = 0; j < numCoveredPoints; ++j) {
220       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
221         numFaces++;
222         f = j;
223       }
224     }
225     if (numFaces != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "%d vertices cannot define more than 1 facet (%d)", dim, numFaces);
226     ierr = DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]);CHKERRQ(ierr);
227     ierr = DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
228   }
229 
230   /* Rebuild cell labels */
231   ierr = DMCreateLabel(*dmNew, rgLabel ? rgLabelName : rgName);CHKERRQ(ierr);
232   ierr = DMGetLabel(*dmNew, rgLabel ? rgLabelName : rgName, &rgLabelNew);CHKERRQ(ierr);
233   for (c = cStart; c < cEnd; ++c) { ierr = DMLabelSetValue(rgLabelNew, c, cellTagsNew[c-cStart]);CHKERRQ(ierr); }
234 
235   /* Clean up */
236   ierr = PetscFree(cells);CHKERRQ(ierr);
237   ierr = PetscFree2(metric, vertices);CHKERRQ(ierr);
238   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
239   ierr = PetscFree2(verTags, cellTags);CHKERRQ(ierr);
240   ierr = PetscFree4(verticesNew, verTagsNew, corners, requiredVer);CHKERRQ(ierr);
241   ierr = PetscFree3(cellsNew, cellTagsNew, requiredCells);CHKERRQ(ierr);
242   ierr = PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces);CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245