xref: /honee/src/setupdm.c (revision f1b4fcbf63a2bebf374a2f2a4f174597200260f4)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Setup DM for Navier-Stokes example using PETSc
6 
7 #include <ceed.h>
8 #include <petscdmplex.h>
9 #include <petscds.h>
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(HoneeOptionsSetValueDefault(NULL, "-dm_plex_simplex", "0"));
30   PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_localize", "0"));  // Delay localization until DMSetupByOrderEnd_FEM
31   PetscCall(DMSetSparseLocalize(*dm, PETSC_FALSE));
32   PetscCall(DMSetBlockingType(*dm, DM_BLOCKING_FIELD_NODE));
33 
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, PetscInt q_extra, SimpleBC bc, Physics phys) {
42   PetscInt num_comp_q = 5;
43   PetscFunctionBeginUser;
44 
45   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
46 
47   {  // Add strong boundary conditions to DM
48     DMLabel label;
49     PetscCall(DMGetLabel(dm, "Face Sets", &label));
50     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label");
51     PetscCall(DMPlexLabelComplete(dm, label));
52 
53     for (PetscInt i = 0; i < problem->num_bc_defs; i++) {
54       BCDefinition    bc_def = problem->bc_defs[i];
55       PetscInt        num_essential_comps, num_label_values;
56       const PetscInt *essential_comps, *label_values;
57       const char     *name;
58 
59       PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps));
60       if (essential_comps > 0) {
61         PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values));
62         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL,
63                                 NULL, NULL));
64       }
65     }
66     {
67       PetscBool use_strongstg = PETSC_FALSE;
68       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
69       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
70     }
71   }
72 
73   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
74 
75   // Empty name for conserved field (because there is only one field)
76   PetscSection section;
77   PetscCall(DMGetLocalSection(dm, &section));
78   PetscCall(PetscSectionSetFieldName(section, 0, ""));
79   switch (phys->state_var) {
80     case STATEVAR_CONSERVATIVE:
81       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
82       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
83       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
84       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
85       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
86       break;
87 
88     case STATEVAR_PRIMITIVE:
89       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
90       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
91       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
92       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
93       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
94       break;
95 
96     case STATEVAR_ENTROPY:
97       PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity"));
98       PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX"));
99       PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY"));
100       PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ"));
101       PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy"));
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 
112   PetscFunctionBeginUser;
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     PetscInt q_order = d + user->app_ctx->q_extra;
124     if (i + 1 == user->app_ctx->viz_refine) d = 1;
125     PetscCall(DMGetVecType(dm, &vec_type));
126     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
127     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
128     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
129     if (!i) user->interp_viz = interp_next;
130     else {
131       Mat C;
132       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
133       PetscCall(MatDestroy(&interp_next));
134       PetscCall(MatDestroy(&user->interp_viz));
135       user->interp_viz = C;
136     }
137   }
138   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
139     PetscCall(DMDestroy(&dm_hierarchy[i]));
140   }
141   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144