xref: /petsc/src/dm/impls/plex/adaptors/parmmg/parmmgadapt.c (revision e831869d00a6b6bf5bd64f1ce9662fca56500707)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <parmmg/libparmmg.h>
3 
4 PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew)
5 {
6   MPI_Comm           comm, tmpComm;
7   const char        *bdName = "_boundary_";
8   DM                 udm, cdm;
9   DMLabel            bdLabelFull;
10   const char        *bdLabelName;
11   IS                 bdIS, globalVertexNum;
12   PetscSection       coordSection;
13   Vec                coordinates;
14   PetscSF            sf;
15   const PetscScalar *coords, *met;
16   const PetscInt    *bdFacesFull;
17   PetscInt          *bdFaces, *bdFaceIds;
18   PetscReal         *vertices, *metric, *verticesNew;
19   PetscInt          *cells;
20   const PetscInt    *gV, *ioffset, *irootloc, *roffset, *rmine, *rremote;
21   PetscInt          *numVerInterfaces, *ngbRanks, *verNgbRank, *interfaces_lv, *interfaces_gv, *intOffset;
22   PetscInt           niranks, nrranks, numNgbRanks, numVerNgbRanksTotal, count, sliceSize;
23   PetscInt          *verTags, *cellTags;
24   PetscInt           dim, Neq, cStart, cEnd, numCells, coff, vStart, vEnd, numVertices;
25   PetscInt           maxConeSize, numBdFaces, bdSize, fStart, fEnd;
26   PetscBool          flg, noInsert, noSwap, noMove;
27   DMLabel            bdLabelNew;
28   PetscReal         *verticesNewLoc, gradationFactor;
29   PetscInt          *verTagsNew, *cellsNew, *cellTagsNew, *corners;
30   PetscInt          *requiredCells, *requiredVer, *facesNew, *faceTagsNew, *ridges, *requiredFaces;
31   PetscInt          *gv_new, *owners, *verticesNewSorted;
32   PetscInt           numCellsNew, numVerticesNew, numCornersNew, numFacesNew, numVerticesNewLoc;
33   PetscInt           i, j, k, p, r, n, v, c, f, off, lv, gv;
34   PetscInt           verbosity, numIter;
35   const PetscMPIInt *iranks, *rranks;
36   PetscMPIInt        numProcs, rank;
37   PMMG_pParMesh      parmesh = NULL;
38   PetscErrorCode     ierr;
39 
40   PetscFunctionBegin;
41   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
42   ierr = MPI_Comm_size(comm, &numProcs);CHKERRMPI(ierr);
43   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
44   ierr = MPI_Comm_dup(comm, &tmpComm);CHKERRMPI(ierr);
45   if (bdLabel) {
46     ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr);
47     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
48     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
49   }
50 
51   /* Get mesh information */
52   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
53   if (dim != 3) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "ParMmg only works in 3D.\n");
54   Neq  = (dim*(dim+1))/2;
55   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
56   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
57   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
58   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
59   numCells    = cEnd - cStart;
60   numVertices = vEnd - vStart;
61 
62   /* Get cell offsets */
63   ierr = PetscMalloc1(numCells*maxConeSize, &cells);CHKERRQ(ierr);
64   for (c = 0, coff = 0; c < numCells; ++c) {
65     const PetscInt *cone;
66     PetscInt        coneSize, cl;
67 
68     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
69     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
70     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart+1;
71   }
72 
73   /* Get vertex coordinate array */
74   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
75   ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
76   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
77   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
78   ierr = PetscMalloc2(numVertices*Neq, &metric, dim*numVertices, &vertices);CHKERRQ(ierr);
79   for (v = 0; v < vEnd-vStart; ++v) {
80     ierr = PetscSectionGetOffset(coordSection, v+vStart, &off);CHKERRQ(ierr);
81     for (i = 0; i < dim; ++i) vertices[dim*v+i] = PetscRealPart(coords[off+i]);
82   }
83   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
84 
85   /* Get boundary mesh */
86   ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull);CHKERRQ(ierr);
87   ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull);CHKERRQ(ierr);
88   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
89   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
90   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
91   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
92     PetscInt *closure = NULL;
93     PetscInt  closureSize, cl;
94 
95     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
96     for (cl = 0; cl < closureSize*2; cl += 2) {
97       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
98     }
99     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
100   }
101   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
102   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
103     PetscInt *closure = NULL;
104     PetscInt  closureSize, cl;
105 
106     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
107     for (cl = 0; cl < closureSize*2; cl += 2) {
108       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart+1;
109     }
110     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
111     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
112     else         {bdFaceIds[f] = 1;}
113   }
114   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
115   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
116 
117   /* Get metric */
118   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
119   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
120   for (v = 0; v < (vEnd-vStart); ++v) {
121     for (i = 0, k = 0; i < dim; ++i) {
122       for (j = i; j < dim; ++j) {
123         metric[Neq*v+k] = PetscRealPart(met[dim*dim*v+dim*i+j]);
124         k++;
125       }
126     }
127   }
128   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
129 
130   /* Build ParMMG communicators: the list of vertices between two partitions  */
131   niranks = nrranks = 0;
132   numNgbRanks = 0;
133   if (numProcs > 1) {
134     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
135     ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
136     ierr = PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc);CHKERRQ(ierr);
137     ierr = PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote);CHKERRQ(ierr);
138     ierr = PetscCalloc1(numProcs, &numVerInterfaces);CHKERRQ(ierr);
139 
140     /* Counting */
141     for (r = 0; r < niranks; ++r) {
142       for (i=ioffset[r], count=0; i<ioffset[r+1]; ++i) {
143         if (irootloc[i] >= vStart && irootloc[i] < vEnd) count++;
144       }
145       numVerInterfaces[iranks[r]] += count;
146     }
147     for (r = 0; r < nrranks; ++r) {
148       for (i=roffset[r], count=0; i<roffset[r+1]; ++i) {
149         if (rmine[i] >= vStart && rmine[i] < vEnd) count++;
150       }
151       numVerInterfaces[rranks[r]] += count;
152     }
153     for (p = 0; p < numProcs; ++p) {
154       if (numVerInterfaces[p]) numNgbRanks++;
155     }
156     ierr = PetscMalloc2(numNgbRanks, &ngbRanks, numNgbRanks, &verNgbRank);CHKERRQ(ierr);
157     for (p = 0, n = 0; p < numProcs; ++p) {
158       if (numVerInterfaces[p]) {
159         ngbRanks[n] = p;
160         verNgbRank[n] = numVerInterfaces[p];
161         n++;
162       }
163     }
164     numVerNgbRanksTotal = 0;
165     for (p = 0; p < numNgbRanks; ++p) numVerNgbRanksTotal += verNgbRank[p];
166 
167     /* For each neighbor, fill in interface arrays */
168     ierr = PetscMalloc3(numVerNgbRanksTotal, &interfaces_lv, numVerNgbRanksTotal, &interfaces_gv,  numNgbRanks+1, &intOffset);CHKERRQ(ierr);
169     intOffset[0] = 0;
170     for (p = 0, r = 0, i = 0; p < numNgbRanks; ++p) {
171       intOffset[p+1] = intOffset[p];
172       if (iranks && iranks[i] == ngbRanks[p]) {
173 
174         /* Add the right slice of irootloc at the right place */
175         sliceSize = ioffset[i+1]-ioffset[i];
176         for (j = 0, count = 0; j < sliceSize; ++j) {
177           if (ioffset[i]+j >= ioffset[niranks]) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
178           v = irootloc[ioffset[i]+j];
179           if (v >= vStart && v < vEnd) {
180             if (intOffset[p+1]+count >= numVerNgbRanksTotal) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
181             interfaces_lv[intOffset[p+1]+count] = v-vStart;
182             count++;
183           }
184         }
185         intOffset[p+1] += count;
186         i++;
187       }
188       if (rranks && rranks[r] == ngbRanks[p]) {
189 
190         /* Add the right slice of rmine at the right place */
191         sliceSize = roffset[r+1]-roffset[r];
192         for (j = 0, count = 0; j < sliceSize; ++j) {
193           if (roffset[r]+j >= roffset[nrranks]) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
194           v = rmine[roffset[r]+j];
195           if (v >= vStart && v < vEnd) {
196             if (intOffset[p+1]+count >= numVerNgbRanksTotal) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Offset out of range");
197             interfaces_lv[intOffset[p+1]+count] = v-vStart;
198             count++;
199           }
200         }
201         intOffset[p+1] += count;
202         r++;
203       }
204       if (intOffset[p+1] != intOffset[p] + verNgbRank[p]) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unequal offsets");
205     }
206     ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
207     ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
208     for (i = 0; i < numVerNgbRanksTotal; ++i) {
209       v = gV[interfaces_lv[i]];
210       interfaces_gv[i] = v < 0 ? -v-1 : v;
211       interfaces_lv[i] += 1;
212       interfaces_gv[i] += 1;
213     }
214     ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
215     ierr = PetscFree(numVerInterfaces);CHKERRQ(ierr);
216   }
217   ierr = DMDestroy(&udm);CHKERRQ(ierr);
218 
219   /* Send the data to ParMmg and remesh */
220   ierr = DMPlexMetricNoInsertion(dm, &noInsert);CHKERRQ(ierr);
221   ierr = DMPlexMetricNoSwapping(dm, &noSwap);CHKERRQ(ierr);
222   ierr = DMPlexMetricNoMovement(dm, &noMove);CHKERRQ(ierr);
223   ierr = DMPlexMetricGetVerbosity(dm, &verbosity);CHKERRQ(ierr);
224   ierr = DMPlexMetricGetNumIterations(dm, &numIter);CHKERRQ(ierr);
225   ierr = DMPlexMetricGetGradationFactor(dm, &gradationFactor);CHKERRQ(ierr);
226   ierr = PetscCalloc2(numVertices, &verTags, numCells, &cellTags);CHKERRQ(ierr);
227   ierr = PMMG_Init_parMesh(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_pMesh, PMMG_ARG_pMet, PMMG_ARG_dim, 3, PMMG_ARG_MPIComm, tmpComm, PMMG_ARG_end);
228   ierr = PMMG_Set_meshSize(parmesh, numVertices, numCells, 0, numBdFaces, 0, 0);
229   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_APImode, PMMG_APIDISTRIB_nodes);
230   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noinsert, noInsert);
231   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_noswap, noSwap);
232   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_nomove, noMove);
233   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_verbose, verbosity);
234   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1);
235   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_niter, numIter);
236   ierr = PMMG_Set_dparameter(parmesh, PMMG_DPARAM_hgrad, gradationFactor);
237   ierr = PMMG_Set_vertices(parmesh, vertices, verTags);
238   for (i=0; i<numCells; i++) ierr = PMMG_Set_tetrahedron(parmesh, cells[4*i+0], cells[4*i+1], cells[4*i+2], cells[4*i+3], 0, i+1);
239   ierr = PMMG_Set_triangles(parmesh, bdFaces, bdFaceIds);
240   ierr = PMMG_Set_metSize(parmesh, MMG5_Vertex, numVertices, MMG5_Tensor);
241   for (i=0; i<numVertices; i++) {
242     PMMG_Set_tensorMet(parmesh, metric[6*i], metric[6*i+1], metric[6*i+2], metric[6*i+3], metric[6*i+4], metric[6*i+5], i+1);
243   }
244   ierr = PMMG_Set_numberOfNodeCommunicators(parmesh, numNgbRanks);
245   for (c = 0; c < numNgbRanks; ++c) {
246     ierr = PMMG_Set_ithNodeCommunicatorSize(parmesh, c, ngbRanks[c], intOffset[c+1]-intOffset[c]);
247     ierr = PMMG_Set_ithNodeCommunicator_nodes(parmesh, c, &interfaces_lv[intOffset[c]], &interfaces_gv[intOffset[c]], 1);
248   }
249   ierr = PMMG_parmmglib_distributed(parmesh);
250   ierr = PetscFree(cells);CHKERRQ(ierr);
251   ierr = PetscFree2(metric, vertices);CHKERRQ(ierr);
252   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
253   ierr = PetscFree2(verTags, cellTags);CHKERRQ(ierr);
254   if (numProcs > 1) {
255     ierr = PetscFree2(ngbRanks, verNgbRank);CHKERRQ(ierr);
256     ierr = PetscFree3(interfaces_lv, interfaces_gv, intOffset);CHKERRQ(ierr);
257   }
258 
259   /* Retrieve mesh from Mmg and create new Plex */
260   numCornersNew = 4;
261   ierr = PMMG_Get_meshSize(parmesh, &numVerticesNew, &numCellsNew, 0, &numFacesNew, 0, 0);
262   ierr = PetscMalloc4(dim*numVerticesNew, &verticesNew, numVerticesNew, &verTagsNew, numVerticesNew, &corners, numVerticesNew, &requiredVer);CHKERRQ(ierr);
263   ierr = PetscMalloc3((dim+1)*numCellsNew, &cellsNew, numCellsNew, &cellTagsNew, numCellsNew, &requiredCells);CHKERRQ(ierr);
264   ierr = PetscMalloc4(dim*numFacesNew, &facesNew, numFacesNew, &faceTagsNew, numFacesNew, &ridges, numFacesNew, &requiredFaces);CHKERRQ(ierr);
265   ierr = PMMG_Get_vertices(parmesh, verticesNew, verTagsNew, corners, requiredVer);
266   ierr = PMMG_Get_tetrahedra(parmesh, cellsNew, cellTagsNew, requiredCells);
267   ierr = PMMG_Get_triangles(parmesh, facesNew, faceTagsNew, requiredFaces);
268   ierr = PetscMalloc2(numVerticesNew, &owners, numVerticesNew, &gv_new);
269   ierr = PMMG_Set_iparameter(parmesh, PMMG_IPARAM_globalNum, 1);
270   ierr = PMMG_Get_verticesGloNum(parmesh, gv_new, owners);
271   for (i = 0; i < dim*numFacesNew; ++i) facesNew[i] -= 1;
272   for (i = 0; i < (dim+1)*numCellsNew; ++i) {
273     cellsNew[i] = gv_new[cellsNew[i]-1]-1;
274   }
275   numVerticesNewLoc = 0;
276   for (i = 0; i < numVerticesNew; ++i) {
277     if (owners[i] == rank) numVerticesNewLoc++;
278   }
279   ierr = PetscMalloc2(numVerticesNewLoc*dim, &verticesNewLoc, numVerticesNew, &verticesNewSorted);CHKERRQ(ierr);
280   for (i = 0, c = 0; i < numVerticesNew; i++) {
281     if (owners[i] == rank) {
282       for (j=0; j<dim; ++j) verticesNewLoc[dim*c+j] = verticesNew[dim*i+j];
283         c++;
284     }
285   }
286   ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNewLoc, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, verticesNewLoc, NULL, &verticesNewSorted, dmNew);CHKERRQ(ierr);
287   ierr = PMMG_Free_all(PMMG_ARG_start, PMMG_ARG_ppParMesh, &parmesh, PMMG_ARG_end);
288 
289   /* Rebuild boundary label*/
290   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
291   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
292   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
293   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
294   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
295   for (i = 0; i < numFacesNew; i++) {
296     PetscInt        numCoveredPoints, numFaces = 0, facePoints[3];
297     const PetscInt *coveredPoints = NULL;
298 
299     for (j = 0; j < dim; ++j) {
300       lv = facesNew[i*dim+j];
301       gv = gv_new[lv]-1;
302       ierr = PetscFindInt(gv, numVerticesNew, verticesNewSorted, &lv);CHKERRQ(ierr);
303       facePoints[j] = lv+vStart;
304     }
305     ierr = DMPlexGetFullJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
306     for (j = 0; j < numCoveredPoints; ++j) {
307       if (coveredPoints[j] >= fStart && coveredPoints[j] < fEnd) {
308         numFaces++;
309         f = j;
310       }
311     }
312     if (numFaces != 1) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "%d vertices cannot define more than 1 facet (%d)", dim, numFaces);
313     ierr = DMLabelSetValue(bdLabelNew, coveredPoints[f], faceTagsNew[i]);CHKERRQ(ierr);
314     ierr = DMPlexRestoreJoin(*dmNew, dim, facePoints, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
315   }
316 
317   /* Clean up */
318   ierr = PetscFree4(verticesNew, verTagsNew, corners, requiredVer);CHKERRQ(ierr);
319   ierr = PetscFree3(cellsNew, cellTagsNew, requiredCells);CHKERRQ(ierr);
320   ierr = PetscFree4(facesNew, faceTagsNew, ridges, requiredFaces);CHKERRQ(ierr);
321   ierr = PetscFree2(owners, gv_new);CHKERRQ(ierr);
322   ierr = PetscFree2(verticesNewLoc, verticesNewSorted);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325