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