xref: /libCEED/examples/fluids/src/setupdm.c (revision 1f97d2f18688a5caa5dd9f004da3151e92048b76)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Setup DM for Navier-Stokes example using PETSc
10 
11 #include <ceed.h>
12 #include <petscdmplex.h>
13 
14 #include "../navierstokes.h"
15 #include "../problems/stg_shur14.h"
16 
17 // Create mesh
18 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) {
19   PetscFunctionBeginUser;
20   // Create DMPLEX
21   PetscCall(DMCreate(comm, dm));
22   PetscCall(DMSetType(*dm, DMPLEX));
23   {
24     PetscBool skip = PETSC_TRUE;
25     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
26     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
27   }
28   PetscCall(DMSetMatType(*dm, mat_type));
29   PetscCall(DMSetVecType(*dm, vec_type));
30 
31   // Set Tensor elements
32   PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"));
33   PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"));
34   // Set CL options
35   PetscCall(DMSetFromOptions(*dm));
36   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
37   PetscFunctionReturn(PETSC_SUCCESS);
38 }
39 
40 // Setup DM
41 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys) {
42   PetscFunctionBeginUser;
43   {
44     // Configure the finite element space and boundary conditions
45     PetscFE  fe;
46     PetscInt num_comp_q = 5;
47     DMLabel  label;
48     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
49     PetscCall(PetscObjectSetName((PetscObject)fe, "Q"));
50     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
51     PetscCall(DMCreateDS(dm));
52     PetscCall(DMGetLabel(dm, "Face Sets", &label));
53     PetscCall(DMPlexLabelComplete(dm, label));
54     // Set wall BCs
55     if (bc->num_wall > 0) {
56       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, NULL, NULL, NULL, NULL));
57     }
58     // Set slip BCs in the x direction
59     if (bc->num_slip[0] > 0) {
60       PetscInt comps[1] = {1};
61       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, NULL, NULL, NULL, NULL));
62     }
63     // Set slip BCs in the y direction
64     if (bc->num_slip[1] > 0) {
65       PetscInt comps[1] = {2};
66       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, NULL, NULL, NULL, NULL));
67     }
68     // Set slip BCs in the z direction
69     if (bc->num_slip[2] > 0) {
70       PetscInt comps[1] = {3};
71       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, NULL, NULL, NULL, NULL));
72     }
73     {
74       PetscBool use_strongstg = PETSC_FALSE;
75       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
76       if (use_strongstg) PetscCall(SetupStrongSTG(dm, bc, problem, phys));
77     }
78 
79     PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
80     PetscCall(PetscFEDestroy(&fe));
81   }
82 
83   // Empty name for conserved field (because there is only one field)
84   PetscSection section;
85   PetscCall(DMGetLocalSection(dm, &section));
86   PetscCall(PetscSectionSetFieldName(section, 0, ""));
87   switch (phys->state_var) {
88     case STATEVAR_CONSERVATIVE:
89       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
90       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Momentum X"));
91       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Momentum Y"));
92       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Momentum Z"));
93       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Energy Density"));
94       break;
95 
96     case STATEVAR_PRIMITIVE:
97       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
98       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Velocity X"));
99       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Velocity Y"));
100       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Velocity Z"));
101       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
102       break;
103   }
104   PetscFunctionReturn(PETSC_SUCCESS);
105 }
106 
107 // Refine DM for high-order viz
108 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
109   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
110   VecType vec_type;
111   PetscFunctionBeginUser;
112 
113   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
114 
115   dm_hierarchy[0] = dm;
116   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
117     Mat interp_next;
118     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
119     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
120     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
121     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
122     d = (d + 1) / 2;
123     if (i + 1 == user->app_ctx->viz_refine) d = 1;
124     PetscCall(DMGetVecType(dm, &vec_type));
125     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
126     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, bc, phys));
127     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
128     if (!i) user->interp_viz = interp_next;
129     else {
130       Mat C;
131       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
132       PetscCall(MatDestroy(&interp_next));
133       PetscCall(MatDestroy(&user->interp_viz));
134       user->interp_viz = C;
135     }
136   }
137   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
138     PetscCall(DMDestroy(&dm_hierarchy[i]));
139   }
140   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
141 
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144