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