1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Setup DM for Navier-Stokes example using PETSc 10 11 #include "../navierstokes.h" 12 #include "../problems/stg_shur14.h" 13 14 // Create mesh 15 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) { 16 PetscFunctionBeginUser; 17 // Create DMPLEX 18 PetscCall(DMCreate(comm, dm)); 19 PetscCall(DMSetType(*dm, DMPLEX)); 20 { 21 PetscBool skip = PETSC_TRUE; 22 PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL)); 23 PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 24 } 25 PetscCall(DMSetMatType(*dm, mat_type)); 26 PetscCall(DMSetVecType(*dm, vec_type)); 27 28 // Set Tensor elements 29 PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0")); 30 PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0")); 31 // Set CL options 32 PetscCall(DMSetFromOptions(*dm)); 33 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 34 PetscFunctionReturn(0); 35 } 36 37 // Setup DM 38 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys) { 39 PetscFunctionBeginUser; 40 { 41 // Configure the finite element space and boundary conditions 42 PetscFE fe; 43 PetscInt num_comp_q = 5; 44 DMLabel label; 45 PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, PETSC_FALSE, degree, PETSC_DECIDE, &fe)); 46 PetscCall(PetscObjectSetName((PetscObject)fe, "Q")); 47 PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 48 PetscCall(DMCreateDS(dm)); 49 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 50 PetscCall(DMPlexLabelComplete(dm, label)); 51 // Set wall BCs 52 if (bc->num_wall > 0) { 53 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, 54 (void (*)(void))problem->bc, NULL, problem->bc_ctx, NULL)); 55 } 56 // Set slip BCs in the x direction 57 if (bc->num_slip[0] > 0) { 58 PetscInt comps[1] = {1}; 59 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, bc->num_slip[0], bc->slips[0], 0, 1, comps, (void (*)(void))NULL, NULL, 60 problem->bc_ctx, NULL)); 61 } 62 // Set slip BCs in the y direction 63 if (bc->num_slip[1] > 0) { 64 PetscInt comps[1] = {2}; 65 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, bc->num_slip[1], bc->slips[1], 0, 1, comps, (void (*)(void))NULL, NULL, 66 problem->bc_ctx, NULL)); 67 } 68 // Set slip BCs in the z direction 69 if (bc->num_slip[2] > 0) { 70 PetscInt comps[1] = {3}; 71 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, bc->num_slip[2], bc->slips[2], 0, 1, comps, (void (*)(void))NULL, NULL, 72 problem->bc_ctx, NULL)); 73 } 74 { 75 PetscBool use_strongstg = PETSC_FALSE; 76 PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 77 78 if (use_strongstg) { 79 PetscCall(SetupStrongSTG(dm, bc, problem, phys)); 80 } 81 } 82 83 PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 84 PetscCall(PetscFEDestroy(&fe)); 85 } 86 87 // Empty name for conserved field (because there is only one field) 88 PetscSection section; 89 PetscCall(DMGetLocalSection(dm, §ion)); 90 PetscCall(PetscSectionSetFieldName(section, 0, "")); 91 switch (phys->state_var) { 92 case STATEVAR_CONSERVATIVE: 93 PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 94 PetscCall(PetscSectionSetComponentName(section, 0, 1, "Momentum X")); 95 PetscCall(PetscSectionSetComponentName(section, 0, 2, "Momentum Y")); 96 PetscCall(PetscSectionSetComponentName(section, 0, 3, "Momentum Z")); 97 PetscCall(PetscSectionSetComponentName(section, 0, 4, "Energy Density")); 98 break; 99 100 case STATEVAR_PRIMITIVE: 101 PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 102 PetscCall(PetscSectionSetComponentName(section, 0, 1, "Velocity X")); 103 PetscCall(PetscSectionSetComponentName(section, 0, 2, "Velocity Y")); 104 PetscCall(PetscSectionSetComponentName(section, 0, 3, "Velocity Z")); 105 PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 106 break; 107 } 108 PetscFunctionReturn(0); 109 } 110 111 // Refine DM for high-order viz 112 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) { 113 DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 114 VecType vec_type; 115 PetscFunctionBeginUser; 116 117 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 118 119 dm_hierarchy[0] = dm; 120 for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) { 121 Mat interp_next; 122 PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 123 PetscCall(DMClearDS(dm_hierarchy[i + 1])); 124 PetscCall(DMClearFields(dm_hierarchy[i + 1])); 125 PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 126 d = (d + 1) / 2; 127 if (i + 1 == user->app_ctx->viz_refine) d = 1; 128 PetscCall(DMGetVecType(dm, &vec_type)); 129 PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 130 PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, bc, phys)); 131 PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 132 if (!i) user->interp_viz = interp_next; 133 else { 134 Mat C; 135 PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 136 PetscCall(MatDestroy(&interp_next)); 137 PetscCall(MatDestroy(&user->interp_viz)); 138 user->interp_viz = C; 139 } 140 } 141 for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) { 142 PetscCall(DMDestroy(&dm_hierarchy[i])); 143 } 144 user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 145 146 PetscFunctionReturn(0); 147 } 148