xref: /petsc/src/dm/impls/plex/adaptors/pragmatic/pragmaticadapt.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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   PetscErrorCode     ierr;
342bc0b47dSJoe Wallwork 
352bc0b47dSJoe Wallwork   PetscFunctionBegin;
362bc0b47dSJoe Wallwork 
372bc0b47dSJoe Wallwork   /* Check for FEM adjacency flags */
382bc0b47dSJoe Wallwork   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
392bc0b47dSJoe Wallwork   ierr = MPI_Comm_size(comm, &numProcs);CHKERRMPI(ierr);
402bc0b47dSJoe Wallwork   if (bdLabel) {
412bc0b47dSJoe Wallwork     ierr = PetscObjectGetName((PetscObject) bdLabel, &bdLabelName);CHKERRQ(ierr);
422bc0b47dSJoe Wallwork     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
43*98921bdaSJacob Faibussowitsch     if (flg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
442bc0b47dSJoe Wallwork   }
459fe9e680SJoe Wallwork   if (rgLabel) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot currently preserve cell tags with Pragmatic");
462bc0b47dSJoe Wallwork #if 0
472bc0b47dSJoe Wallwork   /* Check for overlap by looking for cell in the SF */
482bc0b47dSJoe Wallwork   if (!overlapped) {
492bc0b47dSJoe Wallwork     ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr);
502bc0b47dSJoe Wallwork     if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);}
512bc0b47dSJoe Wallwork   }
522bc0b47dSJoe Wallwork #endif
532bc0b47dSJoe Wallwork 
542bc0b47dSJoe Wallwork   /* Get mesh information */
552bc0b47dSJoe Wallwork   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
562bc0b47dSJoe Wallwork   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
572bc0b47dSJoe Wallwork   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
582bc0b47dSJoe Wallwork   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
592bc0b47dSJoe Wallwork   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
602bc0b47dSJoe Wallwork   numCells = cEnd - cStart;
612bc0b47dSJoe Wallwork   if (numCells == 0) {
622bc0b47dSJoe Wallwork     PetscMPIInt rank;
632bc0b47dSJoe Wallwork 
642bc0b47dSJoe Wallwork     ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
65*98921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot perform mesh adaptation because process %d does not own any cells.", rank);
662bc0b47dSJoe Wallwork   }
672bc0b47dSJoe Wallwork   numVertices = vEnd - vStart;
682bc0b47dSJoe Wallwork   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
692bc0b47dSJoe Wallwork 
702bc0b47dSJoe Wallwork   /* Get cell offsets */
712bc0b47dSJoe Wallwork   for (c = 0, coff = 0; c < numCells; ++c) {
722bc0b47dSJoe Wallwork     const PetscInt *cone;
732bc0b47dSJoe Wallwork     PetscInt        coneSize, cl;
742bc0b47dSJoe Wallwork 
752bc0b47dSJoe Wallwork     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
762bc0b47dSJoe Wallwork     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
772bc0b47dSJoe Wallwork     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
782bc0b47dSJoe Wallwork   }
792bc0b47dSJoe Wallwork 
802bc0b47dSJoe Wallwork   /* Get local-to-global vertex map */
812bc0b47dSJoe Wallwork   ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr);
822bc0b47dSJoe Wallwork   ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
832bc0b47dSJoe Wallwork   ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
842bc0b47dSJoe Wallwork   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
852bc0b47dSJoe Wallwork     if (gV[v] >= 0) ++numLocVertices;
862bc0b47dSJoe Wallwork     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
872bc0b47dSJoe Wallwork   }
882bc0b47dSJoe Wallwork   ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
892bc0b47dSJoe Wallwork   ierr = DMDestroy(&udm);CHKERRQ(ierr);
902bc0b47dSJoe Wallwork 
912bc0b47dSJoe Wallwork   /* Get vertex coordinate arrays */
922bc0b47dSJoe Wallwork   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
932bc0b47dSJoe Wallwork   ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
942bc0b47dSJoe Wallwork   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
952bc0b47dSJoe Wallwork   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
962bc0b47dSJoe Wallwork   for (v = vStart; v < vEnd; ++v) {
972bc0b47dSJoe Wallwork     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
982bc0b47dSJoe Wallwork     x[v-vStart] = PetscRealPart(coords[off+0]);
992bc0b47dSJoe Wallwork     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
1002bc0b47dSJoe Wallwork     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
1012bc0b47dSJoe Wallwork   }
1022bc0b47dSJoe Wallwork   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
1032bc0b47dSJoe Wallwork 
1042bc0b47dSJoe Wallwork   /* Get boundary mesh */
1052bc0b47dSJoe Wallwork   ierr = DMLabelCreate(PETSC_COMM_SELF, bdName, &bdLabelFull);CHKERRQ(ierr);
1062bc0b47dSJoe Wallwork   ierr = DMPlexMarkBoundaryFaces(dm, 1, bdLabelFull);CHKERRQ(ierr);
1072bc0b47dSJoe Wallwork   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
1082bc0b47dSJoe Wallwork   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
1092bc0b47dSJoe Wallwork   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
1102bc0b47dSJoe Wallwork   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1112bc0b47dSJoe Wallwork     PetscInt *closure = NULL;
1122bc0b47dSJoe Wallwork     PetscInt  closureSize, cl;
1132bc0b47dSJoe Wallwork 
1142bc0b47dSJoe Wallwork     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1152bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1162bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
1172bc0b47dSJoe Wallwork     }
1182bc0b47dSJoe Wallwork     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1192bc0b47dSJoe Wallwork   }
1202bc0b47dSJoe Wallwork   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
1212bc0b47dSJoe Wallwork   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1222bc0b47dSJoe Wallwork     PetscInt *closure = NULL;
1232bc0b47dSJoe Wallwork     PetscInt  closureSize, cl;
1242bc0b47dSJoe Wallwork 
1252bc0b47dSJoe Wallwork     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1262bc0b47dSJoe Wallwork     for (cl = 0; cl < closureSize*2; cl += 2) {
1272bc0b47dSJoe Wallwork       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
1282bc0b47dSJoe Wallwork     }
1292bc0b47dSJoe Wallwork     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1302bc0b47dSJoe Wallwork     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
1312bc0b47dSJoe Wallwork     else         {bdFaceIds[f] = 1;}
1322bc0b47dSJoe Wallwork   }
1332bc0b47dSJoe Wallwork   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
1342bc0b47dSJoe Wallwork   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
1352bc0b47dSJoe Wallwork 
1362bc0b47dSJoe Wallwork   /* Get metric */
1372bc0b47dSJoe Wallwork   ierr = VecViewFromOptions(vertexMetric, NULL, "-adapt_metric_view");CHKERRQ(ierr);
1382bc0b47dSJoe Wallwork   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
13993520041SJoe Wallwork   ierr = DMPlexMetricIsIsotropic(dm, &isotropic);CHKERRQ(ierr);
14093520041SJoe Wallwork   ierr = DMPlexMetricIsUniform(dm, &uniform);CHKERRQ(ierr);
14193520041SJoe Wallwork   Nd = PetscSqr(dim);
14293520041SJoe Wallwork   for (v = 0; v < vEnd-vStart; ++v) {
14393520041SJoe Wallwork     for (i = 0; i < dim; ++i) {
14493520041SJoe Wallwork       for (j = 0; j < dim; ++j) {
14593520041SJoe Wallwork         if (isotropic) {
14693520041SJoe Wallwork           if (i == j) {
14793520041SJoe Wallwork             if (uniform) metric[Nd*v+dim*i+j] = PetscRealPart(met[0]);
14893520041SJoe Wallwork             else metric[Nd*v+dim*i+j] = PetscRealPart(met[v]);
14993520041SJoe Wallwork           } else metric[Nd*v+dim*i+j] = 0.0;
15093520041SJoe Wallwork         } else metric[Nd*v+dim*i+j] = PetscRealPart(met[Nd*v+dim*i+j]);
15193520041SJoe Wallwork       }
15293520041SJoe Wallwork     }
15393520041SJoe Wallwork   }
1542bc0b47dSJoe Wallwork   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
1552bc0b47dSJoe Wallwork 
1562bc0b47dSJoe Wallwork #if 0
1572bc0b47dSJoe Wallwork   /* Destroy overlap mesh */
1582bc0b47dSJoe Wallwork   ierr = DMDestroy(&dm);CHKERRQ(ierr);
1592bc0b47dSJoe Wallwork #endif
1602bc0b47dSJoe Wallwork   /* Send to Pragmatic and remesh */
1612bc0b47dSJoe Wallwork   switch (dim) {
1622bc0b47dSJoe Wallwork   case 2:
1632bc0b47dSJoe Wallwork     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);
1642bc0b47dSJoe Wallwork     break;
1652bc0b47dSJoe Wallwork   case 3:
1662bc0b47dSJoe Wallwork     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);
1672bc0b47dSJoe Wallwork     break;
168*98921bdaSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim);
1692bc0b47dSJoe Wallwork   }
1702bc0b47dSJoe Wallwork   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
1712bc0b47dSJoe Wallwork   pragmatic_set_metric(metric);
1722bc0b47dSJoe Wallwork   pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0);
1732bc0b47dSJoe Wallwork   ierr = PetscFree(l2gv);CHKERRQ(ierr);
1742bc0b47dSJoe Wallwork 
1752bc0b47dSJoe Wallwork   /* Retrieve mesh from Pragmatic and create new plex */
1762bc0b47dSJoe Wallwork   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
1772bc0b47dSJoe Wallwork   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
1782bc0b47dSJoe Wallwork   switch (dim) {
1792bc0b47dSJoe Wallwork   case 2:
1802bc0b47dSJoe Wallwork     numCornersNew = 3;
1812bc0b47dSJoe Wallwork     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
1822bc0b47dSJoe Wallwork     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
1832bc0b47dSJoe Wallwork     break;
1842bc0b47dSJoe Wallwork   case 3:
1852bc0b47dSJoe Wallwork     numCornersNew = 4;
1862bc0b47dSJoe Wallwork     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
1872bc0b47dSJoe Wallwork     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
1882bc0b47dSJoe Wallwork     break;
1892bc0b47dSJoe Wallwork   default:
190*98921bdaSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %D", dim);
1912bc0b47dSJoe Wallwork   }
1922bc0b47dSJoe Wallwork   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];}
1932bc0b47dSJoe Wallwork   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
1942bc0b47dSJoe Wallwork   pragmatic_get_elements(cellsNew);
1952bc0b47dSJoe Wallwork   ierr = DMPlexCreateFromCellListParallelPetsc(comm, dim, numCellsNew, numVerticesNew, PETSC_DECIDE, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, NULL, dmNew);CHKERRQ(ierr);
1962bc0b47dSJoe Wallwork 
1972bc0b47dSJoe Wallwork   /* Rebuild boundary label */
1982bc0b47dSJoe Wallwork   pragmatic_get_boundaryTags(&bdTags);
1992bc0b47dSJoe Wallwork   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
2002bc0b47dSJoe Wallwork   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
2012bc0b47dSJoe Wallwork   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
2022bc0b47dSJoe Wallwork   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
2032bc0b47dSJoe Wallwork   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
2042bc0b47dSJoe Wallwork   for (c = cStart; c < cEnd; ++c) {
2052bc0b47dSJoe Wallwork 
2062bc0b47dSJoe Wallwork     /* Only for simplicial meshes */
2072bc0b47dSJoe Wallwork     coff = (c-cStart)*(dim+1);
2082bc0b47dSJoe Wallwork 
2092bc0b47dSJoe Wallwork     /* d is the local cell number of the vertex opposite to the face we are marking */
2102bc0b47dSJoe Wallwork     for (d = 0; d < dim+1; ++d) {
2112bc0b47dSJoe Wallwork       if (bdTags[coff+d]) {
2122bc0b47dSJoe 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 */
2132bc0b47dSJoe Wallwork         const PetscInt *cone;
2142bc0b47dSJoe Wallwork 
2152bc0b47dSJoe Wallwork         /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
2162bc0b47dSJoe Wallwork         ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr);
2172bc0b47dSJoe Wallwork         ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr);
2182bc0b47dSJoe Wallwork       }
2192bc0b47dSJoe Wallwork     }
2202bc0b47dSJoe Wallwork   }
2212bc0b47dSJoe Wallwork 
2222bc0b47dSJoe Wallwork   /* Clean up */
2232bc0b47dSJoe Wallwork   switch (dim) {
2242bc0b47dSJoe Wallwork   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);
2252bc0b47dSJoe Wallwork   break;
2262bc0b47dSJoe Wallwork   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);
2272bc0b47dSJoe Wallwork   break;
2282bc0b47dSJoe Wallwork   }
2292bc0b47dSJoe Wallwork   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
2302bc0b47dSJoe Wallwork   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
2312bc0b47dSJoe Wallwork   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
2322bc0b47dSJoe Wallwork   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
2332bc0b47dSJoe Wallwork   pragmatic_finalize();
2342bc0b47dSJoe Wallwork   PetscFunctionReturn(0);
2352bc0b47dSJoe Wallwork }
236