xref: /petsc/src/dm/impls/plex/plexadapt.c (revision ace5e5346c708ab88b64038dff8fd16ec16aa51f)
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"
104190e864Sbarral PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdyLabelName[], DM *dmCoarsened)
110bbe5d31Sbarral {
120bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
130bbe5d31Sbarral   DM             udm, coordDM;
144190e864Sbarral   DMLabel        bd, bdtags;
150bbe5d31Sbarral   Vec            coordinates;
160bbe5d31Sbarral   PetscSection   coordSection;
170bbe5d31Sbarral   const PetscScalar *coords;
180bbe5d31Sbarral   double        *coarseCoords;
190bbe5d31Sbarral   IS             bdIS;
200bbe5d31Sbarral   PetscReal     *x, *y, *z;
210bbe5d31Sbarral   const PetscInt *faces;
220bbe5d31Sbarral   PetscInt      *cells, *ccells, *bdFaces, *bdFaceIds;
230bbe5d31Sbarral   PetscInt       dim, numCorners, cStart, cEnd, numCells, numCoarseCells, c, vStart, vEnd, numVertices, numCoarseVertices, v, numBdFaces, f, maxConeSize, bdSize, coff;
240bbe5d31Sbarral   const PetscScalar   *metricArray;
250bbe5d31Sbarral   PetscReal     *met;
26*ace5e534Sbarral   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
27f7caa48eSbarral   PetscInt      *boundaryTags;
28f7caa48eSbarral 
290bbe5d31Sbarral #endif
300bbe5d31Sbarral   PetscErrorCode ierr;
310bbe5d31Sbarral   MPI_Comm       comm;
320bbe5d31Sbarral 
330bbe5d31Sbarral   PetscFunctionBegin;
340bbe5d31Sbarral   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
350bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
360bbe5d31Sbarral   if (!dm->coarseMesh) {
370bbe5d31Sbarral     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
380bbe5d31Sbarral     ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
390bbe5d31Sbarral     ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr);
400bbe5d31Sbarral     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
410bbe5d31Sbarral     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
420bbe5d31Sbarral     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
430bbe5d31Sbarral     ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
440bbe5d31Sbarral     ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
450bbe5d31Sbarral     numCells    = cEnd - cStart;
460bbe5d31Sbarral     numVertices = vEnd - vStart;
470bbe5d31Sbarral     ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numCells*maxConeSize, &cells, dim*dim*numVertices, &met);CHKERRQ(ierr);
480bbe5d31Sbarral     ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
490bbe5d31Sbarral     ierr = VecGetArrayRead(metric, &metricArray);CHKERRQ(ierr);
500bbe5d31Sbarral     for (v = vStart; v < vEnd; ++v) {
510bbe5d31Sbarral       PetscInt off;
520bbe5d31Sbarral 
530bbe5d31Sbarral       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
540bbe5d31Sbarral       x[v-vStart] = coords[off+0];
550bbe5d31Sbarral       y[v-vStart] = coords[off+1];
560bbe5d31Sbarral       if (dim > 2) z[v-vStart] = coords[off+2];
570bbe5d31Sbarral 
580bbe5d31Sbarral       ierr = PetscMemcpy(&met[dim*dim*(v-vStart)], &metricArray[dim*off], dim*dim*sizeof(PetscScalar));CHKERRQ(ierr);
590bbe5d31Sbarral     }
600bbe5d31Sbarral     ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
610bbe5d31Sbarral     ierr = VecRestoreArrayRead(metric, &metricArray);CHKERRQ(ierr);
620bbe5d31Sbarral     for (c = 0, coff = 0; c < numCells; ++c) {
630bbe5d31Sbarral       const PetscInt *cone;
640bbe5d31Sbarral       PetscInt        coneSize, cl;
650bbe5d31Sbarral 
660bbe5d31Sbarral       ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
670bbe5d31Sbarral       ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
680bbe5d31Sbarral       for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
690bbe5d31Sbarral     }
700bbe5d31Sbarral     switch (dim) {
710bbe5d31Sbarral     case 2:
720bbe5d31Sbarral       pragmatic_2d_init(&numVertices, &numCells, cells, x, y);
730bbe5d31Sbarral       break;
740bbe5d31Sbarral     case 3:
750bbe5d31Sbarral       pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z);
760bbe5d31Sbarral       break;
770bbe5d31Sbarral     default:
780bbe5d31Sbarral       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
790bbe5d31Sbarral     }
800bbe5d31Sbarral     /* Create boundary mesh */
814190e864Sbarral     if ( bdyLabelName ) {
824190e864Sbarral       ierr = DMGetLabel(dm, bdyLabelName, &bdtags);CHKERRQ(ierr);
834190e864Sbarral     }
844190e864Sbarral     // TODO fix if bdyLabelName = "boundary"
8552b40ca0Sbarral     // TODO a way to use only bdyLabelName would be to loop over its strata, if I can't get all the strata at once ?
860bbe5d31Sbarral     ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr);
870bbe5d31Sbarral     ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr);
880bbe5d31Sbarral     ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr);
890bbe5d31Sbarral     ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr);
900bbe5d31Sbarral     ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr);
9152b40ca0Sbarral     // TODO why not assume that we are considering simplicial meshes, in which case bdSize = dim*numBdFaces ?
920bbe5d31Sbarral     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
930bbe5d31Sbarral       PetscInt *closure = NULL;
940bbe5d31Sbarral       PetscInt  closureSize, cl;
950bbe5d31Sbarral 
960bbe5d31Sbarral       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
970bbe5d31Sbarral       for (cl = 0; cl < closureSize*2; cl += 2) {
980bbe5d31Sbarral         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
990bbe5d31Sbarral       }
1000bbe5d31Sbarral       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1010bbe5d31Sbarral     }
1020bbe5d31Sbarral     ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
1030bbe5d31Sbarral     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1040bbe5d31Sbarral       PetscInt *closure = NULL;
1050bbe5d31Sbarral       PetscInt  closureSize, cl;
1060bbe5d31Sbarral 
1070bbe5d31Sbarral       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1080bbe5d31Sbarral       for (cl = 0; cl < closureSize*2; cl += 2) {
1090bbe5d31Sbarral         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
1100bbe5d31Sbarral       }
1114190e864Sbarral       if ( bdyLabelName )
1124190e864Sbarral         DMLabelGetValue(bdtags, faces[f], &bdFaceIds[f]);
1134190e864Sbarral       else
1140bbe5d31Sbarral         bdFaceIds[f] = 1;
1150bbe5d31Sbarral       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1160bbe5d31Sbarral     }
1170bbe5d31Sbarral     ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
1180bbe5d31Sbarral     ierr = DMLabelDestroy(&bd);CHKERRQ(ierr);
1190bbe5d31Sbarral     pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
1200bbe5d31Sbarral     pragmatic_set_metric(met);
1210bbe5d31Sbarral     pragmatic_adapt();
1220bbe5d31Sbarral     /* Read out mesh */
1230bbe5d31Sbarral     pragmatic_get_info(&numCoarseVertices, &numCoarseCells);
1240bbe5d31Sbarral     ierr = PetscMalloc1(numCoarseVertices*dim, &coarseCoords);CHKERRQ(ierr);
1250bbe5d31Sbarral     switch (dim) {
1260bbe5d31Sbarral     case 2:
1270bbe5d31Sbarral       ierr = PetscMalloc2(numCoarseVertices, &x, numCoarseVertices, &y);CHKERRQ(ierr);
1280bbe5d31Sbarral       pragmatic_get_coords_2d(x, y);
1290bbe5d31Sbarral       numCorners = 3;
1300bbe5d31Sbarral       for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*2+0] = x[v]; coarseCoords[v*2+1] = y[v];}
1310bbe5d31Sbarral       break;
1320bbe5d31Sbarral     case 3:
1330bbe5d31Sbarral       ierr = PetscMalloc3(numCoarseVertices, &x, numCoarseVertices, &y, numCoarseVertices, &z);CHKERRQ(ierr);
1340bbe5d31Sbarral       pragmatic_get_coords_3d(x, y, z);
1350bbe5d31Sbarral       numCorners = 4;
1360bbe5d31Sbarral       for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*3+0] = x[v]; coarseCoords[v*3+1] = y[v]; coarseCoords[v*3+2] = z[v];}
1370bbe5d31Sbarral       break;
1380bbe5d31Sbarral     default:
1390bbe5d31Sbarral       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim);
1400bbe5d31Sbarral     }
1410bbe5d31Sbarral     ierr = PetscMalloc1(numCoarseCells*(dim+1), &ccells);CHKERRQ(ierr); // only for simplicial meshes
1420bbe5d31Sbarral     pragmatic_get_elements(ccells);
1430bbe5d31Sbarral     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCoarseCells, numCoarseVertices, numCorners, PETSC_TRUE, ccells, dim, coarseCoords, &dm->coarseMesh);CHKERRQ(ierr);
144f7caa48eSbarral     /* Read out boundary tags */
145f7caa48eSbarral     pragmatic_get_boundaryTags(&boundaryTags);
146f7caa48eSbarral     if ( !bdyLabelName ) bdyLabelName = "boundary";
14752b40ca0Sbarral     ierr = DMCreateLabel(dm->coarseMesh, bdyLabelName);CHKERRQ(ierr);
14852b40ca0Sbarral     ierr = DMGetLabel(dm->coarseMesh, bdyLabelName, &bd);CHKERRQ(ierr);
149f7caa48eSbarral     ierr = DMPlexGetHeightStratum(dm->coarseMesh, 0, &cStart, &cEnd);CHKERRQ(ierr);
150f7caa48eSbarral     ierr = DMPlexGetDepthStratum(dm->coarseMesh, 0, &vStart, &vEnd);CHKERRQ(ierr);
151f7caa48eSbarral     for (c = cStart; c < cEnd; ++c) {
152f7caa48eSbarral       PetscInt       *cellClosure = NULL;
153f7caa48eSbarral       PetscInt        cellClosureSize, cl;
154f7caa48eSbarral       PetscInt       *facetClosure = NULL;
155f7caa48eSbarral       PetscInt        facetClosureSize, cl2;
156f7caa48eSbarral       const PetscInt *facetList;
157f7caa48eSbarral       PetscInt        facetListSize, f;
158f7caa48eSbarral 
159f7caa48eSbarral       cell = &ccells[(c-cStart)*(dim+1)]; // pointing to the current cell in the ccell table
160f7caa48eSbarral       ierr = DMPlexGetTransitiveClosure(dm->coarseMesh, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr);
161f7caa48eSbarral       // first, encode the permutation of the vertices indices betwenn the cell closure and pragmatic ordering
162f7caa48eSbarral       for (cl = 0; cl < cellClosureSize*2; cl+=2 ) {
163f7caa48eSbarral         if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) continue;
164f7caa48eSbarral         iVer = cellClosure[cl] - vStart;
165f7caa48eSbarral         for (i=0; i<dim+1; ++i) {
166f7caa48eSbarral           if ( iVer == cell[i] ) {
16752b40ca0Sbarral             perm[cl] = i;    // the cl-th element of the closure is the i-th vertex of the cell
168f7caa48eSbarral             break;
169f7caa48eSbarral           }
170f7caa48eSbarral         }
171f7caa48eSbarral       }
172f7caa48eSbarral       ierr = DMPlexGetCone(dm->coarseMesh, c, &facetList);CHKERRQ(ierr);     // list of edges/facets of the cell
173f7caa48eSbarral       ierr = DMPlexGetConeSize(dm->coarseMesh, c, &facetListSize);CHKERRQ(ierr);
174f7caa48eSbarral       // then, for each edge/facet of the cell, find the opposite vertex (ie the one not in the closure of the facet/edge)
175f7caa48eSbarral       for (f=0; f<facetListSize; ++f){
176f7caa48eSbarral         facet = facetList[f];
177f7caa48eSbarral         ierr = DMPlexGetTransitiveClosure(dm->coarseMesh, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr);
178f7caa48eSbarral         // loop over the vertices of the cell closure, and check if each vertex belongs to the edge/facet closure
179f7caa48eSbarral         PetscMemzero(isInFacetClosure, sizeof(isInFacetClosure));
180f7caa48eSbarral         for (cl = 0; cl < cellClosureSize*2; cl+=2) {
181f7caa48eSbarral           if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) isInFacetClosure[cl] = 1;
182f7caa48eSbarral           for (cl2 =0; cl2<facetClosureSize*2; cl2+=2){
183f7caa48eSbarral             if ((facetClosure[cl2] < vStart) || (facetClosure[cl2] >= vEnd)) continue;
184f7caa48eSbarral             if ( cellClosure[cl] == facetClosure[cl2] ) {
185f7caa48eSbarral               isInFacetClosure[cl] = 1;
186f7caa48eSbarral             }
187f7caa48eSbarral           }
188f7caa48eSbarral         }
189f7caa48eSbarral         // 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
190f7caa48eSbarral         for (cl = 0; cl < cellClosureSize*2; cl+=2 ) {
191f7caa48eSbarral           if ( !isInFacetClosure[cl] ) {
192f7caa48eSbarral             idx = (c-cStart)*(dim+1) + perm[cl];
193f7caa48eSbarral             if ( boundaryTags[idx] )
19452b40ca0Sbarral //              ierr = DMSetLabelValue(dm->coarseMesh, bdyLabelName, facet, boundaryTags[idx]);CHKERRQ(ierr);
19552b40ca0Sbarral               ierr = DMLabelSetValue(bd, facet, boundaryTags[idx]);CHKERRQ(ierr);
196f7caa48eSbarral             break;
197f7caa48eSbarral           }
198f7caa48eSbarral         }
199f7caa48eSbarral         ierr = DMPlexRestoreTransitiveClosure(dm->coarseMesh, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr);
200f7caa48eSbarral       }
201f7caa48eSbarral       ierr = DMPlexRestoreTransitiveClosure(dm->coarseMesh, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr);
202f7caa48eSbarral     }
2030bbe5d31Sbarral     pragmatic_finalize();
204*ace5e534Sbarral     ierr = PetscFree5(x, y, z, cells, met);CHKERRQ(ierr);
2050bbe5d31Sbarral     ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
2060bbe5d31Sbarral     ierr = PetscFree(coarseCoords);CHKERRQ(ierr);
207*ace5e534Sbarral     ierr = PetscFree(ccells);CHKERRQ(ierr);
2080bbe5d31Sbarral   }
2090bbe5d31Sbarral #endif
2100bbe5d31Sbarral   ierr = PetscObjectReference((PetscObject) dm->coarseMesh);CHKERRQ(ierr);
2110bbe5d31Sbarral   *dmCoarsened = dm->coarseMesh;
2120bbe5d31Sbarral   PetscFunctionReturn(0);
2130bbe5d31Sbarral }
214