1*7f296bb3SBarry Smith(ch_unstructured)= 2*7f296bb3SBarry Smith 3*7f296bb3SBarry Smith# DMPlex: Unstructured Grids 4*7f296bb3SBarry Smith 5*7f296bb3SBarry SmithThis chapter introduces the `DMPLEX` subclass of `DM`, which allows 6*7f296bb3SBarry Smiththe user to handle unstructured grids (or meshes) using the generic `DM` interface 7*7f296bb3SBarry Smithfor hierarchy and multi-physics. `DMPLEX` was created to remedy a huge 8*7f296bb3SBarry Smithproblem in all current PDE simulation codes, namely that the 9*7f296bb3SBarry Smithdiscretization was so closely tied to the data layout and solver that 10*7f296bb3SBarry Smithswitching discretizations in the same code was not possible. Not only 11*7f296bb3SBarry Smithdoes this preclude the kind of comparison that is necessary for 12*7f296bb3SBarry Smithscientific investigation, but it makes library (as opposed to monolithic 13*7f296bb3SBarry Smithapplication) development impossible. 14*7f296bb3SBarry Smith 15*7f296bb3SBarry Smith## Representing Unstructured Grids 16*7f296bb3SBarry Smith 17*7f296bb3SBarry SmithThe main advantage of `DMPLEX` in representing topology is that it 18*7f296bb3SBarry Smithtreats all the different pieces of a mesh, e.g. cells, faces, edges, and 19*7f296bb3SBarry Smithvertices, in the same way. This allows the interface to be 20*7f296bb3SBarry Smithsmall and simple, while remaining flexible and general. This also allows 21*7f296bb3SBarry Smith“dimension independent programming”, which means that the same algorithm 22*7f296bb3SBarry Smithcan be used unchanged for meshes of different shapes and dimensions. 23*7f296bb3SBarry Smith 24*7f296bb3SBarry SmithAll pieces of the mesh (vertices, edges, faces, and cells) are treated as *points* (or mesh entities), which are each identified by a 25*7f296bb3SBarry Smith`PetscInt`. A mesh is built by relating points to other points, in 26*7f296bb3SBarry Smithparticular specifying a “covering” relation among the points. For 27*7f296bb3SBarry Smithexample, an edge is defined by being covered by two vertices, and a 28*7f296bb3SBarry Smithtriangle can be defined by being covered by three edges (or even by 29*7f296bb3SBarry Smiththree vertices). This structure is known as a [Hasse Diagram](http://en.wikipedia.org/wiki/Hasse_diagram), which is a 30*7f296bb3SBarry SmithDirected Acyclic Graph (DAG) representing a cell complex using the 31*7f296bb3SBarry Smithcovering relation. The graph edges represent the relation, which also 32*7f296bb3SBarry Smithencodes a partially ordered set (poset). 33*7f296bb3SBarry Smith 34*7f296bb3SBarry SmithFor example, we can encode the doublet mesh as in {numref}`fig_doubletMesh`, 35*7f296bb3SBarry Smith 36*7f296bb3SBarry Smith:::{figure} /images/manual/dmplex_doublet_mesh.svg 37*7f296bb3SBarry Smith:name: fig_doubletMesh 38*7f296bb3SBarry Smith 39*7f296bb3SBarry SmithA 2D doublet mesh, two triangles sharing an edge. 40*7f296bb3SBarry Smith::: 41*7f296bb3SBarry Smith 42*7f296bb3SBarry Smithwhich can also be represented as the DAG in 43*7f296bb3SBarry Smith{numref}`fig_doubletDAG`. 44*7f296bb3SBarry Smith 45*7f296bb3SBarry Smith:::{figure} /images/manual/dmplex_doublet_dag.svg 46*7f296bb3SBarry Smith:name: fig_doubletDAG 47*7f296bb3SBarry Smith 48*7f296bb3SBarry SmithThe Hasse diagram for our 2D doublet mesh, expressed as a DAG. 49*7f296bb3SBarry Smith::: 50*7f296bb3SBarry Smith 51*7f296bb3SBarry SmithTo use the PETSc API, we consecutively number the mesh pieces. The 52*7f296bb3SBarry SmithPETSc convention in 3 dimensions is to number first cells, then 53*7f296bb3SBarry Smithvertices, then faces, and then edges. In 2 dimensions the convention is 54*7f296bb3SBarry Smithto number faces, vertices, and then edges. 55*7f296bb3SBarry SmithIn terms of the labels in 56*7f296bb3SBarry Smith{numref}`fig_doubletMesh`, these numberings are 57*7f296bb3SBarry Smith 58*7f296bb3SBarry Smith$$ 59*7f296bb3SBarry Smithf_0 \mapsto \mathtt{0}, f_1 \mapsto \mathtt{1}, \\ v_0 \mapsto \mathtt{2}, v_1 \mapsto \mathtt{3}, v_2 \mapsto \mathtt{4}, v_3 \mapsto \mathtt{5}, \\ e_0 \mapsto \mathtt{6}, e_1 \mapsto \mathtt{7}, e_2 \mapsto \mathtt{8}, e_3 \mapsto \mathtt{9}, e_4 \mapsto \mathtt{10} 60*7f296bb3SBarry Smith$$ 61*7f296bb3SBarry Smith 62*7f296bb3SBarry SmithFirst, we declare the set of points present in a mesh, 63*7f296bb3SBarry Smith 64*7f296bb3SBarry Smith``` 65*7f296bb3SBarry SmithDMPlexSetChart(dm, 0, 11); 66*7f296bb3SBarry Smith``` 67*7f296bb3SBarry Smith 68*7f296bb3SBarry SmithNote that a *chart* here corresponds to a semi-closed interval (e.g 69*7f296bb3SBarry Smith$[0,11) = \{0,1,\ldots,10\}$) specifying the range of indices we’d 70*7f296bb3SBarry Smithlike to use to define points on the current rank. We then define the 71*7f296bb3SBarry Smithcovering relation, which we call the *cone*, which are also the in-edges 72*7f296bb3SBarry Smithin the DAG. In order to preallocate correctly, we first provide sizes, 73*7f296bb3SBarry Smith 74*7f296bb3SBarry Smith``` 75*7f296bb3SBarry Smith/* DMPlexSetConeSize(dm, point, number of points that cover the point); */ 76*7f296bb3SBarry SmithDMPlexSetConeSize(dm, 0, 3); 77*7f296bb3SBarry SmithDMPlexSetConeSize(dm, 1, 3); 78*7f296bb3SBarry SmithDMPlexSetConeSize(dm, 6, 2); 79*7f296bb3SBarry SmithDMPlexSetConeSize(dm, 7, 2); 80*7f296bb3SBarry SmithDMPlexSetConeSize(dm, 8, 2); 81*7f296bb3SBarry SmithDMPlexSetConeSize(dm, 9, 2); 82*7f296bb3SBarry SmithDMPlexSetConeSize(dm, 10, 2); 83*7f296bb3SBarry SmithDMSetUp(dm); 84*7f296bb3SBarry Smith``` 85*7f296bb3SBarry Smith 86*7f296bb3SBarry Smithand then point values (recall each point is an integer that represents a single geometric entity, a cell, face, edge, or vertex), 87*7f296bb3SBarry Smith 88*7f296bb3SBarry Smith``` 89*7f296bb3SBarry Smith/* DMPlexSetCone(dm, point, [points that cover the point]); */ 90*7f296bb3SBarry SmithDMPlexSetCone(dm, 0, [6, 7, 8]); 91*7f296bb3SBarry SmithDMPlexSetCone(dm, 1, [7, 9, 10]); 92*7f296bb3SBarry SmithDMPlexSetCone(dm, 6, [2, 3]); 93*7f296bb3SBarry SmithDMPlexSetCone(dm, 7, [3, 4]); 94*7f296bb3SBarry SmithDMPlexSetCone(dm, 8, [4, 2]); 95*7f296bb3SBarry SmithDMPlexSetCone(dm, 9, [4, 5]); 96*7f296bb3SBarry SmithDMPlexSetCone(dm, 10, [5, 3]); 97*7f296bb3SBarry Smith``` 98*7f296bb3SBarry Smith 99*7f296bb3SBarry SmithThere is also an API for providing the dual relation, using 100*7f296bb3SBarry Smith`DMPlexSetSupportSize()` and `DMPlexSetSupport()`, but this can be 101*7f296bb3SBarry Smithcalculated automatically using the provided `DMPlexSetConeSize()` and `DMPlexSetCone()` information and then calling 102*7f296bb3SBarry Smith 103*7f296bb3SBarry Smith``` 104*7f296bb3SBarry SmithDMPlexSymmetrize(dm); 105*7f296bb3SBarry Smith``` 106*7f296bb3SBarry Smith 107*7f296bb3SBarry SmithThe "symmetrization" is in the sense of the DAG. Each point knows its covering (cone) and each point knows what it covers (support). Note that when using automatic symmetrization, cones will be ordered but supports will not. The user can enforce an ordering on supports by rewriting them after symmetrization using `DMPlexSetSupport()`. 108*7f296bb3SBarry Smith 109*7f296bb3SBarry SmithIn order to support efficient queries, we construct fast 110*7f296bb3SBarry Smithsearch structures and indices for the different types of points using 111*7f296bb3SBarry Smith 112*7f296bb3SBarry Smith``` 113*7f296bb3SBarry SmithDMPlexStratify(dm); 114*7f296bb3SBarry Smith``` 115*7f296bb3SBarry Smith 116*7f296bb3SBarry Smith## Grid Point Orientations 117*7f296bb3SBarry Smith 118*7f296bb3SBarry SmithTODO: fill this out with regard to orientations. 119*7f296bb3SBarry Smith 120*7f296bb3SBarry SmithShould probably reference Hapla and summarize what's there. 121*7f296bb3SBarry Smith 122*7f296bb3SBarry Smith## Dealing with Periodicity 123*7f296bb3SBarry Smith 124*7f296bb3SBarry SmithPlex allows you to represent periodic domains is two ways. Using the default scheme, periodic topology can be represented directly. This ensures that all topological queries can be satisfied, but then care must be taken in representing functions over the mesh, such as the coordinates. The second method is to use a non-periodic topology, but connect certain mesh points using the local-to-global map for that DM. This allows a more general set of mappings to be implemented, such as partial twists, but topological queries on the periodic boundary cease to function. 125*7f296bb3SBarry Smith 126*7f296bb3SBarry SmithFor the default scheme, a call to `DMLocalizeCoordinates()` (which usually happens automatically on mesh creation) creates a second, discontinuous coordinate field. These values can be accessed using `DMGetCellCoordinates()` and `DMGetCellCoordinatesLocal()`. Plex provides a convenience method, `DMPlexGetCellCoordinates()`, that extracts cell coordinates correctly, depending on the periodicity of the mesh. An example of its use is shown below: 127*7f296bb3SBarry Smith 128*7f296bb3SBarry Smith``` 129*7f296bb3SBarry Smithconst PetscScalar *array; 130*7f296bb3SBarry SmithPetscScalar *coords = NULL; 131*7f296bb3SBarry SmithPetscInt numCoords; 132*7f296bb3SBarry SmithPetscBool isDG; 133*7f296bb3SBarry Smith 134*7f296bb3SBarry SmithPetscCall(DMPlexGetCellCoordinates(dm, cell, &isDG, &numCoords, &array, &coords)); 135*7f296bb3SBarry Smithfor (PetscInt cc = 0; cc < numCoords/dim; ++cc) { 136*7f296bb3SBarry Smith if (cc > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " -- ")); 137*7f296bb3SBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF, "(")); 138*7f296bb3SBarry Smith for (PetscInt d = 0; d < dim; ++d) { 139*7f296bb3SBarry Smith if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 140*7f296bb3SBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(coords[cc * dim + d]))); 141*7f296bb3SBarry Smith } 142*7f296bb3SBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF, ")")); 143*7f296bb3SBarry Smith} 144*7f296bb3SBarry SmithPetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 145*7f296bb3SBarry SmithPetscCall(DMPlexRestoreCellCoordinates(dm, cell, &isDG, &numCoords, &array, &coords)); 146*7f296bb3SBarry Smith``` 147*7f296bb3SBarry Smith 148*7f296bb3SBarry Smith## Connecting Grids to Data Using PetscSection: 149*7f296bb3SBarry Smith 150*7f296bb3SBarry SmithA `PetscSection` is used to describe the connection between the grid and data associated with the grid. 151*7f296bb3SBarry SmithSpecifically, it assigns a number of dofs to each mesh entity on the DAG. 152*7f296bb3SBarry SmithSee {any}`ch_petscsection` for more details. 153*7f296bb3SBarry SmithUsing the mesh from {numref}`fig_doubletMesh`, we provide an example of creating a `PetscSection` for a single field. 154*7f296bb3SBarry SmithWe can lay out data for a continuous Galerkin $P_3$ finite element method, 155*7f296bb3SBarry Smith 156*7f296bb3SBarry Smith``` 157*7f296bb3SBarry SmithPetscInt pStart, pEnd, cStart, cEnd, c, vStart, vEnd, v, eStart, eEnd, e; 158*7f296bb3SBarry Smith 159*7f296bb3SBarry SmithDMPlexGetChart(dm, &pStart, &pEnd); 160*7f296bb3SBarry SmithDMPlexGetHeightStratum(dm, 0, &cStart, &cEnd); // cells 161*7f296bb3SBarry SmithDMPlexGetHeightStratum(dm, 1, &eStart, &eEnd); // edges 162*7f296bb3SBarry SmithDMPlexGetHeightStratum(dm, 2, &vStart, &vEnd); // vertices, equivalent to DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd); 163*7f296bb3SBarry SmithPetscSectionSetChart(s, pStart, pEnd); 164*7f296bb3SBarry Smithfor(c = cStart; c < cEnd; ++c) 165*7f296bb3SBarry Smith PetscSectionSetDof(s, c, 1); 166*7f296bb3SBarry Smithfor(v = vStart; v < vEnd; ++v) 167*7f296bb3SBarry Smith PetscSectionSetDof(s, v, 1); 168*7f296bb3SBarry Smithfor(e = eStart; e < eEnd; ++e) 169*7f296bb3SBarry Smith PetscSectionSetDof(s, e, 2); // two dof on each edge 170*7f296bb3SBarry SmithPetscSectionSetUp(s); 171*7f296bb3SBarry Smith``` 172*7f296bb3SBarry Smith 173*7f296bb3SBarry Smith`DMPlexGetHeightStratum()` returns all the points of the requested height 174*7f296bb3SBarry Smithin the DAG. Since this problem is in two dimensions the edges are at 175*7f296bb3SBarry Smithheight 1 and the vertices at height 2 (the cells are always at height 176*7f296bb3SBarry Smith0). One can also use `DMPlexGetDepthStratum()` to use the depth in the 177*7f296bb3SBarry SmithDAG to select the points. `DMPlexGetDepth(dm,&depth)` returns the depth 178*7f296bb3SBarry Smithof the DAG, hence `DMPlexGetDepthStratum(dm,depth-1-h,)` returns the 179*7f296bb3SBarry Smithsame values as `DMPlexGetHeightStratum(dm,h,)`. 180*7f296bb3SBarry Smith 181*7f296bb3SBarry SmithFor $P_3$ elements there is one degree of freedom at each vertex, 2 along 182*7f296bb3SBarry Smitheach edge (resulting in a total of 4 degrees of freedom along each edge 183*7f296bb3SBarry Smithincluding the vertices, thus being able to reproduce a cubic function) 184*7f296bb3SBarry Smithand 1 degree of freedom within the cell (the bubble function which is 185*7f296bb3SBarry Smithzero along all edges). 186*7f296bb3SBarry Smith 187*7f296bb3SBarry SmithNow a PETSc local vector can be created manually using this layout, 188*7f296bb3SBarry Smith 189*7f296bb3SBarry Smith``` 190*7f296bb3SBarry SmithPetscSectionGetStorageSize(s, &n); 191*7f296bb3SBarry SmithVecSetSizes(localVec, n, PETSC_DETERMINE); 192*7f296bb3SBarry SmithVecSetFromOptions(localVec); 193*7f296bb3SBarry Smith``` 194*7f296bb3SBarry Smith 195*7f296bb3SBarry SmithWhen working with `DMPLEX` and `PetscFE` (see below) one can simply get the sections (and related vectors) with 196*7f296bb3SBarry Smith 197*7f296bb3SBarry Smith``` 198*7f296bb3SBarry SmithDMSetLocalSection(dm, s); 199*7f296bb3SBarry SmithDMGetLocalVector(dm, &localVec); 200*7f296bb3SBarry SmithDMGetGlobalVector(dm, &globalVec); 201*7f296bb3SBarry Smith``` 202*7f296bb3SBarry Smith 203*7f296bb3SBarry Smith### DMPlex-specific PetscSection Features: 204*7f296bb3SBarry Smith 205*7f296bb3SBarry SmithThe following features are built into `PetscSection`. 206*7f296bb3SBarry SmithHowever, their usage and purpose is best understood through `DMPLEX`. 207*7f296bb3SBarry Smith 208*7f296bb3SBarry Smith#### Closure: 209*7f296bb3SBarry Smith 210*7f296bb3SBarry SmithClosure information can be attached to a `PetscSection` to allow for more efficient closure information queries. 211*7f296bb3SBarry SmithThis information can either be set directly with `DMPlexCreateClosureIndex()` or generated automatically for a `DMPLEX` via `DMPlexCreateClosureIndex()`. 212*7f296bb3SBarry Smith 213*7f296bb3SBarry Smith#### Symmetries: Accessing data from different orientations 214*7f296bb3SBarry Smith 215*7f296bb3SBarry SmithWhile mesh point orientation information specifies how one mesh point is oriented with respect to another, it does not describe how the dofs associated with that mesh point should be permuted for that orientation. 216*7f296bb3SBarry SmithThis information is supplied via a `PetscSectionSym` object that is attached to the `PetscSection`. 217*7f296bb3SBarry SmithGenerally the setup and usage of this information is handled automatically by PETSc during setup of a Plex using `PetscFE`. 218*7f296bb3SBarry Smith 219*7f296bb3SBarry Smith#### Closure Permutation: 220*7f296bb3SBarry Smith 221*7f296bb3SBarry SmithA permutation of the dof closure of a k-cell may be specified. 222*7f296bb3SBarry SmithThis allows data to be returned in an order that might be more efficiently processed than the default (breadth-first search) ordering. 223*7f296bb3SBarry SmithFor example, for tensor cells such as quadrilaterals, closure data can be permuted to lexicographic order (i.e. a tensor-product ordering). 224*7f296bb3SBarry SmithThis is most commonly done via `DMPlexSetClosurePermutationTensor()`. 225*7f296bb3SBarry SmithCustom permutations can be set using `PetscSectionSetClosurePermutation()`. 226*7f296bb3SBarry Smith 227*7f296bb3SBarry Smith## Data Layout using DMPLEX and PetscFE 228*7f296bb3SBarry Smith 229*7f296bb3SBarry SmithA `DM` can automatically create the local section if given a description of the discretization, for example using a `PetscFE` object. We demonstrate this by creating 230*7f296bb3SBarry Smitha `PetscFE` that can be configured from the command line. It is a single, scalar field, and is added to the `DM` using `DMSetField()`. 231*7f296bb3SBarry SmithWhen a local or global vector is requested, the `DM` builds the local and global sections automatically. 232*7f296bb3SBarry Smith 233*7f296bb3SBarry Smith``` 234*7f296bb3SBarry SmithDMPlexIsSimplex(dm, &simplex); 235*7f296bb3SBarry SmithPetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe); 236*7f296bb3SBarry SmithDMSetField(dm, 0, NULL, (PetscObject) fe); 237*7f296bb3SBarry SmithDMCreateDS(dm); 238*7f296bb3SBarry Smith``` 239*7f296bb3SBarry Smith 240*7f296bb3SBarry SmithHere the call to `DMSetField()` declares the discretization will have one field with the integer label 0 that has one degree of freedom at each point on the `DMPlex`. 241*7f296bb3SBarry SmithTo get the $P_3$ section above, we can either give the option `-petscspace_degree 3`, or call `PetscFECreateLagrange()` and set the degree directly. 242*7f296bb3SBarry Smith 243*7f296bb3SBarry Smith## Partitioning and Ordering 244*7f296bb3SBarry Smith 245*7f296bb3SBarry SmithIn the same way as `MatPartitioning` or 246*7f296bb3SBarry Smith`MatGetOrdering()`, give the results of a partitioning or ordering of a graph defined by a sparse matrix, 247*7f296bb3SBarry Smith`PetscPartitionerDMPlexPartition` or `DMPlexPermute` are encoded in 248*7f296bb3SBarry Smithan `IS`. However, the graph is not the adjacency graph of the matrix 249*7f296bb3SBarry Smithbut the mesh itself. Once the mesh is partitioned and 250*7f296bb3SBarry Smithreordered, the data layout from a `PetscSection` can be used to 251*7f296bb3SBarry Smithautomatically derive a problem partitioning/ordering. 252*7f296bb3SBarry Smith 253*7f296bb3SBarry Smith### Influence of Variables on One Another 254*7f296bb3SBarry Smith 255*7f296bb3SBarry SmithThe Jacobian of a problem represents the influence of some 256*7f296bb3SBarry Smithvariable on other variables in the problem. Very often, this influence 257*7f296bb3SBarry Smithpattern is determined jointly by the computational mesh and 258*7f296bb3SBarry Smithdiscretization. `DMCreateMatrix()` must compute this pattern when it 259*7f296bb3SBarry Smithautomatically creates the properly preallocated Jacobian matrix. In 260*7f296bb3SBarry Smith`DMDA` the influence pattern, or what we will call variable 261*7f296bb3SBarry Smith*adjacency*, depends only on the stencil since the topology is Cartesian 262*7f296bb3SBarry Smithand the discretization is implicitly finite difference. 263*7f296bb3SBarry Smith 264*7f296bb3SBarry SmithIn `DMPLEX`, 265*7f296bb3SBarry Smithwe allow the user to specify the adjacency topologically, while 266*7f296bb3SBarry Smithmaintaining good defaults. The pattern is controlled by two flags. The first flag, `useCone`, 267*7f296bb3SBarry Smithindicates whether variables couple first to their boundary [^boundary-footnote] 268*7f296bb3SBarry Smithand then to 269*7f296bb3SBarry Smithneighboring entities, or the reverse. For example, in finite elements, 270*7f296bb3SBarry Smiththe variables couple to the set of neighboring cells containing the mesh 271*7f296bb3SBarry Smithpoint, and we set the flag to `useCone = PETSC_FALSE`. By contrast, 272*7f296bb3SBarry Smithin finite volumes, cell variables first couple to the cell boundary, and 273*7f296bb3SBarry Smiththen to the neighbors, so we set the flag to `useCone = PETSC_TRUE`. 274*7f296bb3SBarry SmithThe second flag, `useClosure`, indicates whether we consider the 275*7f296bb3SBarry Smithtransitive closure of the neighbor relation above, or just a single 276*7f296bb3SBarry Smithlevel. For example, in finite elements, the entire boundary of any cell 277*7f296bb3SBarry Smithcouples to the interior, and we set the flag to 278*7f296bb3SBarry Smith`useClosure = PETSC_TRUE`. By contrast, in most finite volume methods, 279*7f296bb3SBarry Smithcells couple only across faces, and not through vertices, so we set the 280*7f296bb3SBarry Smithflag to `useClosure = PETSC_FALSE`. However, the power of this method 281*7f296bb3SBarry Smithis its flexibility. If we wanted a finite volume method that coupled all 282*7f296bb3SBarry Smithcells around a vertex, we could easily prescribe that by changing to 283*7f296bb3SBarry Smith`useClosure = PETSC_TRUE`. 284*7f296bb3SBarry Smith 285*7f296bb3SBarry Smith## Evaluating Residuals 286*7f296bb3SBarry Smith 287*7f296bb3SBarry SmithThe evaluation of a residual or Jacobian, for most discretizations has 288*7f296bb3SBarry Smiththe following general form: 289*7f296bb3SBarry Smith 290*7f296bb3SBarry Smith- Traverse the mesh, picking out pieces (which in general overlap), 291*7f296bb3SBarry Smith- Extract some values from the current solution vector, associated with this 292*7f296bb3SBarry Smith piece, 293*7f296bb3SBarry Smith- Calculate some values for the piece, and 294*7f296bb3SBarry Smith- Insert these values into the residual vector 295*7f296bb3SBarry Smith 296*7f296bb3SBarry SmithDMPlex separates these different concerns by passing sets of points from mesh traversal routines to data 297*7f296bb3SBarry Smithextraction routines and back. In this way, the `PetscSection` which 298*7f296bb3SBarry Smithstructures the data inside a `Vec` does not need to know anything 299*7f296bb3SBarry Smithabout the mesh inside a `DMPLEX`. 300*7f296bb3SBarry Smith 301*7f296bb3SBarry SmithThe most common mesh traversal is the transitive closure of a point, 302*7f296bb3SBarry Smithwhich is exactly the transitive closure of a point in the DAG using the 303*7f296bb3SBarry Smithcovering relation. In other words, the transitive closure consists of 304*7f296bb3SBarry Smithall points that cover the given point (generally a cell) plus all points 305*7f296bb3SBarry Smiththat cover those points, etc. So in 2d the transitive closure for a cell 306*7f296bb3SBarry Smithconsists of edges and vertices while in 3d it consists of faces, edges, 307*7f296bb3SBarry Smithand vertices. Note that this closure can be calculated orienting the 308*7f296bb3SBarry Smitharrows in either direction. For example, in a finite element 309*7f296bb3SBarry Smithcalculation, we calculate an integral over each element, and then sum up 310*7f296bb3SBarry Smiththe contributions to the basis function coefficients. The closure of the 311*7f296bb3SBarry Smithelement can be expressed discretely as the transitive closure of the 312*7f296bb3SBarry Smithelement point in the mesh DAG, where each point also has an orientation. 313*7f296bb3SBarry SmithThen we can retrieve the data using `PetscSection` methods, 314*7f296bb3SBarry Smith 315*7f296bb3SBarry Smith``` 316*7f296bb3SBarry SmithPetscScalar *a; 317*7f296bb3SBarry SmithPetscInt numPoints, *points = NULL, p; 318*7f296bb3SBarry Smith 319*7f296bb3SBarry SmithVecGetArrayRead(u,&a); 320*7f296bb3SBarry SmithDMPlexGetTransitiveClosure(dm,cell,PETSC_TRUE,&numPoints,&points); 321*7f296bb3SBarry Smithfor (p = 0; p <= numPoints*2; p += 2) { 322*7f296bb3SBarry Smith PetscInt dof, off, d; 323*7f296bb3SBarry Smith 324*7f296bb3SBarry Smith PetscSectionGetDof(section, points[p], &dof); 325*7f296bb3SBarry Smith PetscSectionGetOffset(section, points[p], &off); 326*7f296bb3SBarry Smith for (d = 0; d <= dof; ++d) { 327*7f296bb3SBarry Smith myfunc(a[off+d]); 328*7f296bb3SBarry Smith } 329*7f296bb3SBarry Smith} 330*7f296bb3SBarry SmithDMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &numPoints, &points); 331*7f296bb3SBarry SmithVecRestoreArrayRead(u, &a); 332*7f296bb3SBarry Smith``` 333*7f296bb3SBarry Smith 334*7f296bb3SBarry SmithThis operation is so common that we have built a convenience method 335*7f296bb3SBarry Smitharound it which returns the values in a contiguous array, correctly 336*7f296bb3SBarry Smithtaking into account the orientations of various mesh points: 337*7f296bb3SBarry Smith 338*7f296bb3SBarry Smith``` 339*7f296bb3SBarry Smithconst PetscScalar *values; 340*7f296bb3SBarry SmithPetscInt csize; 341*7f296bb3SBarry Smith 342*7f296bb3SBarry SmithDMPlexVecGetClosure(dm, section, u, cell, &csize, &values); 343*7f296bb3SBarry Smith// Do integral in quadrature loop putting the result into r[] 344*7f296bb3SBarry SmithDMPlexVecRestoreClosure(dm, section, u, cell, &csize, &values); 345*7f296bb3SBarry SmithDMPlexVecSetClosure(dm, section, residual, cell, &r, ADD_VALUES); 346*7f296bb3SBarry Smith``` 347*7f296bb3SBarry Smith 348*7f296bb3SBarry SmithA simple example of this kind of calculation is in 349*7f296bb3SBarry Smith`DMPlexComputeL2Diff_Plex()` (<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/plex/plexfem.c.html#DMComputeL2Diff_Plex">source</a>). 350*7f296bb3SBarry SmithNote that there is no restriction on the type of cell or dimension of 351*7f296bb3SBarry Smiththe mesh in the code above, so it will work for polyhedral cells, hybrid 352*7f296bb3SBarry Smithmeshes, and meshes of any dimension, without change. We can also reverse 353*7f296bb3SBarry Smiththe covering relation, so that the code works for finite volume methods 354*7f296bb3SBarry Smithwhere we want the data from neighboring cells for each face: 355*7f296bb3SBarry Smith 356*7f296bb3SBarry Smith``` 357*7f296bb3SBarry SmithPetscScalar *a; 358*7f296bb3SBarry SmithPetscInt points[2*2], numPoints, p, dofA, offA, dofB, offB; 359*7f296bb3SBarry Smith 360*7f296bb3SBarry SmithVecGetArray(u, &a); 361*7f296bb3SBarry SmithDMPlexGetTransitiveClosure(dm, cell, PETSC_FALSE, &numPoints, &points); 362*7f296bb3SBarry Smithassert(numPoints == 2); 363*7f296bb3SBarry SmithPetscSectionGetDof(section, points[0*2], &dofA); 364*7f296bb3SBarry SmithPetscSectionGetDof(section, points[1*2], &dofB); 365*7f296bb3SBarry Smithassert(dofA == dofB); 366*7f296bb3SBarry SmithPetscSectionGetOffset(section, points[0*2], &offA); 367*7f296bb3SBarry SmithPetscSectionGetOffset(section, points[1*2], &offB); 368*7f296bb3SBarry Smithmyfunc(a[offA], a[offB]); 369*7f296bb3SBarry SmithVecRestoreArray(u, &a); 370*7f296bb3SBarry Smith``` 371*7f296bb3SBarry Smith 372*7f296bb3SBarry SmithThis kind of calculation is used in 373*7f296bb3SBarry Smith<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex11.c.html">TS Tutorial ex11</a>. 374*7f296bb3SBarry Smith 375*7f296bb3SBarry Smith## Saving and Loading DMPlex Data with HDF5 376*7f296bb3SBarry Smith 377*7f296bb3SBarry SmithPETSc allows users to save/load `DMPLEX`s representing meshes, 378*7f296bb3SBarry Smith`PetscSection`s representing data layouts on the meshes, and 379*7f296bb3SBarry Smith`Vec`s defined on the data layouts to/from an HDF5 file in 380*7f296bb3SBarry Smithparallel, where one can use different number of processes for saving 381*7f296bb3SBarry Smithand for loading. 382*7f296bb3SBarry Smith 383*7f296bb3SBarry Smith### Saving 384*7f296bb3SBarry Smith 385*7f296bb3SBarry SmithThe simplest way to save `DM` data is to use options for configuration. 386*7f296bb3SBarry SmithThis requires only the code 387*7f296bb3SBarry Smith 388*7f296bb3SBarry Smith``` 389*7f296bb3SBarry SmithDMViewFromOptions(dm, NULL, "-dm_view"); 390*7f296bb3SBarry SmithVecViewFromOptions(vec, NULL, "-vec_view") 391*7f296bb3SBarry Smith``` 392*7f296bb3SBarry Smith 393*7f296bb3SBarry Smithalong with the command line options 394*7f296bb3SBarry Smith 395*7f296bb3SBarry Smith```console 396*7f296bb3SBarry Smith$ ./myprog -dm_view hdf5:myprog.h5 -vec_view hdf5:myprog.h5::append 397*7f296bb3SBarry Smith``` 398*7f296bb3SBarry Smith 399*7f296bb3SBarry SmithOptions prefixes can be used to separately control the saving and loading of various fields. 400*7f296bb3SBarry SmithHowever, the user can have finer grained control by explicitly creating the PETSc objects involved. 401*7f296bb3SBarry SmithTo save data to "example.h5" file, we can first create a `PetscViewer` of type `PETSCVIEWERHDF5` in `FILE_MODE_WRITE` mode as: 402*7f296bb3SBarry Smith 403*7f296bb3SBarry Smith``` 404*7f296bb3SBarry SmithPetscViewer viewer; 405*7f296bb3SBarry Smith 406*7f296bb3SBarry SmithPetscViewerHDF5Open(PETSC_COMM_WORLD, "example.h5", FILE_MODE_WRITE, &viewer); 407*7f296bb3SBarry Smith``` 408*7f296bb3SBarry Smith 409*7f296bb3SBarry SmithAs `dm` is a `DMPLEX` object representing a mesh, we first give it a *mesh name*, "plexA", and save it as: 410*7f296bb3SBarry Smith 411*7f296bb3SBarry Smith``` 412*7f296bb3SBarry SmithPetscObjectSetName((PetscObject)dm, "plexA"); 413*7f296bb3SBarry SmithPetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_PETSC); 414*7f296bb3SBarry SmithDMView(dm, viewer); 415*7f296bb3SBarry SmithPetscViewerPopFormat(viewer); 416*7f296bb3SBarry Smith``` 417*7f296bb3SBarry Smith 418*7f296bb3SBarry SmithThe `DMView()` call is shorthand for the following sequence 419*7f296bb3SBarry Smith 420*7f296bb3SBarry Smith``` 421*7f296bb3SBarry SmithDMPlexTopologyView(dm, viewer); 422*7f296bb3SBarry SmithDMPlexCoordinatesView(dm, viewer); 423*7f296bb3SBarry SmithDMPlexLabelsView(dm, viewer); 424*7f296bb3SBarry Smith``` 425*7f296bb3SBarry Smith 426*7f296bb3SBarry SmithIf the *mesh name* is not explicitly set, the default name is used. 427*7f296bb3SBarry SmithIn the above `PETSC_VIEWER_HDF5_PETSC` format was used to save the entire representation of the mesh. 428*7f296bb3SBarry SmithThis format also saves global point numbers attached to the mesh points. 429*7f296bb3SBarry SmithIn this example the set of all global point numbers is $X = [0, 11)$. 430*7f296bb3SBarry Smith 431*7f296bb3SBarry SmithThe data layout, `s`, needs to be wrapped in a `DM` object for it to be saved. 432*7f296bb3SBarry SmithHere, we create the wrapping `DM`, `sdm`, with `DMClone()`, give it a *dm name*, "dmA", attach `s` to `sdm`, and save it as: 433*7f296bb3SBarry Smith 434*7f296bb3SBarry Smith``` 435*7f296bb3SBarry SmithDMClone(dm, &sdm); 436*7f296bb3SBarry SmithPetscObjectSetName((PetscObject)sdm, "dmA"); 437*7f296bb3SBarry SmithDMSetLocalSection(sdm, s); 438*7f296bb3SBarry SmithDMPlexSectionView(dm, viewer, sdm); 439*7f296bb3SBarry Smith``` 440*7f296bb3SBarry Smith 441*7f296bb3SBarry SmithIf the *dm name* is not explicitly set, the default name is to be used. 442*7f296bb3SBarry SmithIn the above, instead of using `DMClone()`, one could also create a new `DMSHELL` object to attach `s` to. 443*7f296bb3SBarry SmithThe first argument of `DMPlexSectionView()` is a `DMPLEX` object that represents the mesh, and the third argument is a `DM` object that carries the data layout that we would like to save. 444*7f296bb3SBarry SmithThey are, in general, two different objects, and the former carries a *mesh name*, while the latter carries a *dm name*. 445*7f296bb3SBarry SmithThese names are used to construct a group structure in the HDF5 file. 446*7f296bb3SBarry SmithNote that the data layout points are associated with the mesh points, so each of them can also be tagged with a global point number in $X$; `DMPlexSectionView()` saves these tags along with the data layout itself, so that, when the mesh and the data layout are loaded separately later, one can associate the points in the former with those in the latter by comparing their global point numbers. 447*7f296bb3SBarry Smith 448*7f296bb3SBarry SmithWe now create a local vector associated with `sdm`, e.g., as: 449*7f296bb3SBarry Smith 450*7f296bb3SBarry Smith``` 451*7f296bb3SBarry SmithVec vec; 452*7f296bb3SBarry Smith 453*7f296bb3SBarry SmithDMGetLocalVector(sdm, &vec); 454*7f296bb3SBarry Smith``` 455*7f296bb3SBarry Smith 456*7f296bb3SBarry SmithAfter setting values of `vec`, we name it "vecA" and save it as: 457*7f296bb3SBarry Smith 458*7f296bb3SBarry Smith``` 459*7f296bb3SBarry SmithPetscObjectSetName((PetscObject)vec, "vecA"); 460*7f296bb3SBarry SmithDMPlexLocalVectorView(dm, viewer, sdm, vec); 461*7f296bb3SBarry Smith``` 462*7f296bb3SBarry Smith 463*7f296bb3SBarry SmithA global vector can be saved in the exact same way with trivial changes. 464*7f296bb3SBarry Smith 465*7f296bb3SBarry SmithAfter saving, we destroy the `PetscViewer` with: 466*7f296bb3SBarry Smith 467*7f296bb3SBarry Smith``` 468*7f296bb3SBarry SmithPetscViewerDestroy(&viewer); 469*7f296bb3SBarry Smith``` 470*7f296bb3SBarry Smith 471*7f296bb3SBarry SmithThe output file "example.h5" now looks like the following: 472*7f296bb3SBarry Smith 473*7f296bb3SBarry Smith``` 474*7f296bb3SBarry Smith$ h5dump --contents example.h5 475*7f296bb3SBarry SmithHDF5 "example.h5" { 476*7f296bb3SBarry SmithFILE_CONTENTS { 477*7f296bb3SBarry Smith group / 478*7f296bb3SBarry Smith group /topologies 479*7f296bb3SBarry Smith group /topologies/plexA 480*7f296bb3SBarry Smith group /topologies/plexA/dms 481*7f296bb3SBarry Smith group /topologies/plexA/dms/dmA 482*7f296bb3SBarry Smith dataset /topologies/plexA/dms/dmA/order 483*7f296bb3SBarry Smith group /topologies/plexA/dms/dmA/section 484*7f296bb3SBarry Smith dataset /topologies/plexA/dms/dmA/section/atlasDof 485*7f296bb3SBarry Smith dataset /topologies/plexA/dms/dmA/section/atlasOff 486*7f296bb3SBarry Smith group /topologies/plexA/dms/dmA/vecs 487*7f296bb3SBarry Smith group /topologies/plexA/dms/dmA/vecs/vecA 488*7f296bb3SBarry Smith dataset /topologies/plexA/dms/dmA/vecs/vecA/vecA 489*7f296bb3SBarry Smith group /topologies/plexA/labels 490*7f296bb3SBarry Smith group /topologies/plexA/topology 491*7f296bb3SBarry Smith dataset /topologies/plexA/topology/cells 492*7f296bb3SBarry Smith dataset /topologies/plexA/topology/cones 493*7f296bb3SBarry Smith dataset /topologies/plexA/topology/order 494*7f296bb3SBarry Smith dataset /topologies/plexA/topology/orientation 495*7f296bb3SBarry Smith } 496*7f296bb3SBarry Smith} 497*7f296bb3SBarry Smith``` 498*7f296bb3SBarry Smith 499*7f296bb3SBarry Smith### Saving in the new parallel HDF5 format 500*7f296bb3SBarry Smith 501*7f296bb3SBarry SmithSince PETSc 3.19, we offer a new format which supports parallel loading. 502*7f296bb3SBarry SmithTo write in this format, you currently need to specify it explicitly using the option 503*7f296bb3SBarry Smith 504*7f296bb3SBarry Smith``` 505*7f296bb3SBarry Smith-dm_plex_view_hdf5_storage_version 3.0.0 506*7f296bb3SBarry Smith``` 507*7f296bb3SBarry Smith 508*7f296bb3SBarry SmithThe storage version is stored in the file and set automatically when loading (described below). 509*7f296bb3SBarry SmithYou can check the storage version of the HDF5 file with 510*7f296bb3SBarry Smith 511*7f296bb3SBarry Smith``` 512*7f296bb3SBarry Smith$ h5dump -a /dmplex_storage_version example.h5 513*7f296bb3SBarry Smith``` 514*7f296bb3SBarry Smith 515*7f296bb3SBarry SmithTo allow for simple and efficient implementation, and good load balancing, the 3.0.0 format changes the way the mesh topology is stored. 516*7f296bb3SBarry SmithDifferent strata (sets of mesh entities with an equal dimension; or vertices, edges, faces, and cells) are now stored separately. 517*7f296bb3SBarry SmithThe new structure of `/topologies/<mesh_name>/topology` is following: 518*7f296bb3SBarry Smith 519*7f296bb3SBarry Smith``` 520*7f296bb3SBarry Smith$ h5dump --contents example.h5 521*7f296bb3SBarry SmithHDF5 "example.h5" { 522*7f296bb3SBarry SmithFILE_CONTENTS { 523*7f296bb3SBarry Smith ... 524*7f296bb3SBarry Smith group /topologies/plexA/topology 525*7f296bb3SBarry Smith dataset /topologies/plexA/topology/permutation 526*7f296bb3SBarry Smith group /topologies/plexA/topology/strata 527*7f296bb3SBarry Smith group /topologies/plexA/topology/strata/0 528*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/0/cone_sizes 529*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/0/cones 530*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/0/orientations 531*7f296bb3SBarry Smith group /topologies/plexA/topology/strata/1 532*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/1/cone_sizes 533*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/1/cones 534*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/1/orientations 535*7f296bb3SBarry Smith group /topologies/plexA/topology/strata/2 536*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/2/cone_sizes 537*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/2/cones 538*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/2/orientations 539*7f296bb3SBarry Smith group /topologies/plexA/topology/strata/3 540*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/3/cone_sizes 541*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/3/cones 542*7f296bb3SBarry Smith dataset /topologies/plexA/topology/strata/3/orientations 543*7f296bb3SBarry Smith } 544*7f296bb3SBarry Smith} 545*7f296bb3SBarry Smith``` 546*7f296bb3SBarry Smith 547*7f296bb3SBarry SmithGroup `/topologies/<mesh_name>/topology/strata` contains a subgroup for each stratum depth (dimension; 0 for vertices up to 3 for cells). 548*7f296bb3SBarry SmithDAG points (mesh entities) have an implicit global numbering, given by the position in `orientations` (or `cone_sizes`) plus the stratum offset. 549*7f296bb3SBarry SmithThe stratum offset is given by a sum of lengths of all previous strata with respect to the order stored in `/topologies/<mesh_name>/topology/permutation`. 550*7f296bb3SBarry SmithThis global numbering is compatible with the explicit numbering in dataset `topology/order` of previous versions. 551*7f296bb3SBarry Smith 552*7f296bb3SBarry SmithFor a DAG point with index `i` at depth `s`, `cone_sizes[i]` gives a size of this point's cone (set of adjacent entities with depth `s-1`). 553*7f296bb3SBarry SmithLet `o = sum(cone_sizes[0:i]])` (in Python syntax). 554*7f296bb3SBarry SmithPoints forming the cone are then given by `cones[o:o+cone_sizes[i]]` (in numbering relative to the depth `s-1`). 555*7f296bb3SBarry SmithThe orientation of the cone with respect to point `i` is then stored in `orientations[i]`. 556*7f296bb3SBarry Smith 557*7f296bb3SBarry Smith### Loading 558*7f296bb3SBarry Smith 559*7f296bb3SBarry SmithTo load data from "example.h5" file, we create a `PetscViewer` 560*7f296bb3SBarry Smithof type `PETSCVIEWERHDF5` in `FILE_MODE_READ` mode as: 561*7f296bb3SBarry Smith 562*7f296bb3SBarry Smith``` 563*7f296bb3SBarry SmithPetscViewerHDF5Open(PETSC_COMM_WORLD, "example.h5", FILE_MODE_READ, &viewer); 564*7f296bb3SBarry Smith``` 565*7f296bb3SBarry Smith 566*7f296bb3SBarry SmithWe then create a `DMPLEX` object, give it a *mesh name*, "plexA", and load 567*7f296bb3SBarry Smiththe mesh as: 568*7f296bb3SBarry Smith 569*7f296bb3SBarry Smith``` 570*7f296bb3SBarry SmithDMCreate(PETSC_COMM_WORLD, &dm); 571*7f296bb3SBarry SmithDMSetType(dm, DMPLEX); 572*7f296bb3SBarry SmithPetscObjectSetName((PetscObject)dm, "plexA"); 573*7f296bb3SBarry SmithPetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_PETSC); 574*7f296bb3SBarry SmithDMLoad(dm, viewer); 575*7f296bb3SBarry SmithPetscViewerPopFormat(viewer); 576*7f296bb3SBarry Smith``` 577*7f296bb3SBarry Smith 578*7f296bb3SBarry Smithwhere `PETSC_VIEWER_HDF5_PETSC` format was again used. The user can have more control by replace the single load call with 579*7f296bb3SBarry Smith 580*7f296bb3SBarry Smith``` 581*7f296bb3SBarry SmithPetscSF sfO; 582*7f296bb3SBarry Smith 583*7f296bb3SBarry SmithDMCreate(PETSC_COMM_WORLD, &dm); 584*7f296bb3SBarry SmithDMSetType(dm, DMPLEX); 585*7f296bb3SBarry SmithPetscObjectSetName((PetscObject)dm, "plexA"); 586*7f296bb3SBarry SmithPetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_PETSC); 587*7f296bb3SBarry SmithDMPlexTopologyLoad(dm, viewer, &sfO); 588*7f296bb3SBarry SmithDMPlexCoordinatesLoad(dm, viewer, sfO); 589*7f296bb3SBarry SmithPetscViewerPopFormat(viewer); 590*7f296bb3SBarry Smith``` 591*7f296bb3SBarry Smith 592*7f296bb3SBarry SmithThe object returned by `DMPlexTopologyLoad()`, `sfO`, is a 593*7f296bb3SBarry Smith`PetscSF` that pushes forward $X$ to the loaded mesh, 594*7f296bb3SBarry Smith`dm`; this `PetscSF` is constructed with the global point 595*7f296bb3SBarry Smithnumber tags that we saved along with the mesh points. 596*7f296bb3SBarry Smith 597*7f296bb3SBarry SmithAs the `DMPLEX` mesh just loaded might not have a desired distribution, 598*7f296bb3SBarry Smithit is common to redistribute the mesh for a better distribution using 599*7f296bb3SBarry Smith`DMPlexDistribute()`, e.g., as: 600*7f296bb3SBarry Smith 601*7f296bb3SBarry Smith``` 602*7f296bb3SBarry SmithDM distributedDM; 603*7f296bb3SBarry SmithPetscInt overlap = 1; 604*7f296bb3SBarry SmithPetscSF sfDist, sf; 605*7f296bb3SBarry Smith 606*7f296bb3SBarry SmithDMPlexDistribute(dm, overlap, &sfDist, &distributedDM); 607*7f296bb3SBarry Smithif (distributedDM) { 608*7f296bb3SBarry Smith DMDestroy(&dm); 609*7f296bb3SBarry Smith dm = distributedDM; 610*7f296bb3SBarry Smith PetscObjectSetName((PetscObject)dm, "plexA"); 611*7f296bb3SBarry Smith} 612*7f296bb3SBarry SmithPetscSFCompose(sfO, sfDist, &sf); 613*7f296bb3SBarry SmithPetscSFDestroy(&sfO); 614*7f296bb3SBarry SmithPetscSFDestroy(&sfDist); 615*7f296bb3SBarry Smith``` 616*7f296bb3SBarry Smith 617*7f296bb3SBarry SmithNote that the new `DMPLEX` does not automatically inherit the *mesh name*, 618*7f296bb3SBarry Smithso we need to name it "plexA" once again. `sfDist` is a `PetscSF` 619*7f296bb3SBarry Smiththat pushes forward the loaded mesh to the redistributed mesh, so, composed 620*7f296bb3SBarry Smithwith `sfO`, it makes the `PetscSF` that pushes forward $X$ 621*7f296bb3SBarry Smithdirectly to the redistributed mesh, which we call `sf`. 622*7f296bb3SBarry Smith 623*7f296bb3SBarry SmithWe then create a new `DM`, `sdm`, with `DMClone()`, give it 624*7f296bb3SBarry Smitha *dm name*, "dmA", and load the on-disk data layout into `sdm` as: 625*7f296bb3SBarry Smith 626*7f296bb3SBarry Smith``` 627*7f296bb3SBarry SmithPetscSF globalDataSF, localDataSF; 628*7f296bb3SBarry Smith 629*7f296bb3SBarry SmithDMClone(dm, &sdm); 630*7f296bb3SBarry SmithPetscObjectSetName((PetscObject)sdm, "dmA"); 631*7f296bb3SBarry SmithDMPlexSectionLoad(dm, viewer, sdm, sf, &globalDataSF, &localDataSF); 632*7f296bb3SBarry Smith``` 633*7f296bb3SBarry Smith 634*7f296bb3SBarry Smithwhere we could also create a new 635*7f296bb3SBarry Smith`DMSHELL` object instead of using `DMClone()`. 636*7f296bb3SBarry SmithEach point in the on-disk data layout being tagged with a global 637*7f296bb3SBarry Smithpoint number in $X$, `DMPlexSectionLoad()` 638*7f296bb3SBarry Smithinternally constructs a `PetscSF` that pushes forward the on-disk 639*7f296bb3SBarry Smithdata layout to $X$. 640*7f296bb3SBarry SmithComposing this with `sf`, `DMPlexSectionLoad()` internally 641*7f296bb3SBarry Smithconstructs another `PetscSF` that pushes forward the on-disk 642*7f296bb3SBarry Smithdata layout directly to the redistributed mesh. It then 643*7f296bb3SBarry Smithreconstructs the data layout `s` on the redistributed mesh and 644*7f296bb3SBarry Smithattaches it to `sdm`. The objects returned by this function, 645*7f296bb3SBarry Smith`globalDataSF` and `localDataSF`, are `PetscSF`s that can 646*7f296bb3SBarry Smithbe used to migrate the on-disk vector data into local and global 647*7f296bb3SBarry Smith`Vec`s defined on `sdm`. 648*7f296bb3SBarry Smith 649*7f296bb3SBarry SmithWe now create a local vector associated with `sdm`, e.g., as: 650*7f296bb3SBarry Smith 651*7f296bb3SBarry Smith``` 652*7f296bb3SBarry SmithVec vec; 653*7f296bb3SBarry Smith 654*7f296bb3SBarry SmithDMGetLocalVector(sdm, &vec); 655*7f296bb3SBarry Smith``` 656*7f296bb3SBarry Smith 657*7f296bb3SBarry SmithWe then name `vec` "vecA" and load the on-disk vector into `vec` as: 658*7f296bb3SBarry Smith 659*7f296bb3SBarry Smith``` 660*7f296bb3SBarry SmithPetscObjectSetName((PetscObject)vec, "vecA"); 661*7f296bb3SBarry SmithDMPlexLocalVectorLoad(dm, viewer, sdm, localDataSF, localVec); 662*7f296bb3SBarry Smith``` 663*7f296bb3SBarry Smith 664*7f296bb3SBarry Smithwhere `localDataSF` knows how to migrate the on-disk vector 665*7f296bb3SBarry Smithdata into a local `Vec` defined on `sdm`. 666*7f296bb3SBarry SmithThe on-disk vector can be loaded into a global vector associated with 667*7f296bb3SBarry Smith`sdm` in the exact same way with trivial changes. 668*7f296bb3SBarry Smith 669*7f296bb3SBarry SmithAfter loading, we destroy the `PetscViewer` with: 670*7f296bb3SBarry Smith 671*7f296bb3SBarry Smith``` 672*7f296bb3SBarry SmithPetscViewerDestroy(&viewer); 673*7f296bb3SBarry Smith``` 674*7f296bb3SBarry Smith 675*7f296bb3SBarry SmithThe above infrastructure works seamlessly in distributed-memory parallel 676*7f296bb3SBarry Smithsettings, in which one can even use different number of processes for 677*7f296bb3SBarry Smithsaving and for loading; a more comprehensive example is found in 678*7f296bb3SBarry Smith<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/plex/tutorials/ex12.c.html">DMPlex Tutorial ex12</a>. 679*7f296bb3SBarry Smith 680*7f296bb3SBarry Smith(dm_adaptor_table)= 681*7f296bb3SBarry Smith 682*7f296bb3SBarry Smith## Metric-based mesh adaptation 683*7f296bb3SBarry Smith 684*7f296bb3SBarry Smith`DMPLEX` supports mesh adaptation using the *Riemannian metric framework*. 685*7f296bb3SBarry SmithThe idea is to use a Riemannian metric space within the mesher. The 686*7f296bb3SBarry Smithmetric space dictates how mesh resolution should be distributed across 687*7f296bb3SBarry Smiththe domain. Using this information, the remesher transforms the mesh such 688*7f296bb3SBarry Smiththat it is a *unit mesh* when viewed in the metric space. That is, the image 689*7f296bb3SBarry Smithof each of its elements under the mapping from Euclidean space into the 690*7f296bb3SBarry Smithmetric space has edges of unit length. 691*7f296bb3SBarry Smith 692*7f296bb3SBarry SmithOne of the main advantages of metric-based mesh adaptation is that it allows 693*7f296bb3SBarry Smithfor fully anisotropic remeshing. That is, it provides a means of controlling 694*7f296bb3SBarry Smiththe shape and orientation of elements in the adapted mesh, as well as their 695*7f296bb3SBarry Smithsize. This can be particularly useful for advection-dominated and 696*7f296bb3SBarry Smithdirectionally-dependent problems. 697*7f296bb3SBarry Smith 698*7f296bb3SBarry SmithSee {cite}`alauzet2010` for further details on metric-based anisotropic mesh 699*7f296bb3SBarry Smithadaptation. 700*7f296bb3SBarry Smith 701*7f296bb3SBarry SmithThe two main ingredients for metric-based mesh adaptation are an input mesh 702*7f296bb3SBarry Smith(i.e. the `DMPLEX`) and a Riemannian metric. The implementation in PETSc assumes 703*7f296bb3SBarry Smiththat the metric is piecewise linear and continuous across elemental boundaries. 704*7f296bb3SBarry SmithSuch an object can be created using the routine 705*7f296bb3SBarry Smith 706*7f296bb3SBarry Smith``` 707*7f296bb3SBarry SmithDMPlexMetricCreate(DM dm, PetscInt field, Vec *metric); 708*7f296bb3SBarry Smith``` 709*7f296bb3SBarry Smith 710*7f296bb3SBarry SmithA metric must be symmetric positive-definite, so that distances may be properly 711*7f296bb3SBarry Smithdefined. This can be checked using 712*7f296bb3SBarry Smith 713*7f296bb3SBarry Smith``` 714*7f296bb3SBarry SmithDMPlexMetricEnforceSPD(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant); 715*7f296bb3SBarry Smith``` 716*7f296bb3SBarry Smith 717*7f296bb3SBarry SmithThis routine may also be used to enforce minimum and maximum tolerated metric 718*7f296bb3SBarry Smithmagnitudes (i.e. cell sizes), as well as maximum anisotropy. These quantities 719*7f296bb3SBarry Smithcan be specified using 720*7f296bb3SBarry Smith 721*7f296bb3SBarry Smith``` 722*7f296bb3SBarry SmithDMPlexMetricSetMinimumMagnitude(DM dm, PetscReal h_min); 723*7f296bb3SBarry SmithDMPlexMetricSetMaximumMagnitude(DM dm, PetscReal h_max); 724*7f296bb3SBarry SmithDMPlexMetricSetMaximumAnisotropy(DM dm, PetscReal a_max); 725*7f296bb3SBarry Smith``` 726*7f296bb3SBarry Smith 727*7f296bb3SBarry Smithor the command line arguments 728*7f296bb3SBarry Smith 729*7f296bb3SBarry Smith``` 730*7f296bb3SBarry Smith-dm_plex_metric_h_min <h_min> 731*7f296bb3SBarry Smith-dm_plex_metric_h_max <h_max> 732*7f296bb3SBarry Smith-dm_plex_metric_a_max <a_max> 733*7f296bb3SBarry Smith``` 734*7f296bb3SBarry Smith 735*7f296bb3SBarry SmithOne simple way to combine two metrics is by simply averaging them entry-by-entry. 736*7f296bb3SBarry SmithAnother is to *intersect* them, which amounts to choosing the greatest level of 737*7f296bb3SBarry Smithrefinement in each direction. These operations are available in PETSc through 738*7f296bb3SBarry Smiththe routines 739*7f296bb3SBarry Smith 740*7f296bb3SBarry Smith``` 741*7f296bb3SBarry SmithDMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec metricAvg); 742*7f296bb3SBarry SmithDMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec metricInt); 743*7f296bb3SBarry Smith``` 744*7f296bb3SBarry Smith 745*7f296bb3SBarry SmithHowever, before combining metrics, it is important that they are scaled in the same 746*7f296bb3SBarry Smithway. Scaling also allows the user to control the number of vertices in the adapted 747*7f296bb3SBarry Smithmesh (in an approximate sense). This is achieved using the $L^p$ normalization 748*7f296bb3SBarry Smithframework, with the routine 749*7f296bb3SBarry Smith 750*7f296bb3SBarry Smith``` 751*7f296bb3SBarry SmithDMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, PetscBool restrictAnisotropy, Vec metricOut, Vec determinant); 752*7f296bb3SBarry Smith``` 753*7f296bb3SBarry Smith 754*7f296bb3SBarry SmithThere are two important parameters for the normalization: the normalization order 755*7f296bb3SBarry Smith$p$ and the target metric complexity, which is analogous to the vertex count. 756*7f296bb3SBarry SmithThey are controlled using 757*7f296bb3SBarry Smith 758*7f296bb3SBarry Smith``` 759*7f296bb3SBarry SmithDMPlexMetricSetNormalizationOrder(DM dm, PetscReal p); 760*7f296bb3SBarry SmithDMPlexMetricSetTargetComplexity(DM dm, PetscReal target); 761*7f296bb3SBarry Smith``` 762*7f296bb3SBarry Smith 763*7f296bb3SBarry Smithor the command line arguments 764*7f296bb3SBarry Smith 765*7f296bb3SBarry Smith```console 766*7f296bb3SBarry Smith-dm_plex_metric_p <p> 767*7f296bb3SBarry Smith-dm_plex_metric_target_complexity <target> 768*7f296bb3SBarry Smith``` 769*7f296bb3SBarry Smith 770*7f296bb3SBarry SmithTwo different metric-based mesh adaptation tools are available in PETSc: 771*7f296bb3SBarry Smith 772*7f296bb3SBarry Smith- [Pragmatic](https://meshadaptation.github.io/); 773*7f296bb3SBarry Smith- [Mmg/ParMmg](https://www.mmgtools.org/). 774*7f296bb3SBarry Smith 775*7f296bb3SBarry SmithMmg is a serial package, whereas ParMmg is the MPI version. 776*7f296bb3SBarry SmithNote that surface meshing is not currently supported and that ParMmg 777*7f296bb3SBarry Smithworks only in three dimensions. Mmg can be used for both two and three dimensional problems. 778*7f296bb3SBarry SmithPragmatic, Mmg and ParMmg may be specified by the command line arguments 779*7f296bb3SBarry Smith 780*7f296bb3SBarry Smith``` 781*7f296bb3SBarry Smith-dm_adaptor pragmatic 782*7f296bb3SBarry Smith-dm_adaptor mmg 783*7f296bb3SBarry Smith-dm_adaptor parmmg 784*7f296bb3SBarry Smith``` 785*7f296bb3SBarry Smith 786*7f296bb3SBarry SmithOnce a metric has been constructed, it can be used to perform metric-based 787*7f296bb3SBarry Smithmesh adaptation using the routine 788*7f296bb3SBarry Smith 789*7f296bb3SBarry Smith``` 790*7f296bb3SBarry SmithDMAdaptMetric(DM dm, Vec metric, DMLabel bdLabel, DMLabel rgLabel, DM dmAdapt); 791*7f296bb3SBarry Smith``` 792*7f296bb3SBarry Smith 793*7f296bb3SBarry Smithwhere `bdLabel` and `rgLabel` are boundary and interior tags to be 794*7f296bb3SBarry Smithpreserved under adaptation, respectively. 795*7f296bb3SBarry Smith 796*7f296bb3SBarry Smith```{rubric} Footnotes 797*7f296bb3SBarry Smith``` 798*7f296bb3SBarry Smith 799*7f296bb3SBarry Smith[^boundary-footnote]: In three dimensions, the boundary of a cell (sometimes called an element) is its faces, the boundary of a face is its edges and the boundary of an edge is the two vertices. 800*7f296bb3SBarry Smith 801*7f296bb3SBarry Smith```{eval-rst} 802*7f296bb3SBarry Smith.. bibliography:: /petsc.bib 803*7f296bb3SBarry Smith :filter: docname in docnames 804*7f296bb3SBarry Smith``` 805