xref: /petsc/src/dm/impls/plex/plexadapt.c (revision f7caa48e5051722277da12d38e1221bf932ec792)
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__
90bbe5d31Sbarral #define __FUNCT__ "DMPlexAdaptdm->coarseMesh"
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*f7caa48eSbarral   PetscInt       *cell, perm[28], isInFacetClosure[28], iVer, i, idx, facet;
27*f7caa48eSbarral   PetscInt      *boundaryTags;
28*f7caa48eSbarral 
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"
854190e864Sbarral     // 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);
910bbe5d31Sbarral     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
920bbe5d31Sbarral       PetscInt *closure = NULL;
930bbe5d31Sbarral       PetscInt  closureSize, cl;
940bbe5d31Sbarral 
950bbe5d31Sbarral       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
960bbe5d31Sbarral       for (cl = 0; cl < closureSize*2; cl += 2) {
970bbe5d31Sbarral         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
980bbe5d31Sbarral       }
990bbe5d31Sbarral       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1000bbe5d31Sbarral     }
1010bbe5d31Sbarral     ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
1020bbe5d31Sbarral     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1030bbe5d31Sbarral       PetscInt *closure = NULL;
1040bbe5d31Sbarral       PetscInt  closureSize, cl;
1050bbe5d31Sbarral 
1060bbe5d31Sbarral       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1070bbe5d31Sbarral       for (cl = 0; cl < closureSize*2; cl += 2) {
1080bbe5d31Sbarral         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
1090bbe5d31Sbarral       }
1104190e864Sbarral       if ( bdyLabelName )
1114190e864Sbarral         DMLabelGetValue(bdtags, faces[f], &bdFaceIds[f]);
1124190e864Sbarral       else
1130bbe5d31Sbarral         bdFaceIds[f] = 1;
1140bbe5d31Sbarral       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1150bbe5d31Sbarral     }
1160bbe5d31Sbarral     ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
1170bbe5d31Sbarral     ierr = DMLabelDestroy(&bd);CHKERRQ(ierr);
1180bbe5d31Sbarral     pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
1190bbe5d31Sbarral     pragmatic_set_metric(met);
1200bbe5d31Sbarral     pragmatic_adapt();
1210bbe5d31Sbarral     /* Read out mesh */
1220bbe5d31Sbarral     pragmatic_get_info(&numCoarseVertices, &numCoarseCells);
1230bbe5d31Sbarral     ierr = PetscMalloc1(numCoarseVertices*dim, &coarseCoords);CHKERRQ(ierr);
1240bbe5d31Sbarral     switch (dim) {
1250bbe5d31Sbarral     case 2:
1260bbe5d31Sbarral       ierr = PetscMalloc2(numCoarseVertices, &x, numCoarseVertices, &y);CHKERRQ(ierr);
1270bbe5d31Sbarral       pragmatic_get_coords_2d(x, y);
1280bbe5d31Sbarral       numCorners = 3;
1290bbe5d31Sbarral       for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*2+0] = x[v]; coarseCoords[v*2+1] = y[v];}
1300bbe5d31Sbarral       break;
1310bbe5d31Sbarral     case 3:
1320bbe5d31Sbarral       ierr = PetscMalloc3(numCoarseVertices, &x, numCoarseVertices, &y, numCoarseVertices, &z);CHKERRQ(ierr);
1330bbe5d31Sbarral       pragmatic_get_coords_3d(x, y, z);
1340bbe5d31Sbarral       numCorners = 4;
1350bbe5d31Sbarral       for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*3+0] = x[v]; coarseCoords[v*3+1] = y[v]; coarseCoords[v*3+2] = z[v];}
1360bbe5d31Sbarral       break;
1370bbe5d31Sbarral     default:
1380bbe5d31Sbarral       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim);
1390bbe5d31Sbarral     }
1400bbe5d31Sbarral     ierr = PetscMalloc1(numCoarseCells*(dim+1), &ccells);CHKERRQ(ierr); // only for simplicial meshes
1410bbe5d31Sbarral     pragmatic_get_elements(ccells);
1420bbe5d31Sbarral     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCoarseCells, numCoarseVertices, numCorners, PETSC_TRUE, ccells, dim, coarseCoords, &dm->coarseMesh);CHKERRQ(ierr);
143*f7caa48eSbarral     /* Read out boundary tags */
144*f7caa48eSbarral     pragmatic_get_boundaryTags(&boundaryTags);
145*f7caa48eSbarral     if ( !bdyLabelName ) bdyLabelName = "boundary";
146*f7caa48eSbarral     ierr = DMCreateLabel(dm->coarseMesh, bdyLabelName);;CHKERRQ(ierr);
147*f7caa48eSbarral     ierr = DMPlexGetHeightStratum(dm->coarseMesh, 0, &cStart, &cEnd);CHKERRQ(ierr);
148*f7caa48eSbarral     ierr = DMPlexGetDepthStratum(dm->coarseMesh, 0, &vStart, &vEnd);CHKERRQ(ierr);
149*f7caa48eSbarral     for (c = cStart; c < cEnd; ++c) {
150*f7caa48eSbarral       PetscInt       *cellClosure = NULL;
151*f7caa48eSbarral       PetscInt        cellClosureSize, cl;
152*f7caa48eSbarral       PetscInt       *facetClosure = NULL;
153*f7caa48eSbarral       PetscInt        facetClosureSize, cl2;
154*f7caa48eSbarral       const PetscInt *facetList;
155*f7caa48eSbarral       PetscInt        facetListSize, f;
156*f7caa48eSbarral 
157*f7caa48eSbarral       cell = &ccells[(c-cStart)*(dim+1)]; // pointing to the current cell in the ccell table
158*f7caa48eSbarral       ierr = DMPlexGetTransitiveClosure(dm->coarseMesh, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr);
159*f7caa48eSbarral       // first, encode the permutation of the vertices indices betwenn the cell closure and pragmatic ordering
160*f7caa48eSbarral       for (cl = 0; cl < cellClosureSize*2; cl+=2 ) {
161*f7caa48eSbarral         if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) continue;
162*f7caa48eSbarral         iVer = cellClosure[cl] - vStart;
163*f7caa48eSbarral         for (i=0; i<dim+1; ++i) {
164*f7caa48eSbarral           if ( iVer == cell[i] ) {
165*f7caa48eSbarral             perm[cl] = i;    // the cl-th vertex of the closure is the i-th vertex of the ccell
166*f7caa48eSbarral             break;
167*f7caa48eSbarral           }
168*f7caa48eSbarral         }
169*f7caa48eSbarral       }
170*f7caa48eSbarral       ierr = DMPlexGetCone(dm->coarseMesh, c, &facetList);CHKERRQ(ierr);     // list of edges/facets of the cell
171*f7caa48eSbarral       ierr = DMPlexGetConeSize(dm->coarseMesh, c, &facetListSize);CHKERRQ(ierr);
172*f7caa48eSbarral       // then, for each edge/facet of the cell, find the opposite vertex (ie the one not in the closure of the facet/edge)
173*f7caa48eSbarral       for (f=0; f<facetListSize; ++f){
174*f7caa48eSbarral         facet = facetList[f];
175*f7caa48eSbarral         ierr = DMPlexGetTransitiveClosure(dm->coarseMesh, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr);
176*f7caa48eSbarral         // loop over the vertices of the cell closure, and check if each vertex belongs to the edge/facet closure
177*f7caa48eSbarral         PetscMemzero(isInFacetClosure, sizeof(isInFacetClosure));
178*f7caa48eSbarral         for (cl = 0; cl < cellClosureSize*2; cl+=2) {
179*f7caa48eSbarral           if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) isInFacetClosure[cl] = 1;
180*f7caa48eSbarral           for (cl2 =0; cl2<facetClosureSize*2; cl2+=2){
181*f7caa48eSbarral             if ((facetClosure[cl2] < vStart) || (facetClosure[cl2] >= vEnd)) continue;
182*f7caa48eSbarral             if ( cellClosure[cl] == facetClosure[cl2] ) {
183*f7caa48eSbarral               isInFacetClosure[cl] = 1;
184*f7caa48eSbarral             }
185*f7caa48eSbarral           }
186*f7caa48eSbarral         }
187*f7caa48eSbarral         // 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
188*f7caa48eSbarral         for (cl = 0; cl < cellClosureSize*2; cl+=2 ) {
189*f7caa48eSbarral           if ( !isInFacetClosure[cl] ) {
190*f7caa48eSbarral             idx = (c-cStart)*(dim+1) + perm[cl];
191*f7caa48eSbarral             if ( boundaryTags[idx] )
192*f7caa48eSbarral               ierr = DMSetLabelValue(dm->coarseMesh, bdyLabelName, facet, boundaryTags[idx]);CHKERRQ(ierr);
193*f7caa48eSbarral             break;
194*f7caa48eSbarral           }
195*f7caa48eSbarral         }
196*f7caa48eSbarral         ierr = DMPlexRestoreTransitiveClosure(dm->coarseMesh, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr);
197*f7caa48eSbarral       }
198*f7caa48eSbarral       ierr = DMPlexRestoreTransitiveClosure(dm->coarseMesh, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr);
199*f7caa48eSbarral     }
2000bbe5d31Sbarral     pragmatic_finalize();
2010bbe5d31Sbarral     ierr = PetscFree4(x, y, z, cells);CHKERRQ(ierr);
2020bbe5d31Sbarral     ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
2030bbe5d31Sbarral     ierr = PetscFree(coarseCoords);CHKERRQ(ierr);
2040bbe5d31Sbarral   }
2050bbe5d31Sbarral #endif
2060bbe5d31Sbarral   ierr = PetscObjectReference((PetscObject) dm->coarseMesh);CHKERRQ(ierr);
2070bbe5d31Sbarral   *dmCoarsened = dm->coarseMesh;
2080bbe5d31Sbarral   PetscFunctionReturn(0);
2090bbe5d31Sbarral }
210