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 <ceed.h> 12 #include <petscdmplex.h> 13 #include <petscds.h> 14 15 #include "../navierstokes.h" 16 #include "../problems/stg_shur14.h" 17 18 // Create mesh 19 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType mat_type, VecType vec_type, DM *dm) { 20 PetscFunctionBeginUser; 21 // Create DMPLEX 22 PetscCall(DMCreate(comm, dm)); 23 PetscCall(DMSetType(*dm, DMPLEX)); 24 { 25 PetscBool skip = PETSC_TRUE; 26 PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_mat_preallocate_skip", &skip, NULL)); 27 PetscCall(DMSetMatrixPreallocateSkip(*dm, skip)); 28 } 29 PetscCall(DMSetMatType(*dm, mat_type)); 30 PetscCall(DMSetVecType(*dm, vec_type)); 31 32 // Set Tensor elements 33 PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_simplex", "0")); 34 PetscCall(PetscOptionsSetValue(NULL, "-dm_sparse_localize", "0")); 35 PetscCall(PetscOptionsSetValue(NULL, "-dm_blocking_type", "field_node")); 36 37 // Set CL options 38 PetscCall(DMSetFromOptions(*dm)); 39 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 40 PetscFunctionReturn(PETSC_SUCCESS); 41 } 42 43 // Setup DM 44 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys) { 45 PetscInt num_comp_q = 5; 46 PetscFunctionBeginUser; 47 48 PetscCall(DMSetupByOrderBegin_FEM(PETSC_TRUE, PETSC_TRUE, degree, PETSC_DECIDE, q_extra, 1, &num_comp_q, dm)); 49 50 { // Add strong boundary conditions to DM 51 DMLabel label; 52 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 53 PetscCall(DMPlexLabelComplete(dm, label)); 54 // Set wall BCs 55 if (bc->num_wall > 0) { 56 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, bc->num_wall, bc->walls, 0, bc->num_comps, bc->wall_comps, NULL, NULL, NULL, NULL)); 57 } 58 // Set symmetry BCs in the x direction 59 if (bc->num_symmetry[0] > 0) { 60 PetscInt comps[1] = {1}; 61 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_x", label, bc->num_symmetry[0], bc->symmetries[0], 0, 1, comps, NULL, NULL, NULL, NULL)); 62 } 63 // Set symmetry BCs in the y direction 64 if (bc->num_symmetry[1] > 0) { 65 PetscInt comps[1] = {2}; 66 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_y", label, bc->num_symmetry[1], bc->symmetries[1], 0, 1, comps, NULL, NULL, NULL, NULL)); 67 } 68 // Set symmetry BCs in the z direction 69 if (bc->num_symmetry[2] > 0) { 70 PetscInt comps[1] = {3}; 71 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "symmetry_z", label, bc->num_symmetry[2], bc->symmetries[2], 0, 1, comps, NULL, NULL, NULL, NULL)); 72 } 73 { 74 PetscBool use_strongstg = PETSC_FALSE; 75 PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 76 if (use_strongstg) PetscCall(SetupStrongStg(dm, bc, problem, phys)); 77 } 78 } 79 80 PetscCall(DMSetupByOrderEnd_FEM(PETSC_TRUE, dm)); 81 82 // Empty name for conserved field (because there is only one field) 83 PetscSection section; 84 PetscCall(DMGetLocalSection(dm, §ion)); 85 PetscCall(PetscSectionSetFieldName(section, 0, "")); 86 switch (phys->state_var) { 87 case STATEVAR_CONSERVATIVE: 88 PetscCall(PetscSectionSetComponentName(section, 0, 0, "Density")); 89 PetscCall(PetscSectionSetComponentName(section, 0, 1, "MomentumX")); 90 PetscCall(PetscSectionSetComponentName(section, 0, 2, "MomentumY")); 91 PetscCall(PetscSectionSetComponentName(section, 0, 3, "MomentumZ")); 92 PetscCall(PetscSectionSetComponentName(section, 0, 4, "TotalEnergy")); 93 break; 94 95 case STATEVAR_PRIMITIVE: 96 PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure")); 97 PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX")); 98 PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY")); 99 PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ")); 100 PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature")); 101 break; 102 } 103 PetscFunctionReturn(PETSC_SUCCESS); 104 } 105 106 // Refine DM for high-order viz 107 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys) { 108 DM dm_hierarchy[user->app_ctx->viz_refine + 1]; 109 VecType vec_type; 110 111 PetscFunctionBeginUser; 112 PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 113 114 dm_hierarchy[0] = dm; 115 for (PetscInt i = 0, d = user->app_ctx->degree; i < user->app_ctx->viz_refine; i++) { 116 Mat interp_next; 117 PetscCall(DMRefine(dm_hierarchy[i], MPI_COMM_NULL, &dm_hierarchy[i + 1])); 118 PetscCall(DMClearDS(dm_hierarchy[i + 1])); 119 PetscCall(DMClearFields(dm_hierarchy[i + 1])); 120 PetscCall(DMSetCoarseDM(dm_hierarchy[i + 1], dm_hierarchy[i])); 121 d = (d + 1) / 2; 122 PetscInt q_order = d + user->app_ctx->q_extra; 123 if (i + 1 == user->app_ctx->viz_refine) d = 1; 124 PetscCall(DMGetVecType(dm, &vec_type)); 125 PetscCall(DMSetVecType(dm_hierarchy[i + 1], vec_type)); 126 PetscCall(SetUpDM(dm_hierarchy[i + 1], problem, d, q_order, bc, phys)); 127 PetscCall(DMCreateInterpolation(dm_hierarchy[i], dm_hierarchy[i + 1], &interp_next, NULL)); 128 if (!i) user->interp_viz = interp_next; 129 else { 130 Mat C; 131 PetscCall(MatMatMult(interp_next, user->interp_viz, MAT_INITIAL_MATRIX, PETSC_DECIDE, &C)); 132 PetscCall(MatDestroy(&interp_next)); 133 PetscCall(MatDestroy(&user->interp_viz)); 134 user->interp_viz = C; 135 } 136 } 137 for (PetscInt i = 1; i < user->app_ctx->viz_refine; i++) { 138 PetscCall(DMDestroy(&dm_hierarchy[i])); 139 } 140 user->dm_viz = dm_hierarchy[user->app_ctx->viz_refine]; 141 PetscFunctionReturn(PETSC_SUCCESS); 142 } 143