xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 6b1afcafdcdda35ab85c942a09161cc59d2f6e07)
10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
20bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
30bbe5d31Sbarral #include <pragmatic/cpragmatic.h>
40bbe5d31Sbarral #endif
50bbe5d31Sbarral 
60bbe5d31Sbarral 
70bbe5d31Sbarral 
80bbe5d31Sbarral #undef __FUNCT__
952b40ca0Sbarral #define __FUNCT__ "DMPlexAdapt"
10*6b1afcafSbarral /* write docstring */
11*6b1afcafSbarral PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdyLabelName[], DM *dmAdapted)
120bbe5d31Sbarral {
130bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
140bbe5d31Sbarral   DM                   udm, coordDM;
15*6b1afcafSbarral   DMLabel              bd, bdTags, bdTagsAdp;
160bbe5d31Sbarral   Vec                  coordinates;
170bbe5d31Sbarral   PetscSection         coordSection;
180bbe5d31Sbarral   const PetscScalar   *coords;
19*6b1afcafSbarral   double              *coordsAdp;
200bbe5d31Sbarral   IS                   bdIS;
21*6b1afcafSbarral   PetscReal           *x, *y, *z, *xAdp, *yAdp, *zAdp;
220bbe5d31Sbarral   const PetscInt      *faces;
23*6b1afcafSbarral   PetscInt            *cells, *cellsAdp, *bdFaces, *bdFaceIds;
24*6b1afcafSbarral   PetscInt             dim, numCornersAdp, cStart, cEnd, numCells, numCellsAdp, c;
25*6b1afcafSbarral   PetscInt             vStart, vEnd, numVertices, numVerticesAdp, v;
26*6b1afcafSbarral   PetscInt             numBdFaces, f, maxConeSize, bdSize, coff;
270bbe5d31Sbarral   const PetscScalar   *metricArray;
280bbe5d31Sbarral   PetscReal           *met;
29ace5e534Sbarral   PetscInt            *cell, perm[30], isInFacetClosure[30], iVer, i, idx, facet; // 30 is twice the max size of a cell closure in 3D for tet meshes
30f7caa48eSbarral   PetscInt            *boundaryTags;
31*6b1afcafSbarral   PetscBool           flag;
32f7caa48eSbarral 
330bbe5d31Sbarral #endif
340bbe5d31Sbarral   PetscErrorCode ierr;
350bbe5d31Sbarral   MPI_Comm       comm;
360bbe5d31Sbarral 
370bbe5d31Sbarral   PetscFunctionBegin;
38*6b1afcafSbarral   printf("DEBUG  PETSc  ======= DMPlexAdapt begin\n");
390bbe5d31Sbarral   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
400bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
410bbe5d31Sbarral   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
420bbe5d31Sbarral   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
430bbe5d31Sbarral   ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr);
440bbe5d31Sbarral   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
450bbe5d31Sbarral   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
460bbe5d31Sbarral   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
470bbe5d31Sbarral   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
480bbe5d31Sbarral   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
490bbe5d31Sbarral   numCells    = cEnd - cStart;
500bbe5d31Sbarral   numVertices = vEnd - vStart;
510bbe5d31Sbarral   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numCells*maxConeSize, &cells, dim*dim*numVertices, &met);CHKERRQ(ierr);
520bbe5d31Sbarral   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
530bbe5d31Sbarral   ierr = VecGetArrayRead(metric, &metricArray);CHKERRQ(ierr);
540bbe5d31Sbarral   for (v = vStart; v < vEnd; ++v) {
550bbe5d31Sbarral     PetscInt off;
560bbe5d31Sbarral 
570bbe5d31Sbarral     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
580bbe5d31Sbarral     x[v-vStart] = coords[off+0];
590bbe5d31Sbarral     y[v-vStart] = coords[off+1];
600bbe5d31Sbarral     if (dim > 2) z[v-vStart] = coords[off+2];
610bbe5d31Sbarral 
620bbe5d31Sbarral     ierr = PetscMemcpy(&met[dim*dim*(v-vStart)], &metricArray[dim*off], dim*dim*sizeof(PetscScalar));CHKERRQ(ierr);
630bbe5d31Sbarral   }
640bbe5d31Sbarral   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
650bbe5d31Sbarral   ierr = VecRestoreArrayRead(metric, &metricArray);CHKERRQ(ierr);
660bbe5d31Sbarral   for (c = 0, coff = 0; c < numCells; ++c) {
670bbe5d31Sbarral     const PetscInt *cone;
680bbe5d31Sbarral     PetscInt        coneSize, cl;
690bbe5d31Sbarral 
700bbe5d31Sbarral     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
710bbe5d31Sbarral     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
720bbe5d31Sbarral     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
730bbe5d31Sbarral   }
740bbe5d31Sbarral   switch (dim) {
750bbe5d31Sbarral   case 2:
760bbe5d31Sbarral     pragmatic_2d_init(&numVertices, &numCells, cells, x, y);
770bbe5d31Sbarral     break;
780bbe5d31Sbarral   case 3:
790bbe5d31Sbarral     pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z);
800bbe5d31Sbarral     break;
810bbe5d31Sbarral   default:
820bbe5d31Sbarral     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
830bbe5d31Sbarral   }
840bbe5d31Sbarral   /* Create boundary mesh */
85*6b1afcafSbarral   printf("DEBUG  PETSc  creating boundary mesh\n");
864190e864Sbarral   if (bdyLabelName) {
87*6b1afcafSbarral     ierr = PetscStrcmp(bdyLabelName, "boundary", &flag);CHKERRQ(ierr);
88*6b1afcafSbarral     if (flag) {
89*6b1afcafSbarral       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdyLabelName);
904190e864Sbarral     }
91*6b1afcafSbarral     ierr = DMGetLabel(dm, bdyLabelName, &bdTags);CHKERRQ(ierr);
92*6b1afcafSbarral   }
93*6b1afcafSbarral   /* TODO a way to use only bdyLabelName would be to loop over its strata, if I can't get all the strata at once ? */
940bbe5d31Sbarral   ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr);
950bbe5d31Sbarral   ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr);
960bbe5d31Sbarral   ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr);
970bbe5d31Sbarral   ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr);
980bbe5d31Sbarral   ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr);
99*6b1afcafSbarral   /* TODO why not assume that we are considering simplicial meshes, in which case bdSize = dim*numBdFaces ? */
1000bbe5d31Sbarral   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1010bbe5d31Sbarral     PetscInt *closure = NULL;
1020bbe5d31Sbarral     PetscInt  closureSize, cl;
1030bbe5d31Sbarral 
1040bbe5d31Sbarral     ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1050bbe5d31Sbarral     for (cl = 0; cl < closureSize*2; cl += 2) {
1060bbe5d31Sbarral       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1070bbe5d31Sbarral     }
1080bbe5d31Sbarral     ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1090bbe5d31Sbarral   }
1100bbe5d31Sbarral   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
1110bbe5d31Sbarral   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1120bbe5d31Sbarral     PetscInt *closure = NULL;
1130bbe5d31Sbarral     PetscInt  closureSize, cl;
1140bbe5d31Sbarral 
1150bbe5d31Sbarral     ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1160bbe5d31Sbarral     for (cl = 0; cl < closureSize*2; cl += 2) {
1170bbe5d31Sbarral       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
1180bbe5d31Sbarral     }
119*6b1afcafSbarral     if (bdyLabelName) {
120*6b1afcafSbarral       ierr = DMLabelGetValue(bdTags, faces[f], &bdFaceIds[f]);CHKERRQ(ierr);
121*6b1afcafSbarral     } else {
1220bbe5d31Sbarral       bdFaceIds[f] = 1;
123*6b1afcafSbarral     }
1240bbe5d31Sbarral     ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1250bbe5d31Sbarral   }
1260bbe5d31Sbarral   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
1270bbe5d31Sbarral   ierr = DMLabelDestroy(&bd);CHKERRQ(ierr);
1280bbe5d31Sbarral   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
1290bbe5d31Sbarral   pragmatic_set_metric(met);
130*6b1afcafSbarral   printf("DEBUG  PETSc  start pragmatic adapt\n");
1310bbe5d31Sbarral   pragmatic_adapt();
132*6b1afcafSbarral   printf("DEBUG  PETSc  end pragmatic adapt\n");
133*6b1afcafSbarral   ierr = PetscFree5(x, y, z, cells, met);CHKERRQ(ierr);
134*6b1afcafSbarral   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
135*6b1afcafSbarral   ierr = PetscFree(coords);CHKERRQ(ierr);
1360bbe5d31Sbarral   /* Read out mesh */
137*6b1afcafSbarral   pragmatic_get_info(&numVerticesAdp, &numCellsAdp);
138*6b1afcafSbarral   ierr = PetscMalloc1(numVerticesAdp*dim, &coordsAdp);CHKERRQ(ierr);
1390bbe5d31Sbarral   switch (dim) {
1400bbe5d31Sbarral   case 2:
141*6b1afcafSbarral     ierr = PetscMalloc2(numVerticesAdp, &xAdp, numVerticesAdp, &yAdp);CHKERRQ(ierr);
142*6b1afcafSbarral     zAdp = NULL;
143*6b1afcafSbarral     pragmatic_get_coords_2d(xAdp, yAdp);
144*6b1afcafSbarral     numCornersAdp = 3;
145*6b1afcafSbarral     for (v = 0; v < numVerticesAdp; ++v) {coordsAdp[v*2+0] = xAdp[v]; coordsAdp[v*2+1] = yAdp[v];}
1460bbe5d31Sbarral     break;
1470bbe5d31Sbarral   case 3:
148*6b1afcafSbarral     ierr = PetscMalloc3(numVerticesAdp, &xAdp, numVerticesAdp, &yAdp, numVerticesAdp, &zAdp);CHKERRQ(ierr);
149*6b1afcafSbarral     pragmatic_get_coords_3d(xAdp, yAdp, zAdp);
150*6b1afcafSbarral     numCornersAdp = 4;
151*6b1afcafSbarral     for (v = 0; v < numVerticesAdp; ++v) {coordsAdp[v*3+0] = xAdp[v]; coordsAdp[v*3+1] = yAdp[v]; coordsAdp[v*3+2] = zAdp[v];}
1520bbe5d31Sbarral     break;
1530bbe5d31Sbarral   default:
154*6b1afcafSbarral     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
1550bbe5d31Sbarral   }
156*6b1afcafSbarral   ierr = PetscMalloc1(numCellsAdp*(dim+1), &cellsAdp);CHKERRQ(ierr); // only for simplicial meshes
157*6b1afcafSbarral   pragmatic_get_elements(cellsAdp);
158*6b1afcafSbarral   ierr = DMPlexCreate(comm, dmAdapted);CHKERRQ(ierr);
159*6b1afcafSbarral   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCellsAdp, numVerticesAdp, numCornersAdp, PETSC_TRUE, cellsAdp, dim, coordsAdp, dmAdapted);CHKERRQ(ierr);
160f7caa48eSbarral   /* Read out boundary tags */
161*6b1afcafSbarral   printf("DEBUG  PETSc  boundary tags reconstruction\n");
162f7caa48eSbarral   pragmatic_get_boundaryTags(&boundaryTags);
163f7caa48eSbarral   if (!bdyLabelName) bdyLabelName = "boundary";
164*6b1afcafSbarral   ierr = DMCreateLabel(*dmAdapted, bdyLabelName);CHKERRQ(ierr);
165*6b1afcafSbarral   ierr = DMGetLabel(*dmAdapted, bdyLabelName, &bdTagsAdp);CHKERRQ(ierr);
166*6b1afcafSbarral   ierr = DMPlexGetHeightStratum(*dmAdapted, 0, &cStart, &cEnd);CHKERRQ(ierr);
167*6b1afcafSbarral   ierr = DMPlexGetDepthStratum(*dmAdapted, 0, &vStart, &vEnd);CHKERRQ(ierr);
168f7caa48eSbarral   for (c = cStart; c < cEnd; ++c) {
169f7caa48eSbarral     PetscInt       *cellClosure = NULL;
170f7caa48eSbarral     PetscInt        cellClosureSize, cl;
171f7caa48eSbarral     PetscInt       *facetClosure = NULL;
172f7caa48eSbarral     PetscInt        facetClosureSize, cl2;
173f7caa48eSbarral     const PetscInt *facetList;
174f7caa48eSbarral     PetscInt        facetListSize, f;
175f7caa48eSbarral 
176*6b1afcafSbarral     // TODO: ne puis-je pas sauter les cells sans aucune frontière ?
177*6b1afcafSbarral     cell = &cellsAdp[(c-cStart)*(dim+1)]; // pointing to the current cell in the ccell table
178*6b1afcafSbarral     ierr = DMPlexGetTransitiveClosure(*dmAdapted, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr);
179*6b1afcafSbarral     /* first, encode the permutation of the vertices indices betwenn the cell closure and pragmatic ordering */
180f7caa48eSbarral     for (cl = 0; cl < cellClosureSize*2; cl+=2) {
181f7caa48eSbarral       if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) continue;
182f7caa48eSbarral       iVer = cellClosure[cl] - vStart;
183f7caa48eSbarral       for (i=0; i<dim+1; ++i) {
184f7caa48eSbarral         if (iVer == cell[i]) {
18552b40ca0Sbarral           perm[cl] = i;    // the cl-th element of the closure is the i-th vertex of the cell
186f7caa48eSbarral           break;
187f7caa48eSbarral         }
188f7caa48eSbarral       }
189f7caa48eSbarral     }
190*6b1afcafSbarral     ierr = DMPlexGetCone(*dmAdapted, c, &facetList);CHKERRQ(ierr);     // list of edges/facets of the cell
191*6b1afcafSbarral     ierr = DMPlexGetConeSize(*dmAdapted, c, &facetListSize);CHKERRQ(ierr);
192*6b1afcafSbarral     /* then, for each edge/facet of the cell, find the opposite vertex (ie the one not in the closure of the facet/edge) */
193f7caa48eSbarral     for (f=0; f<facetListSize; ++f){
194f7caa48eSbarral       facet = facetList[f];
195*6b1afcafSbarral       ierr = DMPlexGetTransitiveClosure(*dmAdapted, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr);
196*6b1afcafSbarral       /* loop over the vertices of the cell closure, and check if each vertex belongs to the edge/facet closure */
197f7caa48eSbarral       PetscMemzero(isInFacetClosure, sizeof(isInFacetClosure));
198f7caa48eSbarral       for (cl = 0; cl < cellClosureSize*2; cl+=2) {
199f7caa48eSbarral         if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) isInFacetClosure[cl] = 1;
200f7caa48eSbarral         for (cl2 =0; cl2<facetClosureSize*2; cl2+=2){
201f7caa48eSbarral           if ((facetClosure[cl2] < vStart) || (facetClosure[cl2] >= vEnd)) continue;
202f7caa48eSbarral           if (cellClosure[cl] == facetClosure[cl2]) {
203f7caa48eSbarral             isInFacetClosure[cl] = 1;
204f7caa48eSbarral           }
205f7caa48eSbarral         }
206f7caa48eSbarral       }
207*6b1afcafSbarral       /* the vertex that was not marked is the vertex opposite to the edge/facet, ie the one giving the edge/facet boundary tag in pragmatic */
208f7caa48eSbarral       for (cl = 0; cl < cellClosureSize*2; cl+=2) {
209f7caa48eSbarral         if (!isInFacetClosure[cl]) {
210f7caa48eSbarral           idx = (c-cStart)*(dim+1) + perm[cl];
211*6b1afcafSbarral           if (boundaryTags[idx]) {
212*6b1afcafSbarral             ierr = DMLabelSetValue(bdTagsAdp, facet, boundaryTags[idx]);CHKERRQ(ierr);
213*6b1afcafSbarral           }
214f7caa48eSbarral           break;
215f7caa48eSbarral         }
216f7caa48eSbarral       }
217*6b1afcafSbarral       ierr = DMPlexRestoreTransitiveClosure(*dmAdapted, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr);
218f7caa48eSbarral     }
219*6b1afcafSbarral     ierr = DMPlexRestoreTransitiveClosure(*dmAdapted, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr);
220f7caa48eSbarral   }
2210bbe5d31Sbarral   pragmatic_finalize();
222*6b1afcafSbarral   ierr = PetscFree5(xAdp, yAdp, zAdp, cellsAdp, coordsAdp);CHKERRQ(ierr);
2230bbe5d31Sbarral #endif
224*6b1afcafSbarral   printf("DEBUG  PETSc  ======= DMPlexAdapt end\n");
2250bbe5d31Sbarral   PetscFunctionReturn(0);
2260bbe5d31Sbarral }
227