xref: /honee/src/setupdm.c (revision b2e5b5b3cbbe5d53fc14e38154896a6dce61003c)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, 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"));
3505a512bdSLeila Ghaffari   // Set CL options
362b916ea7SJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
372b916ea7SJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
38d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
39a515125bSLeila Ghaffari }
40a515125bSLeila Ghaffari 
41a515125bSLeila Ghaffari // Setup DM
4267263decSKenneth E. Jansen PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) {
43a515125bSLeila Ghaffari   PetscInt num_comp_q = 5;
44da4ca0cfSJames Wright   PetscFunctionBeginUser;
45da4ca0cfSJames Wright 
46*b2e5b5b3SJames Wright   PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm));
47da4ca0cfSJames Wright 
48da4ca0cfSJames Wright   {  // Add strong boundary conditions to DM
49047c2946SJed Brown     DMLabel label;
502b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
51c5e9980aSAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
52002797a3SLeila Ghaffari     // Set wall BCs
53002797a3SLeila Ghaffari     if (bc->num_wall > 0) {
54bedfd28dSJames 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));
55002797a3SLeila Ghaffari     }
56002797a3SLeila Ghaffari     // Set slip BCs in the x direction
57002797a3SLeila Ghaffari     if (bc->num_slip[0] > 0) {
58002797a3SLeila Ghaffari       PetscInt comps[1] = {1};
59bedfd28dSJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, NULL, NULL, NULL, NULL));
60002797a3SLeila Ghaffari     }
61002797a3SLeila Ghaffari     // Set slip BCs in the y direction
62002797a3SLeila Ghaffari     if (bc->num_slip[1] > 0) {
63002797a3SLeila Ghaffari       PetscInt comps[1] = {2};
64bedfd28dSJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, NULL, NULL, NULL, NULL));
65002797a3SLeila Ghaffari     }
66002797a3SLeila Ghaffari     // Set slip BCs in the z direction
67002797a3SLeila Ghaffari     if (bc->num_slip[2] > 0) {
68002797a3SLeila Ghaffari       PetscInt comps[1] = {3};
69bedfd28dSJames Wright       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, NULL, NULL, NULL, NULL));
70002797a3SLeila Ghaffari     }
71866f23abSJames Wright     {
72866f23abSJames Wright       PetscBool use_strongstg = PETSC_FALSE;
732b916ea7SJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
7442454adaSJames Wright       if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys));
75866f23abSJames Wright     }
7667263decSKenneth E. Jansen   }
7767263decSKenneth E. Jansen 
78da4ca0cfSJames Wright   PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm));
79cbe60e31SLeila Ghaffari 
80a515125bSLeila Ghaffari   // Empty name for conserved field (because there is only one field)
81a515125bSLeila Ghaffari   PetscSection section;
822b916ea7SJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
832b916ea7SJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
843636f6a4SJames Wright   switch (phys->state_var) {
853636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
862b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
872b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX"));
882b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY"));
892b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ"));
90525869caSJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy"));
913636f6a4SJames Wright       break;
923636f6a4SJames Wright 
933636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
942b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
952b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
962b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
972b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
982b916ea7SJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
993636f6a4SJames Wright       break;
100a515125bSLeila Ghaffari   }
101d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
102a515125bSLeila Ghaffari }
103a515125bSLeila Ghaffari 
104a515125bSLeila Ghaffari // Refine DM for high-order viz
1052b916ea7SJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
106a515125bSLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
107a515125bSLeila Ghaffari   VecType vec_type;
108a515125bSLeila Ghaffari 
10906f41313SJames Wright   PetscFunctionBeginUser;
1102b916ea7SJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
111a515125bSLeila Ghaffari 
112a515125bSLeila Ghaffari   dm_hierarchy[0] = dm;
1132b916ea7SJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
114a515125bSLeila Ghaffari     Mat interp_next;
1152b916ea7SJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1162b916ea7SJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1172b916ea7SJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1182b916ea7SJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
119a515125bSLeila Ghaffari     d                = (d + 1) / 2;
12067263decSKenneth E. Jansen     PetscInt q_order = d + user->app_ctx->q_extra;
121a515125bSLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
1222b916ea7SJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1232b916ea7SJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
12467263decSKenneth E. Jansen     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys));
1252b916ea7SJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
126a515125bSLeila Ghaffari     if (!i) user->interp_viz = interp_next;
127a515125bSLeila Ghaffari     else {
128a515125bSLeila Ghaffari       Mat C;
1292b916ea7SJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1302b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1312b916ea7SJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
132a515125bSLeila Ghaffari       user->interp_viz = C;
133a515125bSLeila Ghaffari     }
134a515125bSLeila Ghaffari   }
135a515125bSLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
1362b916ea7SJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
137a515125bSLeila Ghaffari   }
138a515125bSLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
139d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
140a515125bSLeila Ghaffari }
141