xref: /petsc/src/dm/impls/plex/adaptors/pragmatic/pragmaticadapt.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
12bc0b47dSJoe Wallwork #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
22bc0b47dSJoe Wallwork #include <pragmatic/cpragmatic.h>
32bc0b47dSJoe Wallwork 
49fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric_Pragmatic_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew)
52bc0b47dSJoe Wallwork {
62bc0b47dSJoe Wallwork   MPI_Comm           comm;
72bc0b47dSJoe Wallwork   const char        *bdName = "_boundary_";
82bc0b47dSJoe Wallwork #if 0
92bc0b47dSJoe Wallwork   DM                 odm = dm;
102bc0b47dSJoe Wallwork #endif
112bc0b47dSJoe Wallwork   DM                 udm, cdm;
122bc0b47dSJoe Wallwork   DMLabel            bdLabelFull;
132bc0b47dSJoe Wallwork   const char        *bdLabelName;
142bc0b47dSJoe Wallwork   IS                 bdIS, globalVertexNum;
152bc0b47dSJoe Wallwork   PetscSection       coordSection;
162bc0b47dSJoe Wallwork   Vec                coordinates;
172bc0b47dSJoe Wallwork   const PetscScalar *coords, *met;
182bc0b47dSJoe Wallwork   const PetscInt    *bdFacesFull, *gV;
192bc0b47dSJoe Wallwork   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
202bc0b47dSJoe Wallwork   PetscReal         *x, *y, *z, *metric;
212bc0b47dSJoe Wallwork   PetscInt          *cells;
222bc0b47dSJoe Wallwork   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
2393520041SJoe Wallwork   PetscInt           off, maxConeSize, numBdFaces, f, bdSize, i, j, Nd;
2493520041SJoe Wallwork   PetscBool          flg, isotropic, uniform;
252bc0b47dSJoe Wallwork   DMLabel            bdLabelNew;
262bc0b47dSJoe Wallwork   PetscReal         *coordsNew;
272bc0b47dSJoe Wallwork   PetscInt          *bdTags;
282bc0b47dSJoe Wallwork   PetscReal         *xNew[3] = {NULL, NULL, NULL};
292bc0b47dSJoe Wallwork   PetscInt          *cellsNew;
302bc0b47dSJoe Wallwork   PetscInt           d, numCellsNew, numVerticesNew;
312bc0b47dSJoe Wallwork   PetscInt           numCornersNew, fStart, fEnd;
322bc0b47dSJoe Wallwork   PetscMPIInt        numProcs;
332bc0b47dSJoe Wallwork 
342bc0b47dSJoe Wallwork   PetscFunctionBegin;
352bc0b47dSJoe Wallwork 
362bc0b47dSJoe Wallwork   /* Check for FEM adjacency flags */
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
38*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &numProcs));
392bc0b47dSJoe Wallwork   if (bdLabel) {
40*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetName((PetscObject) bdLabel, &bdLabelName));
41*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(bdLabelName, bdName, &flg));
42*5f80ce2aSJacob Faibussowitsch     PetscCheck(!flg,comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
432bc0b47dSJoe Wallwork   }
44*5f80ce2aSJacob Faibussowitsch   PetscCheck(!rgLabel,comm, PETSC_ERR_ARG_WRONG, "Cannot currently preserve cell tags with Pragmatic");
452bc0b47dSJoe Wallwork #if 0
462bc0b47dSJoe Wallwork   /* Check for overlap by looking for cell in the SF */
472bc0b47dSJoe Wallwork   if (!overlapped) {
48*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexDistributeOverlap(odm, 1, NULL, &dm));
49*5f80ce2aSJacob Faibussowitsch     if (!dm) {dm = odm; CHKERRQ(PetscObjectReference((PetscObject) dm));}
502bc0b47dSJoe Wallwork   }
512bc0b47dSJoe Wallwork #endif
522bc0b47dSJoe Wallwork 
532bc0b47dSJoe Wallwork   /* Get mesh information */
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexUninterpolate(dm, &udm));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetMaxSizes(udm, &maxConeSize, NULL));
592bc0b47dSJoe Wallwork   numCells = cEnd - cStart;
602bc0b47dSJoe Wallwork   if (numCells == 0) {
612bc0b47dSJoe Wallwork     PetscMPIInt rank;
622bc0b47dSJoe Wallwork 
63*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(comm, &rank));
6498921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot perform mesh adaptation because process %d does not own any cells.", rank);
652bc0b47dSJoe Wallwork   }
662bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells));
682bc0b47dSJoe Wallwork 
692bc0b47dSJoe Wallwork   /* Get cell offsets */
702bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
712bc0b47dSJoe Wallwork     const PetscInt *cone;
722bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
732bc0b47dSJoe Wallwork 
74*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetConeSize(udm, c, &coneSize));
75*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetCone(udm, c, &cone));
762bc0b47dSJoe Wallwork     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
772bc0b47dSJoe Wallwork   }
782bc0b47dSJoe Wallwork 
792bc0b47dSJoe Wallwork   /* Get local-to-global vertex map */
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(numVertices, &l2gv));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetVertexNumbering(udm, &globalVertexNum));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(globalVertexNum, &gV));
832bc0b47dSJoe Wallwork   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
842bc0b47dSJoe Wallwork     if (gV[v] >= 0) ++numLocVertices;
852bc0b47dSJoe Wallwork     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
862bc0b47dSJoe Wallwork   }
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(globalVertexNum, &gV));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&udm));
892bc0b47dSJoe Wallwork 
902bc0b47dSJoe Wallwork   /* Get vertex coordinate arrays */
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(dm, &cdm));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(cdm, &coordSection));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(coordinates, &coords));
952bc0b47dSJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
96*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(coordSection, v, &off));
972bc0b47dSJoe Wallwork     x[v-vStart] = PetscRealPart(coords[off+0]);
982bc0b47dSJoe Wallwork     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
992bc0b47dSJoe Wallwork     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
1002bc0b47dSJoe Wallwork   }
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(coordinates, &coords));
1022bc0b47dSJoe Wallwork 
1032bc0b47dSJoe Wallwork   /* Get boundary mesh */
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelGetStratumIS(bdLabelFull, 1, &bdIS));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(bdIS, &bdFacesFull));
1092bc0b47dSJoe Wallwork   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1102bc0b47dSJoe Wallwork     PetscInt *closure = NULL;
1112bc0b47dSJoe Wallwork     PetscInt  closureSize, cl;
1122bc0b47dSJoe Wallwork 
113*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
1142bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1152bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1162bc0b47dSJoe Wallwork     }
117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
1182bc0b47dSJoe Wallwork   }
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds));
1202bc0b47dSJoe Wallwork   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1212bc0b47dSJoe Wallwork     PetscInt *closure = NULL;
1222bc0b47dSJoe Wallwork     PetscInt  closureSize, cl;
1232bc0b47dSJoe Wallwork 
124*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
1252bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1262bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
1272bc0b47dSJoe Wallwork     }
128*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure));
129*5f80ce2aSJacob Faibussowitsch     if (bdLabel) CHKERRQ(DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]));
1302bc0b47dSJoe Wallwork     else         {bdFaceIds[f] = 1;}
1312bc0b47dSJoe Wallwork   }
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&bdIS));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&bdLabelFull));
1342bc0b47dSJoe Wallwork 
1352bc0b47dSJoe Wallwork   /* Get metric */
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view"));
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(vertexMetric, &met));
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricIsIsotropic(dm, &isotropic));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexMetricIsUniform(dm, &uniform));
14093520041SJoe Wallwork   Nd = PetscSqr(dim);
14193520041SJoe Wallwork   for (v = 0; v < vEnd-vStart; ++v) {
14293520041SJoe Wallwork     for (i = 0; i < dim; ++i) {
14393520041SJoe Wallwork       for (j = 0; j < dim; ++j) {
14493520041SJoe Wallwork         if (isotropic) {
14593520041SJoe Wallwork           if (i == j) {
14693520041SJoe Wallwork             if (uniform) metric[Nd*v+dim*i+j] = PetscRealPart(met[0]);
14793520041SJoe Wallwork             else metric[Nd*v+dim*i+j] = PetscRealPart(met[v]);
14893520041SJoe Wallwork           } else metric[Nd*v+dim*i+j] = 0.0;
14993520041SJoe Wallwork         } else metric[Nd*v+dim*i+j] = PetscRealPart(met[Nd*v+dim*i+j]);
15093520041SJoe Wallwork       }
15193520041SJoe Wallwork     }
15293520041SJoe Wallwork   }
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(vertexMetric, &met));
1542bc0b47dSJoe Wallwork 
1552bc0b47dSJoe Wallwork #if 0
1562bc0b47dSJoe Wallwork   /* Destroy overlap mesh */
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
1582bc0b47dSJoe Wallwork #endif
1592bc0b47dSJoe Wallwork   /* Send to Pragmatic and remesh */
1602bc0b47dSJoe Wallwork   switch (dim) {
1612bc0b47dSJoe Wallwork   case 2:
1622bc0b47dSJoe Wallwork     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);
1632bc0b47dSJoe Wallwork     break;
1642bc0b47dSJoe Wallwork   case 3:
1652bc0b47dSJoe Wallwork     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);
1662bc0b47dSJoe Wallwork     break;
16798921bdaSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim);
1682bc0b47dSJoe Wallwork   }
1692bc0b47dSJoe Wallwork   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
1702bc0b47dSJoe Wallwork   pragmatic_set_metric(metric);
1712bc0b47dSJoe Wallwork   pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0);
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(l2gv));
1732bc0b47dSJoe Wallwork 
1742bc0b47dSJoe Wallwork   /* Retrieve mesh from Pragmatic and create new plex */
1752bc0b47dSJoe Wallwork   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numVerticesNew*dim, &coordsNew));
1772bc0b47dSJoe Wallwork   switch (dim) {
1782bc0b47dSJoe Wallwork   case 2:
1792bc0b47dSJoe Wallwork     numCornersNew = 3;
180*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]));
1812bc0b47dSJoe Wallwork     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
1822bc0b47dSJoe Wallwork     break;
1832bc0b47dSJoe Wallwork   case 3:
1842bc0b47dSJoe Wallwork     numCornersNew = 4;
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]));
1862bc0b47dSJoe Wallwork     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
1872bc0b47dSJoe Wallwork     break;
1882bc0b47dSJoe Wallwork   default:
18998921bdaSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim);
1902bc0b47dSJoe Wallwork   }
1912bc0b47dSJoe Wallwork   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];}
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numCellsNew*(dim+1), &cellsNew));
1932bc0b47dSJoe Wallwork   pragmatic_get_elements(cellsNew);
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, NULL, dmNew));
1952bc0b47dSJoe Wallwork 
1962bc0b47dSJoe Wallwork   /* Rebuild boundary label */
1972bc0b47dSJoe Wallwork   pragmatic_get_boundaryTags(&bdTags);
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName));
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
2032bc0b47dSJoe Wallwork   for (c = cStart; c < cEnd; ++c) {
2042bc0b47dSJoe Wallwork 
2052bc0b47dSJoe Wallwork     /* Only for simplicial meshes */
2062bc0b47dSJoe Wallwork     coff = (c-cStart)*(dim+1);
2072bc0b47dSJoe Wallwork 
2082bc0b47dSJoe Wallwork     /* d is the local cell number of the vertex opposite to the face we are marking */
2092bc0b47dSJoe Wallwork     for (d = 0; d < dim+1; ++d) {
2102bc0b47dSJoe Wallwork       if (bdTags[coff+d]) {
2112bc0b47dSJoe Wallwork         const PetscInt  perm[4][4] = {{-1, -1, -1, -1}, {-1, -1, -1, -1}, {1, 2, 0, -1}, {3, 2, 1, 0}}; /* 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() */
215*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexGetCone(*dmNew, c, &cone));
216*5f80ce2aSJacob Faibussowitsch         CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   case 2: CHKERRQ(PetscFree2(xNew[0], xNew[1]));
2242bc0b47dSJoe Wallwork   break;
225*5f80ce2aSJacob Faibussowitsch   case 3: CHKERRQ(PetscFree3(xNew[0], xNew[1], xNew[2]));
2262bc0b47dSJoe Wallwork   break;
2272bc0b47dSJoe Wallwork   }
228*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cellsNew));
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(x, y, z, metric, cells));
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(bdFaces, bdFaceIds));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(coordsNew));
2322bc0b47dSJoe Wallwork   pragmatic_finalize();
2332bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
2342bc0b47dSJoe Wallwork }
235