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 11a515125bSLeila Ghaffari #include "../navierstokes.h" 12866f23abSJames Wright #include "../problems/stg_shur14.h" 13a515125bSLeila Ghaffari 1405a512bdSLeila Ghaffari // Create mesh 15*2b916ea7SJeremy L Thompson PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) { 16a515125bSLeila Ghaffari PetscFunctionBeginUser; 1705a512bdSLeila Ghaffari // Create DMPLEX 18*2b916ea7SJeremy L Thompson PetscCall(DMCreate(comm, dm)); 19*2b916ea7SJeremy L Thompson PetscCall(DMSetType(*dm, DMPLEX)); 20b150a244SJed Brown { 21b150a244SJed Brown PetscBool skip = PETSC_TRUE; 22*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL)); 23b150a244SJed Brown PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 24b150a244SJed Brown } 25*2b916ea7SJeremy L Thompson PetscCall(DMSetMatType(*dm, mat_type)); 26*2b916ea7SJeremy L Thompson PetscCall(DMSetVecType(*dm, vec_type)); 27a9804fe9SJed Brown 2805a512bdSLeila Ghaffari // Set Tensor elements 29*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0")); 30*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0")); 3105a512bdSLeila Ghaffari // Set CL options 32*2b916ea7SJeremy L Thompson PetscCall(DMSetFromOptions(*dm)); 33*2b916ea7SJeremy L Thompson PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 34a515125bSLeila Ghaffari PetscFunctionReturn(0); 35a515125bSLeila Ghaffari } 36a515125bSLeila Ghaffari 37a515125bSLeila Ghaffari // Setup DM 38*2b916ea7SJeremy L Thompson PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys) { 39a515125bSLeila Ghaffari PetscFunctionBeginUser; 40a515125bSLeila Ghaffari { 41a515125bSLeila Ghaffari // Configure the finite element space and boundary conditions 42a515125bSLeila Ghaffari PetscFE fe; 43a515125bSLeila Ghaffari PetscInt num_comp_q = 5; 44047c2946SJed Brown DMLabel label; 45*2b916ea7SJeremy L Thompson PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, PETSC_FALSE, degree, PETSC_DECIDE, &fe)); 46*2b916ea7SJeremy L Thompson PetscCall(PetscObjectSetName((PetscObject)fe, "Q")); 47*2b916ea7SJeremy L Thompson PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 48*2b916ea7SJeremy L Thompson PetscCall(DMCreateDS(dm)); 49*2b916ea7SJeremy L Thompson PetscCall(DMGetLabel(dm, "Face Sets", &label)); 50002797a3SLeila Ghaffari // Set wall BCs 51002797a3SLeila Ghaffari if (bc->num_wall > 0) { 52*2b916ea7SJeremy L Thompson PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, 53*2b916ea7SJeremy L Thompson (void (*)(void))problem->bc, NULL, problem->bc_ctx, NULL)); 54002797a3SLeila Ghaffari } 55002797a3SLeila Ghaffari // Set slip BCs in the x direction 56002797a3SLeila Ghaffari if (bc->num_slip[0] > 0) { 57002797a3SLeila Ghaffari PetscInt comps[1] = {1}; 58*2b916ea7SJeremy 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*2b916ea7SJeremy L Thompson problem->bc_ctx, 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}; 64*2b916ea7SJeremy 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*2b916ea7SJeremy L Thompson problem->bc_ctx, NULL)); 66002797a3SLeila Ghaffari } 67002797a3SLeila Ghaffari // Set slip BCs in the z direction 68002797a3SLeila Ghaffari if (bc->num_slip[2] > 0) { 69002797a3SLeila Ghaffari PetscInt comps[1] = {3}; 70*2b916ea7SJeremy 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*2b916ea7SJeremy L Thompson problem->bc_ctx, NULL)); 72002797a3SLeila Ghaffari } 73866f23abSJames Wright { 74866f23abSJames Wright PetscBool use_strongstg = PETSC_FALSE; 75*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 76866f23abSJames Wright 77866f23abSJames Wright if (use_strongstg) { 78*2b916ea7SJeremy L Thompson PetscCall(SetupStrongSTG(dm, bc, problem, phys)); 79866f23abSJames Wright } 80866f23abSJames Wright } 81866f23abSJames Wright 82*2b916ea7SJeremy L Thompson PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 83*2b916ea7SJeremy L Thompson PetscCall(PetscFEDestroy(&fe)); 84a515125bSLeila Ghaffari } 85cbe60e31SLeila Ghaffari 86a515125bSLeila Ghaffari // Empty name for conserved field (because there is only one field) 87a515125bSLeila Ghaffari PetscSection section; 88*2b916ea7SJeremy L Thompson PetscCall(DMGetLocalSection(dm, §ion)); 89*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetFieldName(section, 0, "")); 903636f6a4SJames Wright switch (phys->state_var) { 913636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 92*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 93*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "Momentum X")); 94*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "Momentum Y")); 95*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "Momentum Z")); 96*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "Energy Density")); 973636f6a4SJames Wright break; 983636f6a4SJames Wright 993636f6a4SJames Wright case STATEVAR_PRIMITIVE: 100*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 101*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 1, "Velocity X")); 102*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 2, "Velocity Y")); 103*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 3, "Velocity Z")); 104*2b916ea7SJeremy L Thompson PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 1053636f6a4SJames Wright break; 106a515125bSLeila Ghaffari } 107a515125bSLeila Ghaffari PetscFunctionReturn(0); 108a515125bSLeila Ghaffari } 109a515125bSLeila Ghaffari 110a515125bSLeila Ghaffari // Refine DM for high-order viz 111*2b916ea7SJeremy L Thompson 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 PetscFunctionBeginUser; 115a515125bSLeila Ghaffari 116*2b916ea7SJeremy L Thompson PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 117a515125bSLeila Ghaffari 118a515125bSLeila Ghaffari dm_hierarchy[0] = dm; 119*2b916ea7SJeremy L Thompson for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) { 120a515125bSLeila Ghaffari Mat interp_next; 121*2b916ea7SJeremy L Thompson PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 122*2b916ea7SJeremy L Thompson PetscCall(DMClearDS(dm_hierarchy[i + 1])); 123*2b916ea7SJeremy L Thompson PetscCall(DMClearFields(dm_hierarchy[i + 1])); 124*2b916ea7SJeremy L Thompson PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 125a515125bSLeila Ghaffari d = (d + 1) / 2; 126a515125bSLeila Ghaffari if (i + 1 == user->app_ctx->viz_refine) d = 1; 127*2b916ea7SJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 128*2b916ea7SJeremy L Thompson PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 129*2b916ea7SJeremy L Thompson PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, bc, phys)); 130*2b916ea7SJeremy L Thompson PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 131a515125bSLeila Ghaffari if (!i) user->interp_viz = interp_next; 132a515125bSLeila Ghaffari else { 133a515125bSLeila Ghaffari Mat C; 134*2b916ea7SJeremy L Thompson PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 135*2b916ea7SJeremy L Thompson PetscCall(MatDestroy(&interp_next)); 136*2b916ea7SJeremy L Thompson PetscCall(MatDestroy(&user->interp_viz)); 137a515125bSLeila Ghaffari user->interp_viz = C; 138a515125bSLeila Ghaffari } 139a515125bSLeila Ghaffari } 140a515125bSLeila Ghaffari for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) { 141*2b916ea7SJeremy L Thompson PetscCall(DMDestroy(&dm_hierarchy[i])); 142a515125bSLeila Ghaffari } 143a515125bSLeila Ghaffari user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 144a515125bSLeila Ghaffari 145a515125bSLeila Ghaffari PetscFunctionReturn(0); 146a515125bSLeila Ghaffari } 147