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