1aa9a5b67SBarry Smith# Meshing for Subsurface Flows in PETSc 2aa9a5b67SBarry Smith 3aa9a5b67SBarry SmithThis tutorials guides users in creating meshes for the TDyCore simulation framework for subsurface flows. The user inputs a surface mesh, a refinement prescription, and an extrusion prescription in order to create the simulation mesh. 4aa9a5b67SBarry Smith 5aa9a5b67SBarry SmithReading the ASCII Output 6aa9a5b67SBarry Smith 7aa9a5b67SBarry SmithFor example, a very simple mesh would start with a square surface mesh divided into two triangles, which is then extruded to form two triangular prisms. This is the first test in the DMPlex tutorial code ex10, 8aa9a5b67SBarry Smith 9aa9a5b67SBarry Smith```console 10*30266083SMatthew G. Knepley$ make -f ./gmakefile test search="dm_impls_plex_tutorials-ex10_0" 11aa9a5b67SBarry Smith``` 12aa9a5b67SBarry Smith 13aa9a5b67SBarry Smithwhich outputs 14aa9a5b67SBarry Smith 15aa9a5b67SBarry Smith```{literalinclude} /../src/dm/impls/plex/tutorials/output/ex10_0.out 16aa9a5b67SBarry Smith``` 17aa9a5b67SBarry Smith 18aa9a5b67SBarry SmithWe can see that there are two 3-cells, meaning three-dimensional cells, and from the `celltype` label we see that those cells have celltype 9, meaning they are triangular prisms. The original surface mesh had 5 edges, so we would expect 10 edges for the two surfaces and four edges connecting those surfaces. This is exactly what we see, since there are 14 1-cells, but 4 of them noted in parentheses are tensor cells created by extrusion. We can see this another way in the celltype label, where there are ten mesh points of type 1, meaning segments, and four mesh points of type 2, meaning tensor products of a vertex and segment. Similarly, there are 9 2-cells, but 5 of them stretch between the two surfaces, meaning they are tensor products of two segments. 19aa9a5b67SBarry Smith 20aa9a5b67SBarry SmithRegular Refinement of Simplex Meshes 21aa9a5b67SBarry Smith 22aa9a5b67SBarry SmithWe can regularly refine the surface before extrusion using `-dm_refine <k>`, where `k` is the number of refinements, 23aa9a5b67SBarry Smith 24aa9a5b67SBarry Smith```console 25*30266083SMatthew G. Knepley$ make -f ./gmakefile test search="dm_impls_plex_tutorials-ex10_1" EXTRA_OPTIONS="-srf_dm_refine 2 -srf_dm_view draw -draw_save $PETSC_DIR/surface.png -draw_save_single_file" 26aa9a5b67SBarry Smith``` 27aa9a5b67SBarry Smith 28aa9a5b67SBarry Smithwhich produces the following surface 29aa9a5b67SBarry Smith 30aa9a5b67SBarry Smith:::{figure} /images/tutorials/meshing/surface.png 31aa9a5b67SBarry Smith:align: center 32aa9a5b67SBarry Smith 33aa9a5b67SBarry Smith**Surface mesh refined twice** 34aa9a5b67SBarry Smith::: 35aa9a5b67SBarry Smith 36aa9a5b67SBarry Smithand the extruded mesh can be visualized using VTK. Here I make the image using Paraview, and give the extrusion 3 layers 37aa9a5b67SBarry Smith 38aa9a5b67SBarry Smith```console 39*30266083SMatthew G. Knepley$ make -f ./gmakefile test search="dm_impls_plex_tutorials-ex10_1" EXTRA_OPTIONS="-dm_view hdf5:$PETSC_DIR/mesh.h5 -dm_extrude 3" 40aa9a5b67SBarry Smith$ $PETSC_DIR/lib/petsc/bin/petsc_gen_xmdf.py mesh.h5 41aa9a5b67SBarry Smith``` 42aa9a5b67SBarry Smith 43aa9a5b67SBarry Smith:::{figure} /images/tutorials/meshing/extrusion.png 44aa9a5b67SBarry Smith:align: center 45aa9a5b67SBarry Smith 46aa9a5b67SBarry Smith**Extruded mesh with refined surface** 47aa9a5b67SBarry Smith::: 48aa9a5b67SBarry Smith 49aa9a5b67SBarry SmithWe can similarly look at this in parallel. Test 2 uses three refinements and three extrusion layers on five processes 50aa9a5b67SBarry Smith 51aa9a5b67SBarry Smith```console 52*30266083SMatthew G. Knepley$ make -f ./gmakefile test search="dm_impls_plex_tutorials-ex10_2" EXTRA_OPTIONS="-dm_view hdf5:$PETSC_DIR/mesh.h5 -dm_partition_view -petscpartitioner_type parmetis" 53aa9a5b67SBarry Smith$ $PETSC_DIR/lib/petsc/bin/petsc_gen_xmdf.py mesh.h5 54aa9a5b67SBarry Smith``` 55aa9a5b67SBarry Smith 56aa9a5b67SBarry Smith:::{figure} /images/tutorials/meshing/extrusionParallel.png 57aa9a5b67SBarry Smith:align: center 58aa9a5b67SBarry Smith 59aa9a5b67SBarry Smith**Parallel extruded mesh with refined surface** 60aa9a5b67SBarry Smith::: 61aa9a5b67SBarry Smith 62aa9a5b67SBarry SmithAdaptive Refinement of Simplex Meshes 63aa9a5b67SBarry Smith 64aa9a5b67SBarry SmithAdaptive refinement of simplicial meshes is somewhat tricky when we demand that the meshes be conforming, as we do in this case. We would like different grid cells to have different levels of refinement, for example headwaters cells in a watershed be refined twice, while river channel cells be refined four times. In order to differentiate between cells, we first mark the cells on the surface using a `DMLabel`. We can do this programmatically, 65aa9a5b67SBarry Smith 66aa9a5b67SBarry Smith```{literalinclude} /../src/dm/impls/plex/tutorials/ex10.c 67aa9a5b67SBarry Smith:append: '}' 68aa9a5b67SBarry Smith:end-at: PetscFunctionReturn(PETSC_SUCCESS); 69aa9a5b67SBarry Smith:start-at: static PetscErrorCode CreateDomainLabel( 70aa9a5b67SBarry Smith``` 71aa9a5b67SBarry Smith 72aa9a5b67SBarry Smithor you can label the mesh using a GUI, such as GMsh, and PETSc will read the label values from the input file. 73aa9a5b67SBarry Smith 74aa9a5b67SBarry SmithWe next create a label marking each cell in the mesh with an action, such as `DM_ADAPT_REFINE` or `DM_ADAPT_COARSEN`. We do this based on a volume constraint, namely that cells with a certain label value should have a certain volume. You could, of course, choose a more complex strategy, but here we just want a clear criterion. We can give volume constraints for label value `v` using the command line argument `-volume_constraint_<v> <vol>`. The mesh is then refined iteratively, checking the volume constraints each time, 75aa9a5b67SBarry Smith 76aa9a5b67SBarry Smith```{literalinclude} /../src/dm/impls/plex/tutorials/ex10.c 77aa9a5b67SBarry Smith:append: '}' 78aa9a5b67SBarry Smith:end-at: PetscCall(DMLabelDestroy(&adaptLabel)); 79aa9a5b67SBarry Smith:start-at: while (adapt) { 80aa9a5b67SBarry Smith``` 81aa9a5b67SBarry Smith 82aa9a5b67SBarry SmithTest 3 from `ex10` constrains the headwater cells (with marker 1) to have volume less than 0.01, and the river channel cells (with marker 2) to be smaller than 0.000625 83aa9a5b67SBarry Smith 84aa9a5b67SBarry Smith```{literalinclude} /../src/dm/impls/plex/tutorials/ex10.c 85aa9a5b67SBarry Smith:lines: 1-3 86aa9a5b67SBarry Smith:start-at: 'suffix: 3' 87aa9a5b67SBarry Smith``` 88aa9a5b67SBarry Smith 89aa9a5b67SBarry SmithWe can look at a parallel run using extra options for the test system 90aa9a5b67SBarry Smith 91aa9a5b67SBarry Smith```console 92*30266083SMatthew G. Knepley$ make -f ./gmakefile test search="dm_impls_plex_tutorials-ex10_3" EXTRA_OPTIONS="-dm_view hdf5:$PETSC_DIR/mesh.h5 -dm_partition_view -petscpartitioner_type parmetis" NP=5 93aa9a5b67SBarry Smith$ $PETSC_DIR/lib/petsc/bin/petsc_gen_xmdf.py mesh.h5 94aa9a5b67SBarry Smith``` 95aa9a5b67SBarry Smith 96aa9a5b67SBarry Smith:::{figure} /images/tutorials/meshing/extrusionAdaptiveParallel.png 97aa9a5b67SBarry Smith:align: center 98aa9a5b67SBarry Smith 99aa9a5b67SBarry Smith**Parallel extruded mesh with adaptively refined surface** 100aa9a5b67SBarry Smith::: 101aa9a5b67SBarry Smith 102aa9a5b67SBarry SmithBy turning on `PetscInfo`, we can see what decisions the refiner is making 103aa9a5b67SBarry Smith 104aa9a5b67SBarry Smith```console 105*30266083SMatthew G. Knepley$ make -f ./gmakefile test search="dm_impls_plex_tutorials-ex10_3" EXTRA_OPTIONS="-info :dm" 106aa9a5b67SBarry Smith# > [0] AdaptMesh(): Adapted mesh, marking 12 cells for refinement, and 0 cells for coarsening 107aa9a5b67SBarry Smith# > [0] AdaptMesh(): Adapted mesh, marking 29 cells for refinement, and 0 cells for coarsening 108aa9a5b67SBarry Smith# > [0] AdaptMesh(): Adapted mesh, marking 84 cells for refinement, and 0 cells for coarsening 109aa9a5b67SBarry Smith# > [0] AdaptMesh(): Adapted mesh, marking 10 cells for refinement, and 0 cells for coarsening 110aa9a5b67SBarry Smith``` 111