1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3af8870a9STimothy Aiken 4af8870a9STimothy Aiken /// @file 5af8870a9STimothy Aiken /// Utility functions for setting up SHOCKTUBE 6af8870a9STimothy Aiken 7af8870a9STimothy Aiken #include "../qfunctions/shocktube.h" 8af8870a9STimothy Aiken 9e419654dSJeremy L Thompson #include <ceed.h> 10e419654dSJeremy L Thompson #include <petscdm.h> 11e419654dSJeremy L Thompson 12149fb536SJames Wright #include <navierstokes.h> 136dd99bedSLeila Ghaffari 149b05e62eSJames Wright static PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx) { 159b05e62eSJames Wright MPI_Comm comm = honee->comm; 169b05e62eSJames Wright 179b05e62eSJames Wright PetscFunctionBeginUser; 189b05e62eSJames Wright PetscCall(PetscPrintf(comm, 199b05e62eSJames Wright " Problem:\n" 209b05e62eSJames Wright " Problem Name : %s\n", 219b05e62eSJames Wright app_ctx->problem_name)); 229b05e62eSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 239b05e62eSJames Wright } 249b05e62eSJames Wright 25f27c2204SJames Wright static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; 26f27c2204SJames Wright 27d3c60affSJames Wright PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx) { 283636f6a4SJames Wright SetupContextShock setup_context; 290c373b74SJames Wright Honee honee = *(Honee *)ctx; 300c373b74SJames Wright MPI_Comm comm = honee->comm; 310c373b74SJames Wright Ceed ceed = honee->ceed; 32af8870a9STimothy Aiken PetscBool implicit; 33af8870a9STimothy Aiken PetscBool yzb; 34af8870a9STimothy Aiken PetscInt stab; 3515a3537eSJed Brown ShockTubeContext shocktube_ctx; 36*f5dc303cSJames Wright CeedQFunctionContext shocktube_qfctx, ics_qfctx; 37af8870a9STimothy Aiken 3815a3537eSJed Brown PetscFunctionBeginUser; 392d898fa6SJames Wright PetscCall(PetscNew(&setup_context)); 402d898fa6SJames Wright PetscCall(PetscNew(&shocktube_ctx)); 41af8870a9STimothy Aiken 42af8870a9STimothy Aiken // ------------------------------------------------------ 43ea615d4cSJames Wright // Create the QFunction context 44af8870a9STimothy Aiken // ------------------------------------------------------ 45af8870a9STimothy Aiken // Driver section initial conditions 46af8870a9STimothy Aiken CeedScalar P_high = 1.0; // Pa 47af8870a9STimothy Aiken CeedScalar rho_high = 1.0; // kg/m^3 48af8870a9STimothy Aiken // Driven section initial conditions 49af8870a9STimothy Aiken CeedScalar P_low = 0.1; // Pa 50af8870a9STimothy Aiken CeedScalar rho_low = 0.125; // kg/m^3 51af8870a9STimothy Aiken // Stabilization parameter 52af8870a9STimothy Aiken CeedScalar c_tau = 0.5; // -, based on Hughes et al (2010) 53af8870a9STimothy Aiken // Tuning parameters for the YZB shock capturing 54af8870a9STimothy Aiken CeedScalar Cyzb = 0.1; // -, used in approximation of (Na),x 55af8870a9STimothy Aiken CeedScalar Byzb = 2.0; // -, 1 for smooth shocks 56af8870a9STimothy Aiken // 2 for sharp shocks 57af8870a9STimothy Aiken PetscReal domain_min[3], domain_max[3], domain_size[3]; 582b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 59493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 60af8870a9STimothy Aiken 61af8870a9STimothy Aiken // ------------------------------------------------------ 62af8870a9STimothy Aiken // Command line Options 63af8870a9STimothy Aiken // ------------------------------------------------------ 64af8870a9STimothy Aiken PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL); 65af8870a9STimothy Aiken 66af8870a9STimothy Aiken // -- Numerical formulation options 672b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 682b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 692b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 702b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL)); 71af8870a9STimothy Aiken 72af8870a9STimothy Aiken // -- Warnings 73af8870a9STimothy Aiken if (stab == STAB_SUPG) { 742b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n")); 75af8870a9STimothy Aiken } 76af8870a9STimothy Aiken if (yzb && implicit) { 772b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n")); 78af8870a9STimothy Aiken } 79af8870a9STimothy Aiken 80af8870a9STimothy Aiken PetscOptionsEnd(); 81af8870a9STimothy Aiken 82af8870a9STimothy Aiken // ------------------------------------------------------ 83ea615d4cSJames Wright // Set up the QFunction context 84af8870a9STimothy Aiken // ------------------------------------------------------ 85af8870a9STimothy Aiken // -- Scale variables to desired units 86c9f37605SMohammed Amin Units units = honee->units; 87493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) { 88c9f37605SMohammed Amin domain_size[i] *= units->meter; 89c9f37605SMohammed Amin domain_min[i] *= units->meter; 90af8870a9STimothy Aiken } 91af8870a9STimothy Aiken CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]); 92af8870a9STimothy Aiken 93af8870a9STimothy Aiken // -- Setup Context 94af8870a9STimothy Aiken setup_context->mid_point = mid_point; 95af8870a9STimothy Aiken setup_context->time = 0.0; 96af8870a9STimothy Aiken setup_context->P_high = P_high; 97af8870a9STimothy Aiken setup_context->rho_high = rho_high; 98af8870a9STimothy Aiken setup_context->P_low = P_low; 99af8870a9STimothy Aiken setup_context->rho_low = rho_low; 100af8870a9STimothy Aiken 101af8870a9STimothy Aiken // -- QFunction Context 1020c373b74SJames Wright honee->phys->implicit = implicit; 10315a3537eSJed Brown shocktube_ctx->implicit = implicit; 10415a3537eSJed Brown shocktube_ctx->stabilization = stab; 10515a3537eSJed Brown shocktube_ctx->yzb = yzb; 10615a3537eSJed Brown shocktube_ctx->Cyzb = Cyzb; 10715a3537eSJed Brown shocktube_ctx->Byzb = Byzb; 10815a3537eSJed Brown shocktube_ctx->c_tau = c_tau; 109af8870a9STimothy Aiken 110*f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx)); 111*f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 112*f5dc303cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 113af8870a9STimothy Aiken 1140c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &shocktube_qfctx)); 115e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(shocktube_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx)); 116e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(shocktube_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 117*f5dc303cSJames Wright 118*f5dc303cSJames Wright // ------------------------------------------------------ 119*f5dc303cSJames Wright // SET UP SHOCKTUBE 120*f5dc303cSJames Wright // ------------------------------------------------------ 121*f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsShockTube, .qf_loc = ICsShockTube_loc, .qfctx = ics_qfctx}; 122*f5dc303cSJames Wright problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = EulerShockTube, .qf_loc = EulerShockTube_loc, .qfctx = shocktube_qfctx}; 123*f5dc303cSJames Wright problem->num_comps_jac_data = 0; 124*f5dc303cSJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 125*f5dc303cSJames Wright problem->print_info = PRINT_SHOCKTUBE; 126*f5dc303cSJames Wright 127*f5dc303cSJames Wright problem->num_components = 5; 128*f5dc303cSJames Wright PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); 129*f5dc303cSJames Wright for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i])); 130d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 131af8870a9STimothy Aiken } 132