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, 16 MatType mat_type, VecType vec_type, 17 DM *dm) { 18 PetscErrorCode ierr; 19 PetscFunctionBeginUser; 20 // Create DMPLEX 21 ierr = DMCreate(comm, dm); CHKERRQ(ierr); 22 ierr = DMSetType(*dm, DMPLEX); CHKERRQ(ierr); 23 ierr = DMSetMatType(*dm, mat_type); CHKERRQ(ierr); 24 ierr = DMSetVecType(*dm, vec_type); CHKERRQ(ierr); 25 26 // Set Tensor elements 27 ierr = PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0"); CHKERRQ(ierr); 28 ierr = PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0"); CHKERRQ(ierr); 29 // Set CL options 30 ierr = DMSetFromOptions(*dm); CHKERRQ(ierr); 31 ierr = DMViewFromOptions(*dm, NULL, "-dm_view"); CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 // Setup DM 36 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, 37 SimpleBC bc, Physics phys) { 38 PetscErrorCode ierr; 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 ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, num_comp_q, 46 PETSC_FALSE, degree, PETSC_DECIDE, 47 &fe); CHKERRQ(ierr); 48 ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 49 ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 50 ierr = DMCreateDS(dm); CHKERRQ(ierr); 51 ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 52 // Set wall BCs 53 if (bc->num_wall > 0) { 54 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 55 bc->num_wall, bc->walls, 0, bc->num_comps, 56 bc->wall_comps, (void(*)(void))problem->bc, 57 NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 58 } 59 // Set slip BCs in the x direction 60 if (bc->num_slip[0] > 0) { 61 PetscInt comps[1] = {1}; 62 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", label, 63 bc->num_slip[0], bc->slips[0], 0, 1, comps, 64 (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 65 } 66 // Set slip BCs in the y direction 67 if (bc->num_slip[1] > 0) { 68 PetscInt comps[1] = {2}; 69 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", label, 70 bc->num_slip[1], bc->slips[1], 0, 1, comps, 71 (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 72 } 73 // Set slip BCs in the z direction 74 if (bc->num_slip[2] > 0) { 75 PetscInt comps[1] = {3}; 76 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", label, 77 bc->num_slip[2], bc->slips[2], 0, 1, comps, 78 (void(*)(void))NULL, NULL, problem->bc_ctx, NULL); CHKERRQ(ierr); 79 } 80 { 81 PetscBool use_strongstg = PETSC_FALSE; 82 ierr = PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL); 83 CHKERRQ(ierr); 84 85 if (use_strongstg) { 86 ierr = SetupStrongSTG(dm, bc, problem); CHKERRQ(ierr); 87 } 88 } 89 90 ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 91 CHKERRQ(ierr); 92 ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 93 } 94 { 95 // Empty name for conserved field (because there is only one field) 96 PetscSection section; 97 ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 98 ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 99 ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 100 CHKERRQ(ierr); 101 ierr = PetscSectionSetComponentName(section, 0, 1, "Momentum X"); 102 CHKERRQ(ierr); 103 ierr = PetscSectionSetComponentName(section, 0, 2, "Momentum Y"); 104 CHKERRQ(ierr); 105 ierr = PetscSectionSetComponentName(section, 0, 3, "Momentum Z"); 106 CHKERRQ(ierr); 107 ierr = PetscSectionSetComponentName(section, 0, 4, "Energy Density"); 108 CHKERRQ(ierr); 109 } 110 PetscFunctionReturn(0); 111 } 112 113 // Refine DM for high-order viz 114 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, 115 SimpleBC bc, Physics phys) { 116 PetscErrorCode ierr; 117 DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 118 VecType vec_type; 119 PetscFunctionBeginUser; 120 121 ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 122 123 dm_hierarchy[0] = dm; 124 for (PetscInt i = 0, d = user->app_ctx->degree; 125 i < user->app_ctx->viz_refine; i++) { 126 Mat interp_next; 127 ierr = DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i+1]); 128 CHKERRQ(ierr); 129 ierr = DMClearDS(dm_hierarchy[i+1]); CHKERRQ(ierr); 130 ierr = DMClearFields(dm_hierarchy[i+1]); CHKERRQ(ierr); 131 ierr = DMSetCoarseDM(dm_hierarchy[i+1], dm_hierarchy[i]); CHKERRQ(ierr); 132 d = (d + 1) / 2; 133 if (i + 1 == user->app_ctx->viz_refine) d = 1; 134 ierr = DMGetVecType(dm, &vec_type); CHKERRQ(ierr); 135 ierr = DMSetVecType(dm_hierarchy[i+1], vec_type); CHKERRQ(ierr); 136 ierr = SetUpDM(dm_hierarchy[i+1], problem, d, bc, phys); 137 CHKERRQ(ierr); 138 ierr = DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i+1], &interp_next, 139 NULL); CHKERRQ(ierr); 140 if (!i) user->interp_viz = interp_next; 141 else { 142 Mat C; 143 ierr = MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, 144 PETSC_DECIDE, &C); CHKERRQ(ierr); 145 ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 146 ierr = MatDestroy(&user->interp_viz); CHKERRQ(ierr); 147 user->interp_viz = C; 148 } 149 } 150 for (PetscInt i=1; i<user->app_ctx->viz_refine; i++) { 151 ierr = DMDestroy(&dm_hierarchy[i]); CHKERRQ(ierr); 152 } 153 user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 154 155 PetscFunctionReturn(0); 156 } 157