xref: /libCEED/examples/fluids/src/setupdm.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Setup DM for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #include "../navierstokes.h"
12363b60e1SJames Wright #include "../problems/stg_shur14.h"
1377841947SLeila Ghaffari 
141864f1c2SLeila Ghaffari // Create mesh
15*2b730f8bSJeremy L Thompson PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) {
1677841947SLeila Ghaffari   PetscFunctionBeginUser;
171864f1c2SLeila Ghaffari   // Create DMPLEX
18*2b730f8bSJeremy L Thompson   PetscCall(DMCreate(comm, dm));
19*2b730f8bSJeremy L Thompson   PetscCall(DMSetType(*dm, DMPLEX));
20c32b0260SJed Brown   {
21c32b0260SJed Brown     PetscBool skip = PETSC_TRUE;
22*2b730f8bSJeremy L Thompson     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
23c32b0260SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
24c32b0260SJed Brown   }
25*2b730f8bSJeremy L Thompson   PetscCall(DMSetMatType(*dm, mat_type));
26*2b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(*dm, vec_type));
274ea65e7bSJed Brown 
281864f1c2SLeila Ghaffari   // Set Tensor elements
29*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"));
30*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"));
311864f1c2SLeila Ghaffari   // Set CL options
32*2b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
33*2b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3477841947SLeila Ghaffari   PetscFunctionReturn(0);
3577841947SLeila Ghaffari }
3677841947SLeila Ghaffari 
3777841947SLeila Ghaffari // Setup DM
38*2b730f8bSJeremy L Thompson PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys) {
3977841947SLeila Ghaffari   PetscFunctionBeginUser;
4077841947SLeila Ghaffari   {
4177841947SLeila Ghaffari     // Configure the finite element space and boundary conditions
4277841947SLeila Ghaffari     PetscFE  fe;
4377841947SLeila Ghaffari     PetscInt num_comp_q = 5;
44496f2382SJed Brown     DMLabel  label;
45*2b730f8bSJeremy L Thompson     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
46*2b730f8bSJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)fe, "Q"));
47*2b730f8bSJeremy L Thompson     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
48*2b730f8bSJeremy L Thompson     PetscCall(DMCreateDS(dm));
49*2b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
502fe7aee7SLeila Ghaffari     // Set wall BCs
512fe7aee7SLeila Ghaffari     if (bc->num_wall > 0) {
52*2b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps,
53*2b730f8bSJeremy L Thompson                               (void (*)(void))problem->bc, NULL, problem->bc_ctx, NULL));
542fe7aee7SLeila Ghaffari     }
552fe7aee7SLeila Ghaffari     // Set slip BCs in the x direction
562fe7aee7SLeila Ghaffari     if (bc->num_slip[0] > 0) {
572fe7aee7SLeila Ghaffari       PetscInt comps[1] = {1};
58*2b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, (void (*)(void))NULL, NULL,
59*2b730f8bSJeremy L Thompson                               problem->bc_ctx, NULL));
602fe7aee7SLeila Ghaffari     }
612fe7aee7SLeila Ghaffari     // Set slip BCs in the y direction
622fe7aee7SLeila Ghaffari     if (bc->num_slip[1] > 0) {
632fe7aee7SLeila Ghaffari       PetscInt comps[1] = {2};
64*2b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, (void (*)(void))NULL, NULL,
65*2b730f8bSJeremy L Thompson                               problem->bc_ctx, NULL));
662fe7aee7SLeila Ghaffari     }
672fe7aee7SLeila Ghaffari     // Set slip BCs in the z direction
682fe7aee7SLeila Ghaffari     if (bc->num_slip[2] > 0) {
692fe7aee7SLeila Ghaffari       PetscInt comps[1] = {3};
70*2b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, (void (*)(void))NULL, NULL,
71*2b730f8bSJeremy L Thompson                               problem->bc_ctx, NULL));
722fe7aee7SLeila Ghaffari     }
73363b60e1SJames Wright     {
74363b60e1SJames Wright       PetscBool use_strongstg = PETSC_FALSE;
75*2b730f8bSJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
76363b60e1SJames Wright 
77363b60e1SJames Wright       if (use_strongstg) {
78*2b730f8bSJeremy L Thompson         PetscCall(SetupStrongSTG(dm, bc, problem, phys));
79363b60e1SJames Wright       }
80363b60e1SJames Wright     }
81363b60e1SJames Wright 
82*2b730f8bSJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
83*2b730f8bSJeremy L Thompson     PetscCall(PetscFEDestroy(&fe));
8477841947SLeila Ghaffari   }
85dc805cc4SLeila Ghaffari 
8677841947SLeila Ghaffari   // Empty name for conserved field (because there is only one field)
8777841947SLeila Ghaffari   PetscSection section;
88*2b730f8bSJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
89*2b730f8bSJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
9097baf651SJames Wright   switch (phys->state_var) {
9197baf651SJames Wright     case STATEVAR_CONSERVATIVE:
92*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
93*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Momentum X"));
94*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Momentum Y"));
95*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Momentum Z"));
96*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Energy Density"));
9797baf651SJames Wright       break;
9897baf651SJames Wright 
9997baf651SJames Wright     case STATEVAR_PRIMITIVE:
100*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
101*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Velocity X"));
102*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Velocity Y"));
103*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Velocity Z"));
104*2b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
10597baf651SJames Wright       break;
10677841947SLeila Ghaffari   }
10777841947SLeila Ghaffari   PetscFunctionReturn(0);
10877841947SLeila Ghaffari }
10977841947SLeila Ghaffari 
11077841947SLeila Ghaffari // Refine DM for high-order viz
111*2b730f8bSJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
11277841947SLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
11377841947SLeila Ghaffari   VecType vec_type;
11477841947SLeila Ghaffari   PetscFunctionBeginUser;
11577841947SLeila Ghaffari 
116*2b730f8bSJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
11777841947SLeila Ghaffari 
11877841947SLeila Ghaffari   dm_hierarchy[0] = dm;
119*2b730f8bSJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
12077841947SLeila Ghaffari     Mat interp_next;
121*2b730f8bSJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
122*2b730f8bSJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
123*2b730f8bSJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
124*2b730f8bSJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
12577841947SLeila Ghaffari     d = (d + 1) / 2;
12677841947SLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
127*2b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
128*2b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
129*2b730f8bSJeremy L Thompson     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, bc, phys));
130*2b730f8bSJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
13177841947SLeila Ghaffari     if (!i) user->interp_viz = interp_next;
13277841947SLeila Ghaffari     else {
13377841947SLeila Ghaffari       Mat C;
134*2b730f8bSJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
135*2b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
136*2b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
13777841947SLeila Ghaffari       user->interp_viz = C;
13877841947SLeila Ghaffari     }
13977841947SLeila Ghaffari   }
14077841947SLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
141*2b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
14277841947SLeila Ghaffari   }
14377841947SLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
14477841947SLeila Ghaffari 
14577841947SLeila Ghaffari   PetscFunctionReturn(0);
14677841947SLeila Ghaffari }
147