xref: /honee/src/setupdm.c (revision eb9b4fe1a6357cb57da823758ee8a00861ee7fff)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdmplex.h>
967263decSKenneth E. Jansen #include <petscds.h>
10e419654dSJeremy L Thompson 
11149fb536SJames Wright #include <navierstokes.h>
12866f23abSJames Wright #include "../problems/stg_shur14.h"
13a515125bSLeila Ghaffari 
1405a512bdSLeila Ghaffari // Create mesh
15991aef52SJames Wright PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType mat_type, VecType vec_type, DM *dm) {
16*eb9b4fe1SJames Wright   PetscBool isCGNS;
17*eb9b4fe1SJames Wright   char      filename[PETSC_MAX_PATH_LEN] = "";
18*eb9b4fe1SJames Wright 
19a515125bSLeila Ghaffari   PetscFunctionBeginUser;
20*eb9b4fe1SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL));
21*eb9b4fe1SJames Wright   PetscCall(HoneeCheckFilenameExtension(comm, filename, ".cgns", &isCGNS));
22*eb9b4fe1SJames Wright 
23*eb9b4fe1SJames Wright   if (isCGNS) {
24*eb9b4fe1SJames Wright     // Must create from file to keep the sfNatural field in tact
25*eb9b4fe1SJames Wright     PetscCall(DMPlexCreateFromFile(comm, filename, "HONEE", PETSC_TRUE, dm));
26*eb9b4fe1SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_filename", ""));
27*eb9b4fe1SJames Wright   } else {
282b916ea7SJeremy L Thompson     PetscCall(DMCreate(comm, dm));
292b916ea7SJeremy L Thompson     PetscCall(DMSetType(*dm, DMPLEX));
30*eb9b4fe1SJames Wright   }
31*eb9b4fe1SJames Wright 
32b150a244SJed Brown   {
33b150a244SJed Brown     PetscBool skip = PETSC_TRUE;
342b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
35b150a244SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
36b150a244SJed Brown   }
372b916ea7SJeremy L Thompson   PetscCall(DMSetMatType(*dm, mat_type));
382b916ea7SJeremy L Thompson   PetscCall(DMSetVecType(*dm, vec_type));
3959572072SJames Wright   PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_plex_simplex", "0"));
40e747eef9SJames Wright   PetscCall(HoneeOptionsSetValueDefault(NULL, "-dm_localize", "0"));  // Delay localization until DMSetupByOrderEnd_FEM
41e747eef9SJames Wright   PetscCall(DMSetSparseLocalize(*dm, PETSC_FALSE));
42e747eef9SJames Wright   PetscCall(DMSetBlockingType(*dm, DM_BLOCKING_FIELD_NODE));
43ebfabadfSJames Wright 
4405a512bdSLeila Ghaffari   // Set CL options
452b916ea7SJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
462b916ea7SJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
47*eb9b4fe1SJames Wright   {  // Force isoperiodic point SF to be created to update sfNatural.
48*eb9b4fe1SJames Wright     // Needs to be done before removing the field corresponding to sfNatural, but after IsoperiodicFaces have been set
49*eb9b4fe1SJames Wright     PetscSection dummy;
50*eb9b4fe1SJames Wright     PetscCall(DMGetGlobalSection(*dm, &dummy));
51*eb9b4fe1SJames Wright   }
52d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
53a515125bSLeila Ghaffari }
54a515125bSLeila Ghaffari 
55a515125bSLeila Ghaffari // Setup DM
56991aef52SJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
57a515125bSLeila Ghaffari   PetscInt num_comp_q = 5;
58da4ca0cfSJames Wright 
59*eb9b4fe1SJames Wright   PetscFunctionBeginUser;
60*eb9b4fe1SJames Wright   PetscCall(DMClearFields(dm));
61*eb9b4fe1SJames Wright   PetscCall(DMSetLocalSection(dm, NULL));
62b2e5b5b3SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
63da4ca0cfSJames Wright 
64da4ca0cfSJames Wright   {  // Add strong boundary conditions to DM
65047c2946SJed Brown     DMLabel label;
662b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
675892ddc8SJames Wright     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label");
68c5e9980aSAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
69487a3b6eSJames Wright 
70487a3b6eSJames Wright     for (PetscInt i = 0; i < problem->num_bc_defs; i++) {
71487a3b6eSJames Wright       BCDefinition    bc_def = problem->bc_defs[i];
72487a3b6eSJames Wright       PetscInt        num_essential_comps, num_label_values;
73487a3b6eSJames Wright       const PetscInt *essential_comps, *label_values;
74487a3b6eSJames Wright       const char     *name;
75487a3b6eSJames Wright 
76487a3b6eSJames Wright       PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps));
77487a3b6eSJames Wright       if (essential_comps > 0) {
78487a3b6eSJames Wright         PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values));
79487a3b6eSJames Wright         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL,
80487a3b6eSJames Wright                                 NULL, NULL));
81002797a3SLeila Ghaffari       }
82002797a3SLeila Ghaffari     }
83866f23abSJames Wright     {
84866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
852b916ea7SJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
8642454adaSJames Wright       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
87866f23abSJames Wright     }
8867263decSKenneth E. Jansen   }
8967263decSKenneth E. Jansen 
90da4ca0cfSJames Wright   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
91cbe60e31SLeila Ghaffari 
92a515125bSLeila Ghaffari   // Empty name for conserved field (because there is only one field)
93a515125bSLeila Ghaffari   PetscSection section;
942b916ea7SJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
952b916ea7SJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
963636f6a4SJames Wright   switch (phys->state_var) {
973636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
982b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
992b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
1002b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
1012b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
102525869caSJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
1033636f6a4SJames Wright       break;
1043636f6a4SJames Wright 
1053636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
1062b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
1072b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
1082b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
1092b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
1102b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
1113636f6a4SJames Wright       break;
1129b103f75SJames Wright 
1139b103f75SJames Wright     case STATEVAR_ENTROPY:
1149b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity"));
1159b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX"));
1169b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY"));
1179b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ"));
1189b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy"));
1199b103f75SJames Wright       break;
120a515125bSLeila Ghaffari   }
121d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
122a515125bSLeila Ghaffari }
123a515125bSLeila Ghaffari 
124a515125bSLeila Ghaffari // Refine DM for high-order viz
1250c373b74SJames Wright PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys) {
1260c373b74SJames Wright   DM      dm_hierarchy[honee->app_ctx->viz_refine + 1];
127a515125bSLeila Ghaffari   VecType vec_type;
128a515125bSLeila Ghaffari 
12906f41313SJames Wright   PetscFunctionBeginUser;
1302b916ea7SJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
131a515125bSLeila Ghaffari 
132a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
1330c373b74SJames Wright   for (PetscInt i = 0, d = honee->app_ctx->degree; i < honee->app_ctx->viz_refine; i++) {
134a515125bSLeila Ghaffari     Mat interp_next;
1352b916ea7SJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1362b916ea7SJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1372b916ea7SJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1382b916ea7SJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
139a515125bSLeila Ghaffari     d                = (d + 1) / 2;
1400c373b74SJames Wright     PetscInt q_order = d + honee->app_ctx->q_extra;
1410c373b74SJames Wright     if (i + 1 == honee->app_ctx->viz_refine) d = 1;
1422b916ea7SJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1432b916ea7SJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
14467263decSKenneth E. Jansen     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
1452b916ea7SJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
1460c373b74SJames Wright     if (!i) honee->interp_viz = interp_next;
147a515125bSLeila Ghaffari     else {
148a515125bSLeila Ghaffari       Mat C;
1490c373b74SJames Wright       PetscCall(MatMatMult(interp_next, honee->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1502b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1510c373b74SJames Wright       PetscCall(MatDestroy(&honee->interp_viz));
1520c373b74SJames Wright       honee->interp_viz = C;
153a515125bSLeila Ghaffari     }
154a515125bSLeila Ghaffari   }
1550c373b74SJames Wright   for (PetscInt i = 1; i < honee->app_ctx->viz_refine; i++) {
1562b916ea7SJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
157a515125bSLeila Ghaffari   }
1580c373b74SJames Wright   honee->dm_viz = dm_hierarchy[honee->app_ctx->viz_refine];
159d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
160a515125bSLeila Ghaffari }
161