xref: /libCEED/examples/fluids/src/setupdm.c (revision 49aac155e7a09736f56fb3abac0f57dab29f7cbf)
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 
11*49aac155SJeremy L Thompson #include <ceed.h>
12*49aac155SJeremy L Thompson #include <petscdmplex.h>
13*49aac155SJeremy L Thompson 
1477841947SLeila Ghaffari #include "../navierstokes.h"
15363b60e1SJames Wright #include "../problems/stg_shur14.h"
1677841947SLeila Ghaffari 
171864f1c2SLeila Ghaffari // Create mesh
182b730f8bSJeremy L Thompson PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) {
1977841947SLeila Ghaffari   PetscFunctionBeginUser;
201864f1c2SLeila Ghaffari   // Create DMPLEX
212b730f8bSJeremy L Thompson   PetscCall(DMCreate(comm, dm));
222b730f8bSJeremy L Thompson   PetscCall(DMSetType(*dm, DMPLEX));
23c32b0260SJed Brown   {
24c32b0260SJed Brown     PetscBool skip = PETSC_TRUE;
252b730f8bSJeremy L Thompson     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL));
26c32b0260SJed Brown     PetscCall(DMSetMatrixPreallocateSkip(*dm, skip));
27c32b0260SJed Brown   }
282b730f8bSJeremy L Thompson   PetscCall(DMSetMatType(*dm, mat_type));
292b730f8bSJeremy L Thompson   PetscCall(DMSetVecType(*dm, vec_type));
304ea65e7bSJed Brown 
311864f1c2SLeila Ghaffari   // Set Tensor elements
322b730f8bSJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"));
332b730f8bSJeremy L Thompson   PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"));
341864f1c2SLeila Ghaffari   // Set CL options
352b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(*dm));
362b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3777841947SLeila Ghaffari   PetscFunctionReturn(0);
3877841947SLeila Ghaffari }
3977841947SLeila Ghaffari 
4077841947SLeila Ghaffari // Setup DM
412b730f8bSJeremy L Thompson PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys) {
4277841947SLeila Ghaffari   PetscFunctionBeginUser;
4377841947SLeila Ghaffari   {
4477841947SLeila Ghaffari     // Configure the finite element space and boundary conditions
4577841947SLeila Ghaffari     PetscFE  fe;
4677841947SLeila Ghaffari     PetscInt num_comp_q = 5;
47496f2382SJed Brown     DMLabel  label;
482b730f8bSJeremy L Thompson     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
492b730f8bSJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)fe, "Q"));
502b730f8bSJeremy L Thompson     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
512b730f8bSJeremy L Thompson     PetscCall(DMCreateDS(dm));
522b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &label));
53ca69d878SAdeleke O. Bankole     PetscCall(DMPlexLabelComplete(dm, label));
542fe7aee7SLeila Ghaffari     // Set wall BCs
552fe7aee7SLeila Ghaffari     if (bc->num_wall > 0) {
562b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps,
572b730f8bSJeremy L Thompson                               (void (*)(void))problem->bc, NULL, problem->bc_ctx, NULL));
582fe7aee7SLeila Ghaffari     }
592fe7aee7SLeila Ghaffari     // Set slip BCs in the x direction
602fe7aee7SLeila Ghaffari     if (bc->num_slip[0] > 0) {
612fe7aee7SLeila Ghaffari       PetscInt comps[1] = {1};
622b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, (void (*)(void))NULL, NULL,
632b730f8bSJeremy L Thompson                               problem->bc_ctx, NULL));
642fe7aee7SLeila Ghaffari     }
652fe7aee7SLeila Ghaffari     // Set slip BCs in the y direction
662fe7aee7SLeila Ghaffari     if (bc->num_slip[1] > 0) {
672fe7aee7SLeila Ghaffari       PetscInt comps[1] = {2};
682b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, (void (*)(void))NULL, NULL,
692b730f8bSJeremy L Thompson                               problem->bc_ctx, NULL));
702fe7aee7SLeila Ghaffari     }
712fe7aee7SLeila Ghaffari     // Set slip BCs in the z direction
722fe7aee7SLeila Ghaffari     if (bc->num_slip[2] > 0) {
732fe7aee7SLeila Ghaffari       PetscInt comps[1] = {3};
742b730f8bSJeremy L Thompson       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, (void (*)(void))NULL, NULL,
752b730f8bSJeremy L Thompson                               problem->bc_ctx, NULL));
762fe7aee7SLeila Ghaffari     }
77363b60e1SJames Wright     {
78363b60e1SJames Wright       PetscBool use_strongstg = PETSC_FALSE;
792b730f8bSJeremy L Thompson       PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
80363b60e1SJames Wright 
81363b60e1SJames Wright       if (use_strongstg) {
822b730f8bSJeremy L Thompson         PetscCall(SetupStrongSTG(dm, bc, problem, phys));
83363b60e1SJames Wright       }
84363b60e1SJames Wright     }
85363b60e1SJames Wright 
862b730f8bSJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
872b730f8bSJeremy L Thompson     PetscCall(PetscFEDestroy(&fe));
8877841947SLeila Ghaffari   }
89dc805cc4SLeila Ghaffari 
9077841947SLeila Ghaffari   // Empty name for conserved field (because there is only one field)
9177841947SLeila Ghaffari   PetscSection section;
922b730f8bSJeremy L Thompson   PetscCall(DMGetLocalSection(dm, &section));
932b730f8bSJeremy L Thompson   PetscCall(PetscSectionSetFieldName(section, 0, ""));
9497baf651SJames Wright   switch (phys->state_var) {
9597baf651SJames Wright     case STATEVAR_CONSERVATIVE:
962b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density"));
972b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Momentum X"));
982b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Momentum Y"));
992b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Momentum Z"));
1002b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Energy Density"));
10197baf651SJames Wright       break;
10297baf651SJames Wright 
10397baf651SJames Wright     case STATEVAR_PRIMITIVE:
1042b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
1052b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 1, "Velocity X"));
1062b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 2, "Velocity Y"));
1072b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 3, "Velocity Z"));
1082b730f8bSJeremy L Thompson       PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
10997baf651SJames Wright       break;
11077841947SLeila Ghaffari   }
11177841947SLeila Ghaffari   PetscFunctionReturn(0);
11277841947SLeila Ghaffari }
11377841947SLeila Ghaffari 
11477841947SLeila Ghaffari // Refine DM for high-order viz
1152b730f8bSJeremy L Thompson PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) {
11677841947SLeila Ghaffari   DM      dm_hierarchy[user->app_ctx->viz_refine + 1];
11777841947SLeila Ghaffari   VecType vec_type;
11877841947SLeila Ghaffari   PetscFunctionBeginUser;
11977841947SLeila Ghaffari 
1202b730f8bSJeremy L Thompson   PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
12177841947SLeila Ghaffari 
12277841947SLeila Ghaffari   dm_hierarchy[0] = dm;
1232b730f8bSJeremy L Thompson   for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) {
12477841947SLeila Ghaffari     Mat interp_next;
1252b730f8bSJeremy L Thompson     PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1]));
1262b730f8bSJeremy L Thompson     PetscCall(DMClearDS(dm_hierarchy[i + 1]));
1272b730f8bSJeremy L Thompson     PetscCall(DMClearFields(dm_hierarchy[i + 1]));
1282b730f8bSJeremy L Thompson     PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i]));
12977841947SLeila Ghaffari     d = (d + 1) / 2;
13077841947SLeila Ghaffari     if (i + 1 == user->app_ctx->viz_refine) d = 1;
1312b730f8bSJeremy L Thompson     PetscCall(DMGetVecType(dm, &vec_type));
1322b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type));
1332b730f8bSJeremy L Thompson     PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, bc, phys));
1342b730f8bSJeremy L Thompson     PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL));
13577841947SLeila Ghaffari     if (!i) user->interp_viz = interp_next;
13677841947SLeila Ghaffari     else {
13777841947SLeila Ghaffari       Mat C;
1382b730f8bSJeremy L Thompson       PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C));
1392b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&interp_next));
1402b730f8bSJeremy L Thompson       PetscCall(MatDestroy(&user->interp_viz));
14177841947SLeila Ghaffari       user->interp_viz = C;
14277841947SLeila Ghaffari     }
14377841947SLeila Ghaffari   }
14477841947SLeila Ghaffari   for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) {
1452b730f8bSJeremy L Thompson     PetscCall(DMDestroy(&dm_hierarchy[i]));
14677841947SLeila Ghaffari   }
14777841947SLeila Ghaffari   user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine];
14877841947SLeila Ghaffari 
14977841947SLeila Ghaffari   PetscFunctionReturn(0);
15077841947SLeila Ghaffari }
151