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__ 952b40ca0Sbarral #define __FUNCT__ "DMPlexAdapt" 100f7e110dSbarral /*@ 110f7e110dSbarral DMPlexAdapt - Generates a mesh adapted to the specified metric field using the pragmatic library. 120f7e110dSbarral 130f7e110dSbarral Input Parameters: 140f7e110dSbarral + dm - The DM object 150f7e110dSbarral . metric - The metric to which the mesh is adapted, defined vertex-wise. 16*379fadc5Sbarral - bdyLabelName - Label name for boundary tags. These will be preserved in the output mesh. bdyLabelName should be "" (empty string) if there is no such label, and should be different from "boundary". 170f7e110dSbarral 180f7e110dSbarral Output Parameter: 190f7e110dSbarral . dmAdapted - Pointer to the DM object containing the adapted mesh 200f7e110dSbarral 210f7e110dSbarral Level: advanced 220f7e110dSbarral 230f7e110dSbarral .seealso: DMCoarsen(), DMRefine() 240f7e110dSbarral @*/ 256b1afcafSbarral PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdyLabelName[], DM *dmAdapted) 260bbe5d31Sbarral { 270bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC 280bbe5d31Sbarral DM udm, coordDM; 296b1afcafSbarral DMLabel bd, bdTags, bdTagsAdp; 300bbe5d31Sbarral Vec coordinates; 310bbe5d31Sbarral PetscSection coordSection; 320bbe5d31Sbarral IS bdIS; 330f7e110dSbarral double *coordsAdp; 340f7e110dSbarral const PetscScalar *coords; 356b1afcafSbarral PetscReal *x, *y, *z, *xAdp, *yAdp, *zAdp; 360f7e110dSbarral const PetscScalar *metricArray; 370f7e110dSbarral PetscReal *met; 380bbe5d31Sbarral const PetscInt *faces; 396b1afcafSbarral PetscInt *cells, *cellsAdp, *bdFaces, *bdFaceIds; 400f7e110dSbarral PetscInt *boundaryTags; 416b1afcafSbarral PetscInt dim, numCornersAdp, cStart, cEnd, numCells, numCellsAdp, c; 426b1afcafSbarral PetscInt vStart, vEnd, numVertices, numVerticesAdp, v; 436b1afcafSbarral PetscInt numBdFaces, f, maxConeSize, bdSize, coff; 440f7e110dSbarral PetscInt *cell, iVer, cellOff, i, j, idx, facet; 450f7e110dSbarral PetscInt perm[4], isInFacetClosure[4]; // 4 = max number of facets for an element in 2D & 3D. Only for simplicial meshes 46*379fadc5Sbarral PetscBool flag, bdyLabelExt, hasBdyFacet; 470bbe5d31Sbarral #endif 480bbe5d31Sbarral PetscErrorCode ierr; 490bbe5d31Sbarral MPI_Comm comm; 500bbe5d31Sbarral 510bbe5d31Sbarral PetscFunctionBegin; 520bbe5d31Sbarral ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 530bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC 540bbe5d31Sbarral ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 550bbe5d31Sbarral ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 560bbe5d31Sbarral ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr); 570bbe5d31Sbarral ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 580bbe5d31Sbarral ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 590bbe5d31Sbarral ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 600bbe5d31Sbarral ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 610bbe5d31Sbarral ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 620bbe5d31Sbarral numCells = cEnd - cStart; 630bbe5d31Sbarral numVertices = vEnd - vStart; 640bbe5d31Sbarral ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numCells*maxConeSize, &cells, dim*dim*numVertices, &met);CHKERRQ(ierr); 650bbe5d31Sbarral ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 660bbe5d31Sbarral ierr = VecGetArrayRead(metric, &metricArray);CHKERRQ(ierr); 670bbe5d31Sbarral for (v = vStart; v < vEnd; ++v) { 680bbe5d31Sbarral PetscInt off; 690bbe5d31Sbarral 700bbe5d31Sbarral ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 710bbe5d31Sbarral x[v-vStart] = coords[off+0]; 720bbe5d31Sbarral y[v-vStart] = coords[off+1]; 730bbe5d31Sbarral if (dim > 2) z[v-vStart] = coords[off+2]; 740bbe5d31Sbarral ierr = PetscMemcpy(&met[dim*dim*(v-vStart)], &metricArray[dim*off], dim*dim*sizeof(PetscScalar));CHKERRQ(ierr); 750bbe5d31Sbarral } 760bbe5d31Sbarral ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 770bbe5d31Sbarral ierr = VecRestoreArrayRead(metric, &metricArray);CHKERRQ(ierr); 780bbe5d31Sbarral for (c = 0, coff = 0; c < numCells; ++c) { 790bbe5d31Sbarral const PetscInt *cone; 800bbe5d31Sbarral PetscInt coneSize, cl; 810bbe5d31Sbarral 820bbe5d31Sbarral ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 830bbe5d31Sbarral ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 840bbe5d31Sbarral for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 850bbe5d31Sbarral } 860bbe5d31Sbarral switch (dim) { 870bbe5d31Sbarral case 2: 880bbe5d31Sbarral pragmatic_2d_init(&numVertices, &numCells, cells, x, y); 890bbe5d31Sbarral break; 900bbe5d31Sbarral case 3: 910bbe5d31Sbarral pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z); 920bbe5d31Sbarral break; 930bbe5d31Sbarral default: 940bbe5d31Sbarral SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 950bbe5d31Sbarral } 960bbe5d31Sbarral /* Create boundary mesh */ 97*379fadc5Sbarral if (!bdyLabelName) { 98*379fadc5Sbarral SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_NULL, "Argument bdyLabelName cannot be NULL", dim); 99*379fadc5Sbarral } else { 100*379fadc5Sbarral ierr = PetscStrcmp(bdyLabelName, "", &flag);CHKERRQ(ierr); 101*379fadc5Sbarral if (flag) bdyLabelExt = PETSC_FALSE; 102*379fadc5Sbarral else bdyLabelExt = PETSC_TRUE; 103*379fadc5Sbarral } 104*379fadc5Sbarral if (bdyLabelExt) { 1056b1afcafSbarral ierr = PetscStrcmp(bdyLabelName, "boundary", &flag);CHKERRQ(ierr); 1066b1afcafSbarral if (flag) { 1076b1afcafSbarral SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdyLabelName); 1084190e864Sbarral } 1096b1afcafSbarral ierr = DMGetLabel(dm, bdyLabelName, &bdTags);CHKERRQ(ierr); 110*379fadc5Sbarral if (!bdTags) { 111*379fadc5Sbarral SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Label %s does not exist in DM", bdyLabelName); 112*379fadc5Sbarral } 1136b1afcafSbarral } 1140f7e110dSbarral /* TODO To avoid marking bdy facets again if bdyLabelName exists, could/should we loop over bdyLabelName stratas ? */ 1150bbe5d31Sbarral ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr); 1160bbe5d31Sbarral ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr); 1170bbe5d31Sbarral ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr); 1180bbe5d31Sbarral ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr); 1190bbe5d31Sbarral ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr); 1206b1afcafSbarral /* TODO why not assume that we are considering simplicial meshes, in which case bdSize = dim*numBdFaces ? */ 1210bbe5d31Sbarral for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 1220bbe5d31Sbarral PetscInt *closure = NULL; 1230bbe5d31Sbarral PetscInt closureSize, cl; 1240bbe5d31Sbarral 1250bbe5d31Sbarral ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1260bbe5d31Sbarral for (cl = 0; cl < closureSize*2; cl += 2) { 1270bbe5d31Sbarral if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 1280bbe5d31Sbarral } 1290bbe5d31Sbarral ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1300bbe5d31Sbarral } 1310bbe5d31Sbarral ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 1320bbe5d31Sbarral for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 1330bbe5d31Sbarral PetscInt *closure = NULL; 1340bbe5d31Sbarral PetscInt closureSize, cl; 1350bbe5d31Sbarral 1360bbe5d31Sbarral ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1370bbe5d31Sbarral for (cl = 0; cl < closureSize*2; cl += 2) { 1380bbe5d31Sbarral if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 1390bbe5d31Sbarral } 140*379fadc5Sbarral if (bdyLabelExt) { 1416b1afcafSbarral ierr = DMLabelGetValue(bdTags, faces[f], &bdFaceIds[f]);CHKERRQ(ierr); 1426b1afcafSbarral } else { 1430bbe5d31Sbarral bdFaceIds[f] = 1; 1446b1afcafSbarral } 1450bbe5d31Sbarral ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1460bbe5d31Sbarral } 1470bbe5d31Sbarral ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 1480bbe5d31Sbarral ierr = DMLabelDestroy(&bd);CHKERRQ(ierr); 1490bbe5d31Sbarral pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 1500bbe5d31Sbarral pragmatic_set_metric(met); 1510bbe5d31Sbarral pragmatic_adapt(); 1526b1afcafSbarral ierr = PetscFree5(x, y, z, cells, met);CHKERRQ(ierr); 1536b1afcafSbarral ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 1546b1afcafSbarral ierr = PetscFree(coords);CHKERRQ(ierr); 1550bbe5d31Sbarral /* Read out mesh */ 1566b1afcafSbarral pragmatic_get_info(&numVerticesAdp, &numCellsAdp); 1576b1afcafSbarral ierr = PetscMalloc1(numVerticesAdp*dim, &coordsAdp);CHKERRQ(ierr); 1580bbe5d31Sbarral switch (dim) { 1590bbe5d31Sbarral case 2: 1606b1afcafSbarral ierr = PetscMalloc2(numVerticesAdp, &xAdp, numVerticesAdp, &yAdp);CHKERRQ(ierr); 1616b1afcafSbarral zAdp = NULL; 1626b1afcafSbarral pragmatic_get_coords_2d(xAdp, yAdp); 1636b1afcafSbarral numCornersAdp = 3; 1646b1afcafSbarral for (v = 0; v < numVerticesAdp; ++v) {coordsAdp[v*2+0] = xAdp[v]; coordsAdp[v*2+1] = yAdp[v];} 1650bbe5d31Sbarral break; 1660bbe5d31Sbarral case 3: 1676b1afcafSbarral ierr = PetscMalloc3(numVerticesAdp, &xAdp, numVerticesAdp, &yAdp, numVerticesAdp, &zAdp);CHKERRQ(ierr); 1686b1afcafSbarral pragmatic_get_coords_3d(xAdp, yAdp, zAdp); 1696b1afcafSbarral numCornersAdp = 4; 1706b1afcafSbarral for (v = 0; v < numVerticesAdp; ++v) {coordsAdp[v*3+0] = xAdp[v]; coordsAdp[v*3+1] = yAdp[v]; coordsAdp[v*3+2] = zAdp[v];} 1710bbe5d31Sbarral break; 1720bbe5d31Sbarral default: 1736b1afcafSbarral SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 1740bbe5d31Sbarral } 1756b1afcafSbarral ierr = PetscMalloc1(numCellsAdp*(dim+1), &cellsAdp);CHKERRQ(ierr); // only for simplicial meshes 1766b1afcafSbarral pragmatic_get_elements(cellsAdp); 1776b1afcafSbarral ierr = DMPlexCreate(comm, dmAdapted);CHKERRQ(ierr); 1786b1afcafSbarral ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCellsAdp, numVerticesAdp, numCornersAdp, PETSC_TRUE, cellsAdp, dim, coordsAdp, dmAdapted);CHKERRQ(ierr); 179f7caa48eSbarral /* Read out boundary tags */ 180f7caa48eSbarral pragmatic_get_boundaryTags(&boundaryTags); 181*379fadc5Sbarral if (!bdyLabelExt) bdyLabelName = "boundary"; 1826b1afcafSbarral ierr = DMCreateLabel(*dmAdapted, bdyLabelName);CHKERRQ(ierr); 1836b1afcafSbarral ierr = DMGetLabel(*dmAdapted, bdyLabelName, &bdTagsAdp);CHKERRQ(ierr); 1846b1afcafSbarral ierr = DMPlexGetHeightStratum(*dmAdapted, 0, &cStart, &cEnd);CHKERRQ(ierr); 1856b1afcafSbarral ierr = DMPlexGetDepthStratum(*dmAdapted, 0, &vStart, &vEnd);CHKERRQ(ierr); 186f7caa48eSbarral for (c = cStart; c < cEnd; ++c) { 187f7caa48eSbarral PetscInt *cellClosure = NULL; 188f7caa48eSbarral PetscInt cellClosureSize, cl; 189f7caa48eSbarral PetscInt *facetClosure = NULL; 190f7caa48eSbarral PetscInt facetClosureSize, cl2; 191f7caa48eSbarral const PetscInt *facetList; 192f7caa48eSbarral PetscInt facetListSize, f; 193f7caa48eSbarral 1940f7e110dSbarral cellOff = (c-cStart)*(dim+1); // gives the offset corresponding to the cell in pragmatic boundary arrays. Only for simplicial meshes 195*379fadc5Sbarral hasBdyFacet = PETSC_FALSE; 1960f7e110dSbarral for (i = 0; i < dim+1; ++i) { // only for simplicial meshes 1972499d69cSbarral if (boundaryTags[cellOff+i]) { 198*379fadc5Sbarral hasBdyFacet = PETSC_TRUE; 1992499d69cSbarral break; 2002499d69cSbarral } 2012499d69cSbarral } 2020f7e110dSbarral if (!hasBdyFacet) continue; // The cell has no boundary edge/facet => no boundary tagging 2032499d69cSbarral 2040f7e110dSbarral cell = &cellsAdp[cellOff]; // pointing to the current cell in the cellsAdp table 2056b1afcafSbarral ierr = DMPlexGetTransitiveClosure(*dmAdapted, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr); 2066b1afcafSbarral /* first, encode the permutation of the vertices indices betwenn the cell closure and pragmatic ordering */ 207fc167ce7Sbarral j = 0; 208f7caa48eSbarral for (cl = 0; cl < cellClosureSize*2; cl += 2) { 209f7caa48eSbarral if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) continue; 210f7caa48eSbarral iVer = cellClosure[cl] - vStart; 2110f7e110dSbarral for (i = 0; i < dim+1; ++i) { // only for simplicial meshes 212f7caa48eSbarral if (iVer == cell[i]) { 213fc167ce7Sbarral perm[j] = i; // the cl-th element of the closure is the i-th vertex of the cell 214f7caa48eSbarral break; 215f7caa48eSbarral } 216f7caa48eSbarral } 217fc167ce7Sbarral j++; 218f7caa48eSbarral } 2196b1afcafSbarral ierr = DMPlexGetCone(*dmAdapted, c, &facetList);CHKERRQ(ierr); // list of edges/facets of the cell 2206b1afcafSbarral ierr = DMPlexGetConeSize(*dmAdapted, c, &facetListSize);CHKERRQ(ierr); 2210f7e110dSbarral /* then, for each edge/facet of the cell, find the opposite vertex (ie the one in the closure of the cell and not in the closure of the facet/edge) */ 222f7caa48eSbarral for (f = 0; f < facetListSize; ++f){ 223f7caa48eSbarral facet = facetList[f]; 2246b1afcafSbarral ierr = DMPlexGetTransitiveClosure(*dmAdapted, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr); 2256b1afcafSbarral /* loop over the vertices of the cell closure, and check if each vertex belongs to the edge/facet closure */ 2260f7e110dSbarral /* TODO should we use bitmaps to mark vertices instead of a static array ? */ 227f7caa48eSbarral PetscMemzero(isInFacetClosure, sizeof(isInFacetClosure)); 228fc167ce7Sbarral j = 0; 229f7caa48eSbarral for (cl = 0; cl < cellClosureSize*2; cl += 2) { 2300f7e110dSbarral if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) continue; 231f7caa48eSbarral for (cl2 = 0; cl2 < facetClosureSize*2; cl2 += 2){ 232f7caa48eSbarral if ((facetClosure[cl2] < vStart) || (facetClosure[cl2] >= vEnd)) continue; 233f7caa48eSbarral if (cellClosure[cl] == facetClosure[cl2]) { 234fc167ce7Sbarral isInFacetClosure[j] = 1; 235f7caa48eSbarral } 236f7caa48eSbarral } 237fc167ce7Sbarral j++; 238f7caa48eSbarral } 2396b1afcafSbarral /* 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 */ 240fc167ce7Sbarral j = 0; 241f7caa48eSbarral for (cl = 0; cl < cellClosureSize*2; cl += 2) { 242fc167ce7Sbarral if ((cellClosure[cl] < vStart) || (cellClosure[cl] >= vEnd)) continue; 243fc167ce7Sbarral if (!isInFacetClosure[j]) { 244fc167ce7Sbarral idx = cellOff + perm[j]; 2456b1afcafSbarral if (boundaryTags[idx]) { 2466b1afcafSbarral ierr = DMLabelSetValue(bdTagsAdp, facet, boundaryTags[idx]);CHKERRQ(ierr); 2476b1afcafSbarral } 248f7caa48eSbarral break; 249f7caa48eSbarral } 250fc167ce7Sbarral j++; 251f7caa48eSbarral } 2526b1afcafSbarral ierr = DMPlexRestoreTransitiveClosure(*dmAdapted, facet, PETSC_TRUE, &facetClosureSize, &facetClosure);CHKERRQ(ierr); 253f7caa48eSbarral } 2546b1afcafSbarral ierr = DMPlexRestoreTransitiveClosure(*dmAdapted, c, PETSC_TRUE, &cellClosureSize, &cellClosure);CHKERRQ(ierr); 255f7caa48eSbarral } 2560bbe5d31Sbarral pragmatic_finalize(); 2576b1afcafSbarral ierr = PetscFree5(xAdp, yAdp, zAdp, cellsAdp, coordsAdp);CHKERRQ(ierr); 2580bbe5d31Sbarral #endif 2590bbe5d31Sbarral PetscFunctionReturn(0); 2600bbe5d31Sbarral } 261