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