xref: /honee/src/setupdm.c (revision 149fb5361f5198e41f87e8815a7e9dbfee84a96a)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdmplex.h>
1367263decSKenneth E. Jansen #include <petscds.h>
14e419654dSJeremy L Thompson 
15*149fb536SJames Wright #include <navierstokes.h>
16866f23abSJames Wright #include "../problems/stg_shur14.h"
17a515125bSLeila Ghaffari 
1805a512bdSLeila Ghaffari // Create mesh
19991aef52SJames Wright PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType mat_type, VecType vec_type, DM *dm) {
20a515125bSLeila Ghaffari   PetscFunctionBeginUser;
2105a512bdSLeila Ghaffari   // Create DMPLEX
222b916ea7SJeremy L Thompson   PetscCall(DMCreate(comm, dm));
232b916ea7SJeremy L Thompson   PetscCall(DMSetType(*dm, DMPLEX));
24b150a244SJed Brown   {
25b150a244SJed Brown     PetscBool skip = PETSC_TRUE;
262b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
27b150a244SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
28b150a244SJed Brown   }
292b916ea7SJeremy L Thompson   PetscCall(DMSetMatType(*dm, mat_type));
302b916ea7SJeremy L Thompson   PetscCall(DMSetVecType(*dm, vec_type));
31a9804fe9SJed Brown 
3205a512bdSLeila Ghaffari   // Set Tensor elements
332b916ea7SJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"));
342b916ea7SJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"));
359bfedf88SJames Wright   PetscCall(PetscOptionsSetValue(NULL, "-dm_localize", "0"));  // Localization done in DMSetupByOrderEnd_FEM
36ebfabadfSJames Wright   PetscCall(PetscOptionsSetValue(NULL, "-dm_blocking_type", "field_node"));
37ebfabadfSJames Wright 
3805a512bdSLeila Ghaffari   // Set CL options
392b916ea7SJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
402b916ea7SJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
41d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
42a515125bSLeila Ghaffari }
43a515125bSLeila Ghaffari 
44a515125bSLeila Ghaffari // Setup DM
45991aef52SJames Wright PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
46a515125bSLeila Ghaffari   PetscInt num_comp_q = 5;
47da4ca0cfSJames Wright   PetscFunctionBeginUser;
48da4ca0cfSJames Wright 
49b2e5b5b3SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
50da4ca0cfSJames Wright 
51da4ca0cfSJames Wright   {  // Add strong boundary conditions to DM
52047c2946SJed Brown     DMLabel label;
532b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
54c5e9980aSAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
55487a3b6eSJames Wright 
56487a3b6eSJames Wright     for (PetscInt i = 0; i < problem->num_bc_defs; i++) {
57487a3b6eSJames Wright       BCDefinition    bc_def = problem->bc_defs[i];
58487a3b6eSJames Wright       PetscInt        num_essential_comps, num_label_values;
59487a3b6eSJames Wright       const PetscInt *essential_comps, *label_values;
60487a3b6eSJames Wright       const char     *name;
61487a3b6eSJames Wright 
62487a3b6eSJames Wright       PetscCall(BCDefinitionGetEssential(bc_def, &num_essential_comps, &essential_comps));
63487a3b6eSJames Wright       if (essential_comps > 0) {
64487a3b6eSJames Wright         PetscCall(BCDefinitionGetInfo(bc_def, &name, &num_label_values, &label_values));
65487a3b6eSJames Wright         PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, name, label, num_label_values, label_values, 0, num_essential_comps, essential_comps, NULL, NULL,
66487a3b6eSJames Wright                                 NULL, NULL));
67002797a3SLeila Ghaffari       }
68002797a3SLeila Ghaffari     }
69866f23abSJames Wright     {
70866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
712b916ea7SJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
7242454adaSJames Wright       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
73866f23abSJames Wright     }
7467263decSKenneth E. Jansen   }
7567263decSKenneth E. Jansen 
76da4ca0cfSJames Wright   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
77cbe60e31SLeila Ghaffari 
78a515125bSLeila Ghaffari   // Empty name for conserved field (because there is only one field)
79a515125bSLeila Ghaffari   PetscSection section;
802b916ea7SJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
812b916ea7SJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
823636f6a4SJames Wright   switch (phys->state_var) {
833636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
842b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
852b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
862b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
872b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
88525869caSJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
893636f6a4SJames Wright       break;
903636f6a4SJames Wright 
913636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
922b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
932b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
942b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
952b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
962b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
973636f6a4SJames Wright       break;
989b103f75SJames Wright 
999b103f75SJames Wright     case STATEVAR_ENTROPY:
1009b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 0, "EntropyDensity"));
1019b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 1, "EntropyMomentumX"));
1029b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 2, "EntropyMomentumY"));
1039b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 3, "EntropyMomentumZ"));
1049b103f75SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "EntropyTotalEnergy"));
1059b103f75SJames Wright       break;
106a515125bSLeila Ghaffari   }
107d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
108a515125bSLeila Ghaffari }
109a515125bSLeila Ghaffari 
110a515125bSLeila Ghaffari // Refine DM for high-order viz
111991aef52SJames Wright PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys) {
112a515125bSLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
113a515125bSLeila Ghaffari   VecType vec_type;
114a515125bSLeila Ghaffari 
11506f41313SJames Wright   PetscFunctionBeginUser;
1162b916ea7SJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
117a515125bSLeila Ghaffari 
118a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
1192b916ea7SJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
120a515125bSLeila Ghaffari     Mat interp_next;
1212b916ea7SJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1222b916ea7SJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1232b916ea7SJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1242b916ea7SJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
125a515125bSLeila Ghaffari     d                = (d + 1) / 2;
12667263decSKenneth E. Jansen     PetscInt q_order = d + user->app_ctx->q_extra;
127a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
1282b916ea7SJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1292b916ea7SJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
13067263decSKenneth E. Jansen     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
1312b916ea7SJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
132a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
133a515125bSLeila Ghaffari     else {
134a515125bSLeila Ghaffari       Mat C;
1352b916ea7SJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1362b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1372b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
138a515125bSLeila Ghaffari       user->interp_viz = C;
139a515125bSLeila Ghaffari     }
140a515125bSLeila Ghaffari   }
141a515125bSLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
1422b916ea7SJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
143a515125bSLeila Ghaffari   }
144a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
145d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
146a515125bSLeila Ghaffari }
147