12bc0b47dSJoe Wallwork #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
22bc0b47dSJoe Wallwork #include <pragmatic/cpragmatic.h>
32bc0b47dSJoe Wallwork
DMAdaptMetric_Pragmatic_Plex(DM dm,Vec vertexMetric,DMLabel bdLabel,DMLabel rgLabel,DM * dmNew)4d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMAdaptMetric_Pragmatic_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DMLabel rgLabel, DM *dmNew)
5d71ae5a4SJacob Faibussowitsch {
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 /* 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]));
129ad540459SPierre Jolivet 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 */
159*d84739d8SStefano Zampini PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1602bc0b47dSJoe Wallwork switch (dim) {
161d71ae5a4SJacob Faibussowitsch case 2:
162*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_2d_mpi_init", pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm));
163d71ae5a4SJacob Faibussowitsch break;
164d71ae5a4SJacob Faibussowitsch case 3:
165*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_3d_mpi_init", pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm));
166d71ae5a4SJacob Faibussowitsch break;
167d71ae5a4SJacob Faibussowitsch default:
168d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %" PetscInt_FMT, dim);
1692bc0b47dSJoe Wallwork }
170*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_set_boundary", pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds));
171*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_set_metric", pragmatic_set_metric(metric));
172*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_adapt", pragmatic_adapt(((DM_Plex *)dm->data)->remeshBd ? 1 : 0));
1739566063dSJacob Faibussowitsch PetscCall(PetscFree(l2gv));
174*d84739d8SStefano Zampini PetscCall(PetscFPTrapPop());
1752bc0b47dSJoe Wallwork
1762bc0b47dSJoe Wallwork /* Retrieve mesh from Pragmatic and create new plex */
177*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_get_info_mpi", pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew));
1789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numVerticesNew * dim, &coordsNew));
1792bc0b47dSJoe Wallwork switch (dim) {
1802bc0b47dSJoe Wallwork case 2:
1812bc0b47dSJoe Wallwork numCornersNew = 3;
1829566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]));
183*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_get_coords_2d_mpi", pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]));
1842bc0b47dSJoe Wallwork break;
1852bc0b47dSJoe Wallwork case 3:
1862bc0b47dSJoe Wallwork numCornersNew = 4;
1879566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]));
188*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_get_coords_3d_mpi", pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]));
1892bc0b47dSJoe Wallwork break;
190d71ae5a4SJacob Faibussowitsch default:
191d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %" PetscInt_FMT, dim);
1922bc0b47dSJoe Wallwork }
1939371c9d4SSatish Balay for (v = 0; v < numVerticesNew; ++v) {
1949371c9d4SSatish Balay for (d = 0; d < dim; ++d) coordsNew[v * dim + d] = xNew[d][v];
1959371c9d4SSatish Balay }
1969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCellsNew * (dim + 1), &cellsNew));
197*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_get_elements", pragmatic_get_elements(cellsNew));
1989566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, NULL, dmNew));
1992bc0b47dSJoe Wallwork
2002bc0b47dSJoe Wallwork /* Rebuild boundary label */
201*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_get_boundaryTags", pragmatic_get_boundaryTags(&bdTags));
2029566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName));
2039566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew));
2049566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd));
2059566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd));
2069566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd));
2072bc0b47dSJoe Wallwork for (c = cStart; c < cEnd; ++c) {
2082bc0b47dSJoe Wallwork /* Only for simplicial meshes */
2092bc0b47dSJoe Wallwork coff = (c - cStart) * (dim + 1);
2102bc0b47dSJoe Wallwork
2112bc0b47dSJoe Wallwork /* d is the local cell number of the vertex opposite to the face we are marking */
2122bc0b47dSJoe Wallwork for (d = 0; d < dim + 1; ++d) {
2132bc0b47dSJoe Wallwork if (bdTags[coff + d]) {
2149371c9d4SSatish Balay const PetscInt perm[4][4] = {
2159371c9d4SSatish Balay {-1, -1, -1, -1},
2169371c9d4SSatish Balay {-1, -1, -1, -1},
2179371c9d4SSatish Balay {1, 2, 0, -1},
2189371c9d4SSatish Balay {3, 2, 1, 0 }
2199371c9d4SSatish Balay }; /* perm[d] = face opposite */
2202bc0b47dSJoe Wallwork const PetscInt *cone;
2212bc0b47dSJoe Wallwork
2222bc0b47dSJoe Wallwork /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
2239566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(*dmNew, c, &cone));
2249566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff + d]));
2252bc0b47dSJoe Wallwork }
2262bc0b47dSJoe Wallwork }
2272bc0b47dSJoe Wallwork }
2282bc0b47dSJoe Wallwork
2292bc0b47dSJoe Wallwork /* Clean up */
2302bc0b47dSJoe Wallwork switch (dim) {
231d71ae5a4SJacob Faibussowitsch case 2:
232d71ae5a4SJacob Faibussowitsch PetscCall(PetscFree2(xNew[0], xNew[1]));
233d71ae5a4SJacob Faibussowitsch break;
234d71ae5a4SJacob Faibussowitsch case 3:
235d71ae5a4SJacob Faibussowitsch PetscCall(PetscFree3(xNew[0], xNew[1], xNew[2]));
236d71ae5a4SJacob Faibussowitsch break;
2372bc0b47dSJoe Wallwork }
2389566063dSJacob Faibussowitsch PetscCall(PetscFree(cellsNew));
2399566063dSJacob Faibussowitsch PetscCall(PetscFree5(x, y, z, metric, cells));
2409566063dSJacob Faibussowitsch PetscCall(PetscFree2(bdFaces, bdFaceIds));
2419566063dSJacob Faibussowitsch PetscCall(PetscFree(coordsNew));
242*d84739d8SStefano Zampini PetscStackCallExternalVoid("pragmatic_finalize", pragmatic_finalize());
2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2442bc0b47dSJoe Wallwork }
245