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