xref: /honee/src/setupdm.c (revision ea615d4cc464aa6ad650c06fae6d120cc2465bc4)
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
5*ea615d4cSJames Wright /// Setup DM for HONEE
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) {
16eb9b4fe1SJames Wright   PetscBool isCGNS;
17eb9b4fe1SJames Wright   char      filename[PETSC_MAX_PATH_LEN] = "";
18eb9b4fe1SJames Wright 
19a515125bSLeila Ghaffari   PetscFunctionBeginUser;
20eb9b4fe1SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL));
21eb9b4fe1SJames Wright   PetscCall(HoneeCheckFilenameExtension(comm, filename, ".cgns", &isCGNS));
22eb9b4fe1SJames Wright 
23eb9b4fe1SJames Wright   if (isCGNS) {
24eb9b4fe1SJames Wright     // Must create from file to keep the sfNatural field in tact
25eb9b4fe1SJames Wright     PetscCall(DMPlexCreateFromFile(comm, filename, "HONEE", PETSC_TRUE, dm));
26eb9b4fe1SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_filename", ""));
27eb9b4fe1SJames Wright   } else {
282b916ea7SJeremy L Thompson     PetscCall(DMCreate(comm, dm));
292b916ea7SJeremy L Thompson     PetscCall(DMSetType(*dm, DMPLEX));
30eb9b4fe1SJames Wright   }
31eb9b4fe1SJames 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));
43e7934809SJames Wright   PetscCall(DMPlexSetPartitionBalance(*dm, PETSC_TRUE));
44ebfabadfSJames Wright 
4505a512bdSLeila Ghaffari   // Set CL options
462b916ea7SJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
472b916ea7SJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
48eb9b4fe1SJames Wright   {  // Force isoperiodic point SF to be created to update sfNatural.
49eb9b4fe1SJames Wright     // Needs to be done before removing the field corresponding to sfNatural, but after IsoperiodicFaces have been set
50c4a0f6c7SJames Wright     PetscInt num_fields;
51c4a0f6c7SJames Wright     PetscCall(DMGetNumFields(*dm, &num_fields));
52c4a0f6c7SJames Wright     if (num_fields) {
53eb9b4fe1SJames Wright       PetscSection dummy;
54eb9b4fe1SJames Wright       PetscCall(DMGetGlobalSection(*dm, &dummy));
55eb9b4fe1SJames Wright     }
56c4a0f6c7SJames Wright   }
57d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
58a515125bSLeila Ghaffari }
59a515125bSLeila Ghaffari 
60a515125bSLeila Ghaffari // Setup DM
61991aef52SJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
62a515125bSLeila Ghaffari   PetscInt num_comp_q = 5;
63da4ca0cfSJames Wright 
64eb9b4fe1SJames Wright   PetscFunctionBeginUser;
65eb9b4fe1SJames Wright   PetscCall(DMClearFields(dm));
66eb9b4fe1SJames Wright   PetscCall(DMSetLocalSection(dm, NULL));
67b2e5b5b3SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
68da4ca0cfSJames Wright 
69da4ca0cfSJames Wright   {  // Add strong boundary conditions to DM
70047c2946SJed Brown     DMLabel label;
712b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
725892ddc8SJames Wright     PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DM does not have 'Face Sets' label");
73c5e9980aSAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
74487a3b6eSJames Wright 
75487a3b6eSJames Wright     for (PetscInt i = 0; i < problem->num_bc_defs; i++) {
76487a3b6eSJames Wright       BCDefinition    bc_def = problem->bc_defs[i];
77487a3b6eSJames Wright       PetscInt        num_essential_comps, num_label_values;
78487a3b6eSJames Wright       const PetscInt *essential_comps, *label_values;
79487a3b6eSJames Wright       const char     *name;
80487a3b6eSJames Wright 
81487a3b6eSJames Wright       PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps));
82487a3b6eSJames Wright       if (essential_comps > 0) {
83487a3b6eSJames Wright         PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values));
84487a3b6eSJames Wright         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL,
85487a3b6eSJames Wright                                 NULL, NULL));
86002797a3SLeila Ghaffari       }
87002797a3SLeila Ghaffari     }
88866f23abSJames Wright     {
89866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
902b916ea7SJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
9142454adaSJames Wright       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
92866f23abSJames Wright     }
9367263decSKenneth E. Jansen   }
9467263decSKenneth E. Jansen 
95da4ca0cfSJames Wright   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
96cbe60e31SLeila Ghaffari 
97a515125bSLeila Ghaffari   // Empty name for conserved field (because there is only one field)
98a515125bSLeila Ghaffari   PetscSection section;
992b916ea7SJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
1002b916ea7SJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
1013636f6a4SJames Wright   switch (phys->state_var) {
1023636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
1032b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
1042b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
1052b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
1062b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
107525869caSJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
1083636f6a4SJames Wright       break;
1093636f6a4SJames Wright 
1103636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
1112b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
1122b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
1132b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
1142b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
1152b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
1163636f6a4SJames Wright       break;
1179b103f75SJames Wright 
1189b103f75SJames Wright     case STATEVAR_ENTROPY:
1199b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity"));
1209b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX"));
1219b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY"));
1229b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ"));
1239b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy"));
1249b103f75SJames Wright       break;
125a515125bSLeila Ghaffari   }
126d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
127a515125bSLeila Ghaffari }
128a515125bSLeila Ghaffari 
129a515125bSLeila Ghaffari // Refine DM for high-order viz
1300c373b74SJames Wright PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys) {
1310c373b74SJames Wright   DM      dm_hierarchy[honee->app_ctx->viz_refine + 1];
132a515125bSLeila Ghaffari   VecType vec_type;
133a515125bSLeila Ghaffari 
13406f41313SJames Wright   PetscFunctionBeginUser;
1352b916ea7SJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
136a515125bSLeila Ghaffari 
137a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
1380c373b74SJames Wright   for (PetscInt i = 0, d = honee->app_ctx->degree; i < honee->app_ctx->viz_refine; i++) {
139a515125bSLeila Ghaffari     Mat interp_next;
1402b916ea7SJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1412b916ea7SJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1422b916ea7SJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1432b916ea7SJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
144a515125bSLeila Ghaffari     d                = (d + 1) / 2;
1450c373b74SJames Wright     PetscInt q_order = d + honee->app_ctx->q_extra;
1460c373b74SJames Wright     if (i + 1 == honee->app_ctx->viz_refine) d = 1;
1472b916ea7SJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1482b916ea7SJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
14967263decSKenneth E. Jansen     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
1502b916ea7SJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
1510c373b74SJames Wright     if (!i) honee->interp_viz = interp_next;
152a515125bSLeila Ghaffari     else {
153a515125bSLeila Ghaffari       Mat C;
1540c373b74SJames Wright       PetscCall(MatMatMult(interp_next, honee->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1552b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1560c373b74SJames Wright       PetscCall(MatDestroy(&honee->interp_viz));
1570c373b74SJames Wright       honee->interp_viz = C;
158a515125bSLeila Ghaffari     }
159a515125bSLeila Ghaffari   }
1600c373b74SJames Wright   for (PetscInt i = 1; i < honee->app_ctx->viz_refine; i++) {
1612b916ea7SJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
162a515125bSLeila Ghaffari   }
1630c373b74SJames Wright   honee->dm_viz = dm_hierarchy[honee->app_ctx->viz_refine];
164d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
165a515125bSLeila Ghaffari }
166