xref: /honee/src/setupdm.c (revision dc936754fc0ae21fa21dd641901ba00fc24c3769)
1*dc936754SJeremy 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 
15a515125bSLeila Ghaffari #include "../navierstokes.h"
16866f23abSJames Wright #include "../problems/stg_shur14.h"
17a515125bSLeila Ghaffari 
1805a512bdSLeila Ghaffari // Create mesh
192b916ea7SJeremy L Thompson 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"));
35ebfabadfSJames Wright   PetscCall(PetscOptionsSetValue(NULL, "-dm_blocking_type", "field_node"));
36ebfabadfSJames Wright 
3705a512bdSLeila Ghaffari   // Set CL options
382b916ea7SJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
392b916ea7SJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
40d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
41a515125bSLeila Ghaffari }
42a515125bSLeila Ghaffari 
43a515125bSLeila Ghaffari // Setup DM
4467263decSKenneth E. Jansen PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
45a515125bSLeila Ghaffari   PetscInt num_comp_q = 5;
46da4ca0cfSJames Wright   PetscFunctionBeginUser;
47da4ca0cfSJames Wright 
48b2e5b5b3SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
49da4ca0cfSJames Wright 
50da4ca0cfSJames Wright   {  // Add strong boundary conditions to DM
51047c2946SJed Brown     DMLabel label;
522b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
53c5e9980aSAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
54002797a3SLeila Ghaffari     // Set wall BCs
55002797a3SLeila Ghaffari     if (bc->num_wall > 0) {
56bedfd28dSJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, NULL, NULL, NULL, NULL));
57002797a3SLeila Ghaffari     }
58fce2147eSJames Wright     // Set symmetry BCs in the x direction
59fce2147eSJames Wright     if (bc->num_symmetry[0] > 0) {
60002797a3SLeila Ghaffari       PetscInt comps[1] = {1};
61fce2147eSJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_x", label, bc->num_symmetry[0], bc->symmetries[0], 0, 1, comps, NULL, NULL, NULL, NULL));
62002797a3SLeila Ghaffari     }
63fce2147eSJames Wright     // Set symmetry BCs in the y direction
64fce2147eSJames Wright     if (bc->num_symmetry[1] > 0) {
65002797a3SLeila Ghaffari       PetscInt comps[1] = {2};
66fce2147eSJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_y", label, bc->num_symmetry[1], bc->symmetries[1], 0, 1, comps, NULL, NULL, NULL, NULL));
67002797a3SLeila Ghaffari     }
68fce2147eSJames Wright     // Set symmetry BCs in the z direction
69fce2147eSJames Wright     if (bc->num_symmetry[2] > 0) {
70002797a3SLeila Ghaffari       PetscInt comps[1] = {3};
71fce2147eSJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_z", label, bc->num_symmetry[2], bc->symmetries[2], 0, 1, comps, NULL, NULL, NULL, NULL));
72002797a3SLeila Ghaffari     }
73866f23abSJames Wright     {
74866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
752b916ea7SJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
7642454adaSJames Wright       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
77866f23abSJames Wright     }
7867263decSKenneth E. Jansen   }
7967263decSKenneth E. Jansen 
80da4ca0cfSJames Wright   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
81cbe60e31SLeila Ghaffari 
82a515125bSLeila Ghaffari   // Empty name for conserved field (because there is only one field)
83a515125bSLeila Ghaffari   PetscSection section;
842b916ea7SJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
852b916ea7SJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
863636f6a4SJames Wright   switch (phys->state_var) {
873636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
882b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
892b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
902b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
912b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
92525869caSJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
933636f6a4SJames Wright       break;
943636f6a4SJames Wright 
953636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
962b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
972b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
982b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
992b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
1002b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
1013636f6a4SJames Wright       break;
102a515125bSLeila Ghaffari   }
103d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
104a515125bSLeila Ghaffari }
105a515125bSLeila Ghaffari 
106a515125bSLeila Ghaffari // Refine DM for high-order viz
1072b916ea7SJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
108a515125bSLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
109a515125bSLeila Ghaffari   VecType vec_type;
110a515125bSLeila Ghaffari 
11106f41313SJames Wright   PetscFunctionBeginUser;
1122b916ea7SJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
113a515125bSLeila Ghaffari 
114a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
1152b916ea7SJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
116a515125bSLeila Ghaffari     Mat interp_next;
1172b916ea7SJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1182b916ea7SJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1192b916ea7SJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1202b916ea7SJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
121a515125bSLeila Ghaffari     d                = (d + 1) / 2;
12267263decSKenneth E. Jansen     PetscInt q_order = d + user->app_ctx->q_extra;
123a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
1242b916ea7SJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1252b916ea7SJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
12667263decSKenneth E. Jansen     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
1272b916ea7SJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
128a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
129a515125bSLeila Ghaffari     else {
130a515125bSLeila Ghaffari       Mat C;
1312b916ea7SJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1322b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1332b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
134a515125bSLeila Ghaffari       user->interp_viz = C;
135a515125bSLeila Ghaffari     }
136a515125bSLeila Ghaffari   }
137a515125bSLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
1382b916ea7SJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
139a515125bSLeila Ghaffari   }
140a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
141d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
142a515125bSLeila Ghaffari }
143