xref: /petsc/src/dm/impls/plex/adaptors/pragmatic/pragmaticadapt.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
12bc0b47dSJoe Wallwork #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
22bc0b47dSJoe Wallwork #include <pragmatic/cpragmatic.h>
32bc0b47dSJoe Wallwork 
4*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMAdaptMetric_Pragmatic_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew) {
52bc0b47dSJoe Wallwork   MPI_Comm    comm;
62bc0b47dSJoe Wallwork   const char *bdName = "_boundary_";
72bc0b47dSJoe Wallwork #if 0
82bc0b47dSJoe Wallwork   DM                 odm = dm;
92bc0b47dSJoe Wallwork #endif
102bc0b47dSJoe Wallwork   DM                 udm, cdm;
112bc0b47dSJoe Wallwork   DMLabel            bdLabelFull;
122bc0b47dSJoe Wallwork   const char        *bdLabelName;
132bc0b47dSJoe Wallwork   IS                 bdIS, globalVertexNum;
142bc0b47dSJoe Wallwork   PetscSection       coordSection;
152bc0b47dSJoe Wallwork   Vec                coordinates;
162bc0b47dSJoe Wallwork   const PetscScalar *coords, *met;
172bc0b47dSJoe Wallwork   const PetscInt    *bdFacesFull, *gV;
182bc0b47dSJoe Wallwork   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
192bc0b47dSJoe Wallwork   PetscReal         *x, *y, *z, *metric;
202bc0b47dSJoe Wallwork   PetscInt          *cells;
212bc0b47dSJoe Wallwork   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
2293520041SJoe Wallwork   PetscInt           off, maxConeSize, numBdFaces, f, bdSize, i, j, Nd;
2393520041SJoe Wallwork   PetscBool          flg, isotropic, uniform;
242bc0b47dSJoe Wallwork   DMLabel            bdLabelNew;
252bc0b47dSJoe Wallwork   PetscReal         *coordsNew;
262bc0b47dSJoe Wallwork   PetscInt          *bdTags;
272bc0b47dSJoe Wallwork   PetscReal         *xNew[3] = {NULL, NULL, NULL};
282bc0b47dSJoe Wallwork   PetscInt          *cellsNew;
292bc0b47dSJoe Wallwork   PetscInt           d, numCellsNew, numVerticesNew;
302bc0b47dSJoe Wallwork   PetscInt           numCornersNew, fStart, fEnd;
312bc0b47dSJoe Wallwork   PetscMPIInt        numProcs;
322bc0b47dSJoe Wallwork 
332bc0b47dSJoe Wallwork   PetscFunctionBegin;
342bc0b47dSJoe Wallwork 
352bc0b47dSJoe Wallwork   /* Check for FEM adjacency flags */
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &numProcs));
382bc0b47dSJoe Wallwork   if (bdLabel) {
399566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)bdLabel, &bdLabelName));
409566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(bdLabelName, bdName, &flg));
415f80ce2aSJacob Faibussowitsch     PetscCheck(!flg, comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
422bc0b47dSJoe Wallwork   }
435f80ce2aSJacob Faibussowitsch   PetscCheck(!rgLabel, comm, PETSC_ERR_ARG_WRONG, "Cannot currently preserve cell tags with Pragmatic");
442bc0b47dSJoe Wallwork #if 0
452bc0b47dSJoe Wallwork   /* Check for overlap by looking for cell in the SF */
462bc0b47dSJoe Wallwork   if (!overlapped) {
479566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeOverlap(odm, 1, NULL, &dm));
489566063dSJacob Faibussowitsch     if (!dm) {dm = odm; PetscCall(PetscObjectReference((PetscObject) dm));}
492bc0b47dSJoe Wallwork   }
502bc0b47dSJoe Wallwork #endif
512bc0b47dSJoe Wallwork 
522bc0b47dSJoe Wallwork   /* Get mesh information */
539566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
549566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
559566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
569566063dSJacob Faibussowitsch   PetscCall(DMPlexUninterpolate(dm, &udm));
579566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
582bc0b47dSJoe Wallwork   numCells = cEnd - cStart;
592bc0b47dSJoe Wallwork   if (numCells == 0) {
602bc0b47dSJoe Wallwork     PetscMPIInt rank;
612bc0b47dSJoe Wallwork 
629566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &rank));
6398921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot perform mesh adaptation because process %d does not own any cells.", rank);
642bc0b47dSJoe Wallwork   }
652bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
669566063dSJacob Faibussowitsch   PetscCall(PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices * PetscSqr(dim), &metric, numCells * maxConeSize, &cells));
672bc0b47dSJoe Wallwork 
682bc0b47dSJoe Wallwork   /* Get cell offsets */
692bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
702bc0b47dSJoe Wallwork     const PetscInt *cone;
712bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
722bc0b47dSJoe Wallwork 
739566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(udm, c, &coneSize));
749566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(udm, c, &cone));
752bc0b47dSJoe Wallwork     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
762bc0b47dSJoe Wallwork   }
772bc0b47dSJoe Wallwork 
782bc0b47dSJoe Wallwork   /* Get local-to-global vertex map */
799566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(numVertices, &l2gv));
809566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVertexNumbering(udm, &globalVertexNum));
819566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(globalVertexNum, &gV));
822bc0b47dSJoe Wallwork   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
832bc0b47dSJoe Wallwork     if (gV[v] >= 0) ++numLocVertices;
842bc0b47dSJoe Wallwork     l2gv[v] = gV[v] < 0 ? -(gV[v] + 1) : gV[v];
852bc0b47dSJoe Wallwork   }
869566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(globalVertexNum, &gV));
879566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&udm));
882bc0b47dSJoe Wallwork 
892bc0b47dSJoe Wallwork   /* Get vertex coordinate arrays */
909566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
919566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(cdm, &coordSection));
929566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
942bc0b47dSJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
959566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(coordSection, v, &off));
962bc0b47dSJoe Wallwork     x[v - vStart] = PetscRealPart(coords[off + 0]);
972bc0b47dSJoe Wallwork     if (dim > 1) y[v - vStart] = PetscRealPart(coords[off + 1]);
982bc0b47dSJoe Wallwork     if (dim > 2) z[v - vStart] = PetscRealPart(coords[off + 2]);
992bc0b47dSJoe Wallwork   }
1009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
1012bc0b47dSJoe Wallwork 
1022bc0b47dSJoe Wallwork   /* Get boundary mesh */
1039566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull));
1049566063dSJacob Faibussowitsch   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull));
1059566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(bdLabelFull, 1, &bdIS));
1069566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces));
1079566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(bdIS, &bdFacesFull));
1082bc0b47dSJoe Wallwork   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1092bc0b47dSJoe Wallwork     PetscInt *closure = NULL;
1102bc0b47dSJoe Wallwork     PetscInt  closureSize, cl;
1112bc0b47dSJoe Wallwork 
1129566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
1132bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize * 2; cl += 2) {
1142bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1152bc0b47dSJoe Wallwork     }
1169566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
1172bc0b47dSJoe Wallwork   }
1189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds));
1192bc0b47dSJoe Wallwork   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1202bc0b47dSJoe Wallwork     PetscInt *closure = NULL;
1212bc0b47dSJoe Wallwork     PetscInt  closureSize, cl;
1222bc0b47dSJoe Wallwork 
1239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
1242bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize * 2; cl += 2) {
1252bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
1262bc0b47dSJoe Wallwork     }
1279566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
1289566063dSJacob Faibussowitsch     if (bdLabel) PetscCall(DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]));
1292bc0b47dSJoe Wallwork     else { bdFaceIds[f] = 1; }
1302bc0b47dSJoe Wallwork   }
1319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&bdIS));
1329566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&bdLabelFull));
1332bc0b47dSJoe Wallwork 
1342bc0b47dSJoe Wallwork   /* Get metric */
1359566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
1369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(vertexMetric, &met));
1379566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
1389566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
13993520041SJoe Wallwork   Nd = PetscSqr(dim);
14093520041SJoe Wallwork   for (v = 0; v < vEnd - vStart; ++v) {
14193520041SJoe Wallwork     for (i = 0; i < dim; ++i) {
14293520041SJoe Wallwork       for (j = 0; j < dim; ++j) {
14393520041SJoe Wallwork         if (isotropic) {
14493520041SJoe Wallwork           if (i == j) {
14593520041SJoe Wallwork             if (uniform) metric[Nd * v + dim * i + j] = PetscRealPart(met[0]);
14693520041SJoe Wallwork             else metric[Nd * v + dim * i + j] = PetscRealPart(met[v]);
14793520041SJoe Wallwork           } else metric[Nd * v + dim * i + j] = 0.0;
14893520041SJoe Wallwork         } else metric[Nd * v + dim * i + j] = PetscRealPart(met[Nd * v + dim * i + j]);
14993520041SJoe Wallwork       }
15093520041SJoe Wallwork     }
15193520041SJoe Wallwork   }
1529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(vertexMetric, &met));
1532bc0b47dSJoe Wallwork 
1542bc0b47dSJoe Wallwork #if 0
1552bc0b47dSJoe Wallwork   /* Destroy overlap mesh */
1569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1572bc0b47dSJoe Wallwork #endif
1582bc0b47dSJoe Wallwork   /* Send to Pragmatic and remesh */
1592bc0b47dSJoe Wallwork   switch (dim) {
160*9371c9d4SSatish Balay   case 2: pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm); break;
161*9371c9d4SSatish Balay   case 3: pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm); break;
16263a3b9bcSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %" PetscInt_FMT, dim);
1632bc0b47dSJoe Wallwork   }
1642bc0b47dSJoe Wallwork   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
1652bc0b47dSJoe Wallwork   pragmatic_set_metric(metric);
1662bc0b47dSJoe Wallwork   pragmatic_adapt(((DM_Plex *)dm->data)->remeshBd ? 1 : 0);
1679566063dSJacob Faibussowitsch   PetscCall(PetscFree(l2gv));
1682bc0b47dSJoe Wallwork 
1692bc0b47dSJoe Wallwork   /* Retrieve mesh from Pragmatic and create new plex */
1702bc0b47dSJoe Wallwork   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
1719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numVerticesNew * dim, &coordsNew));
1722bc0b47dSJoe Wallwork   switch (dim) {
1732bc0b47dSJoe Wallwork   case 2:
1742bc0b47dSJoe Wallwork     numCornersNew = 3;
1759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]));
1762bc0b47dSJoe Wallwork     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
1772bc0b47dSJoe Wallwork     break;
1782bc0b47dSJoe Wallwork   case 3:
1792bc0b47dSJoe Wallwork     numCornersNew = 4;
1809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]));
1812bc0b47dSJoe Wallwork     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
1822bc0b47dSJoe Wallwork     break;
183*9371c9d4SSatish Balay   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %" PetscInt_FMT, dim);
1842bc0b47dSJoe Wallwork   }
185*9371c9d4SSatish Balay   for (v = 0; v < numVerticesNew; ++v) {
186*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) coordsNew[v * dim + d] = xNew[d][v];
187*9371c9d4SSatish Balay   }
1889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numCellsNew * (dim + 1), &cellsNew));
1892bc0b47dSJoe Wallwork   pragmatic_get_elements(cellsNew);
1909566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, NULL, dmNew));
1912bc0b47dSJoe Wallwork 
1922bc0b47dSJoe Wallwork   /* Rebuild boundary label */
1932bc0b47dSJoe Wallwork   pragmatic_get_boundaryTags(&bdTags);
1949566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName));
1959566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew));
1969566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
1979566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
1989566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
1992bc0b47dSJoe Wallwork   for (c = cStart; c < cEnd; ++c) {
2002bc0b47dSJoe Wallwork     /* Only for simplicial meshes */
2012bc0b47dSJoe Wallwork     coff = (c - cStart) * (dim + 1);
2022bc0b47dSJoe Wallwork 
2032bc0b47dSJoe Wallwork     /* d is the local cell number of the vertex opposite to the face we are marking */
2042bc0b47dSJoe Wallwork     for (d = 0; d < dim + 1; ++d) {
2052bc0b47dSJoe Wallwork       if (bdTags[coff + d]) {
206*9371c9d4SSatish Balay         const PetscInt perm[4][4] = {
207*9371c9d4SSatish Balay           {-1, -1, -1, -1},
208*9371c9d4SSatish Balay           {-1, -1, -1, -1},
209*9371c9d4SSatish Balay           {1,  2,  0,  -1},
210*9371c9d4SSatish Balay           {3,  2,  1,  0 }
211*9371c9d4SSatish Balay         }; /* perm[d] = face opposite */
2122bc0b47dSJoe Wallwork         const PetscInt *cone;
2132bc0b47dSJoe Wallwork 
2142bc0b47dSJoe Wallwork         /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
2159566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCone(*dmNew, c, &cone));
2169566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff + d]));
2172bc0b47dSJoe Wallwork       }
2182bc0b47dSJoe Wallwork     }
2192bc0b47dSJoe Wallwork   }
2202bc0b47dSJoe Wallwork 
2212bc0b47dSJoe Wallwork   /* Clean up */
2222bc0b47dSJoe Wallwork   switch (dim) {
223*9371c9d4SSatish Balay   case 2: PetscCall(PetscFree2(xNew[0], xNew[1])); break;
224*9371c9d4SSatish Balay   case 3: PetscCall(PetscFree3(xNew[0], xNew[1], xNew[2])); break;
2252bc0b47dSJoe Wallwork   }
2269566063dSJacob Faibussowitsch   PetscCall(PetscFree(cellsNew));
2279566063dSJacob Faibussowitsch   PetscCall(PetscFree5(x, y, z, metric, cells));
2289566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bdFaces, bdFaceIds));
2299566063dSJacob Faibussowitsch   PetscCall(PetscFree(coordsNew));
2302bc0b47dSJoe Wallwork   pragmatic_finalize();
2312bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
2322bc0b47dSJoe Wallwork }
233