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