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 /// Utility functions for setting up SHOCKTUBE 10 11 #include "../qfunctions/shocktube.h" 12 13 #include <ceed.h> 14 #include <petscdm.h> 15 16 #include "../navierstokes.h" 17 #include "../qfunctions/setupgeo.h" 18 19 PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 20 SetupContextShock setup_context; 21 User user = *(User *)ctx; 22 MPI_Comm comm = PETSC_COMM_WORLD; 23 PetscBool implicit; 24 PetscBool yzb; 25 PetscInt stab; 26 PetscBool has_curr_time = PETSC_FALSE; 27 ShockTubeContext shocktube_ctx; 28 CeedQFunctionContext shocktube_context; 29 30 PetscFunctionBeginUser; 31 PetscCall(PetscCalloc1(1, &setup_context)); 32 PetscCall(PetscCalloc1(1, &shocktube_ctx)); 33 34 // ------------------------------------------------------ 35 // SET UP SHOCKTUBE 36 // ------------------------------------------------------ 37 problem->dim = 3; 38 problem->q_data_size_vol = 10; 39 problem->q_data_size_sur = 4; 40 problem->setup_vol.qfunction = Setup; 41 problem->setup_vol.qfunction_loc = Setup_loc; 42 problem->setup_sur.qfunction = SetupBoundary; 43 problem->setup_sur.qfunction_loc = SetupBoundary_loc; 44 problem->ics.qfunction = ICsShockTube; 45 problem->ics.qfunction_loc = ICsShockTube_loc; 46 problem->apply_vol_rhs.qfunction = EulerShockTube; 47 problem->apply_vol_rhs.qfunction_loc = EulerShockTube_loc; 48 problem->apply_vol_ifunction.qfunction = NULL; 49 problem->apply_vol_ifunction.qfunction_loc = NULL; 50 problem->bc = Exact_ShockTube; 51 problem->bc_ctx = setup_context; 52 problem->non_zero_time = PETSC_FALSE; 53 problem->print_info = PRINT_SHOCKTUBE; 54 55 // ------------------------------------------------------ 56 // Create the libCEED context 57 // ------------------------------------------------------ 58 // Driver section initial conditions 59 CeedScalar P_high = 1.0; // Pa 60 CeedScalar rho_high = 1.0; // kg/m^3 61 // Driven section initial conditions 62 CeedScalar P_low = 0.1; // Pa 63 CeedScalar rho_low = 0.125; // kg/m^3 64 // Stabilization parameter 65 CeedScalar c_tau = 0.5; // -, based on Hughes et al (2010) 66 // Tuning parameters for the YZB shock capturing 67 CeedScalar Cyzb = 0.1; // -, used in approximation of (Na),x 68 CeedScalar Byzb = 2.0; // -, 1 for smooth shocks 69 // 2 for sharp shocks 70 PetscReal domain_min[3], domain_max[3], domain_size[3]; 71 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 72 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 73 74 // ------------------------------------------------------ 75 // Create the PETSc context 76 // ------------------------------------------------------ 77 PetscScalar meter = 1e-2; // 1 meter in scaled length units 78 PetscScalar second = 1e-2; // 1 second in scaled time units 79 80 // ------------------------------------------------------ 81 // Command line Options 82 // ------------------------------------------------------ 83 PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL); 84 85 // -- Numerical formulation options 86 PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 87 PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 88 PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 89 PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL)); 90 91 // -- Units 92 PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 93 meter = fabs(meter); 94 PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 95 second = fabs(second); 96 97 // -- Warnings 98 if (stab == STAB_SUPG) { 99 PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n")); 100 } 101 if (yzb && implicit) { 102 PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n")); 103 } 104 105 PetscOptionsEnd(); 106 107 // ------------------------------------------------------ 108 // Set up the PETSc context 109 // ------------------------------------------------------ 110 user->units->meter = meter; 111 user->units->second = second; 112 113 // ------------------------------------------------------ 114 // Set up the libCEED context 115 // ------------------------------------------------------ 116 // -- Scale variables to desired units 117 for (PetscInt i = 0; i < 3; i++) { 118 domain_size[i] *= meter; 119 domain_min[i] *= meter; 120 } 121 problem->dm_scale = meter; 122 CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]); 123 124 // -- Setup Context 125 setup_context->mid_point = mid_point; 126 setup_context->time = 0.0; 127 setup_context->P_high = P_high; 128 setup_context->rho_high = rho_high; 129 setup_context->P_low = P_low; 130 setup_context->rho_low = rho_low; 131 132 // -- QFunction Context 133 user->phys->implicit = implicit; 134 user->phys->has_curr_time = has_curr_time; 135 shocktube_ctx->implicit = implicit; 136 shocktube_ctx->stabilization = stab; 137 shocktube_ctx->yzb = yzb; 138 shocktube_ctx->Cyzb = Cyzb; 139 shocktube_ctx->Byzb = Byzb; 140 shocktube_ctx->c_tau = c_tau; 141 142 CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context); 143 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context); 144 CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc); 145 146 CeedQFunctionContextCreate(user->ceed, &shocktube_context); 147 CeedQFunctionContextSetData(shocktube_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx); 148 CeedQFunctionContextSetDataDestroy(shocktube_context, CEED_MEM_HOST, FreeContextPetsc); 149 problem->apply_vol_rhs.qfunction_context = shocktube_context; 150 PetscFunctionReturn(0); 151 } 152 153 PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, AppCtx app_ctx) { 154 MPI_Comm comm = PETSC_COMM_WORLD; 155 PetscFunctionBeginUser; 156 157 PetscCall(PetscPrintf(comm, 158 " Problem:\n" 159 " Problem Name : %s\n", 160 app_ctx->problem_name)); 161 162 PetscFunctionReturn(0); 163 } 164