xref: /honee/src/setupdm.c (revision 78a2675014c975c848215d7d366517830c162f9d)
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   PetscBool isCGNS;
17   char      filename[PETSC_MAX_PATH_LEN] = "";
18 
19   PetscFunctionBeginUser;
20   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL));
21   PetscCall(HoneeCheckFilenameExtension(comm, filename, ".cgns", &isCGNS));
22 
23   if (isCGNS) {
24     // Must create from file to keep the sfNatural field in tact
25     PetscCall(DMPlexCreateFromFile(comm, filename, "HONEE", PETSC_TRUE, dm));
26     PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_filename", ""));
27   } else {
28     PetscCall(DMCreate(comm, dm));
29     PetscCall(DMSetType(*dm, DMPLEX));
30   }
31 
32   {
33     PetscBool skip = PETSC_TRUE;
34     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
35     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
36   }
37   PetscCall(DMSetMatType(*dm, mat_type));
38   PetscCall(DMSetVecType(*dm, vec_type));
39   PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_simplex", "0"));
40   PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_localize", "0"));  // Delay localization until DMSetupByOrderEnd_FEM
41   PetscCall(DMSetSparseLocalize(*dm, PETSC_FALSE));
42   PetscCall(DMSetBlockingType(*dm, DM_BLOCKING_FIELD_NODE));
43 
44   // Set CL options
45   PetscCall(DMSetFromOptions(*dm));
46   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
47   {  // Force isoperiodic point SF to be created to update sfNatural.
48     // Needs to be done before removing the field corresponding to sfNatural, but after IsoperiodicFaces have been set
49     PetscInt num_fields;
50     PetscCall(DMGetNumFields(*dm, &num_fields));
51     if (num_fields) {
52       PetscSection dummy;
53       PetscCall(DMGetGlobalSection(*dm, &dummy));
54     }
55   }
56   PetscFunctionReturn(PETSC_SUCCESS);
57 }
58 
59 // Setup DM
60 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
61   PetscInt num_comp_q = 5;
62 
63   PetscFunctionBeginUser;
64   PetscCall(DMClearFields(dm));
65   PetscCall(DMSetLocalSection(dm, NULL));
66   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
67 
68   {  // Add strong boundary conditions to DM
69     DMLabel label;
70     PetscCall(DMGetLabel(dm, "Face Sets", &label));
71     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label");
72     PetscCall(DMPlexLabelComplete(dm, label));
73 
74     for (PetscInt i = 0; i < problem->num_bc_defs; i++) {
75       BCDefinition    bc_def = problem->bc_defs[i];
76       PetscInt        num_essential_comps, num_label_values;
77       const PetscInt *essential_comps, *label_values;
78       const char     *name;
79 
80       PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps));
81       if (essential_comps > 0) {
82         PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values));
83         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL,
84                                 NULL, NULL));
85       }
86     }
87     {
88       PetscBool use_strongstg = PETSC_FALSE;
89       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
90       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
91     }
92   }
93 
94   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
95 
96   // Empty name for conserved field (because there is only one field)
97   PetscSection section;
98   PetscCall(DMGetLocalSection(dm, &section));
99   PetscCall(PetscSectionSetFieldName(section, 0, ""));
100   switch (phys->state_var) {
101     case STATEVAR_CONSERVATIVE:
102       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
103       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
104       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
105       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
106       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
107       break;
108 
109     case STATEVAR_PRIMITIVE:
110       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
111       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
112       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
113       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
114       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
115       break;
116 
117     case STATEVAR_ENTROPY:
118       PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity"));
119       PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX"));
120       PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY"));
121       PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ"));
122       PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy"));
123       break;
124   }
125   PetscFunctionReturn(PETSC_SUCCESS);
126 }
127 
128 // Refine DM for high-order viz
129 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys) {
130   DM      dm_hierarchy[honee->app_ctx->viz_refine + 1];
131   VecType vec_type;
132 
133   PetscFunctionBeginUser;
134   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
135 
136   dm_hierarchy[0] = dm;
137   for (PetscInt i = 0, d = honee->app_ctx->degree; i < honee->app_ctx->viz_refine; i++) {
138     Mat interp_next;
139     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
140     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
141     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
142     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
143     d                = (d + 1) / 2;
144     PetscInt q_order = d + honee->app_ctx->q_extra;
145     if (i + 1 == honee->app_ctx->viz_refine) d = 1;
146     PetscCall(DMGetVecType(dm, &vec_type));
147     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
148     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
149     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
150     if (!i) honee->interp_viz = interp_next;
151     else {
152       Mat C;
153       PetscCall(MatMatMult(interp_next, honee->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
154       PetscCall(MatDestroy(&interp_next));
155       PetscCall(MatDestroy(&honee->interp_viz));
156       honee->interp_viz = C;
157     }
158   }
159   for (PetscInt i = 1; i < honee->app_ctx->viz_refine; i++) {
160     PetscCall(DMDestroy(&dm_hierarchy[i]));
161   }
162   honee->dm_viz = dm_hierarchy[honee->app_ctx->viz_refine];
163   PetscFunctionReturn(PETSC_SUCCESS);
164 }
165