// SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause /// @file /// Utility functions for setting up SHOCKTUBE #include "../qfunctions/shocktube.h" #include #include #include static PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx) { MPI_Comm comm = honee->comm; PetscFunctionBeginUser; PetscCall(PetscPrintf(comm, " Problem:\n" " Problem Name : %s\n", app_ctx->problem_name)); PetscFunctionReturn(PETSC_SUCCESS); } static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"}; PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx) { SetupContextShock setup_context; Honee honee = *(Honee *)ctx; MPI_Comm comm = honee->comm; Ceed ceed = honee->ceed; PetscBool implicit; PetscBool yzb; PetscInt stab; ShockTubeContext shocktube_ctx; CeedQFunctionContext shocktube_qfctx, ics_qfctx; PetscFunctionBeginUser; PetscCall(PetscNew(&setup_context)); PetscCall(PetscNew(&shocktube_ctx)); // ------------------------------------------------------ // Create the QFunction context // ------------------------------------------------------ // Driver section initial conditions CeedScalar P_high = 1.0; // Pa CeedScalar rho_high = 1.0; // kg/m^3 // Driven section initial conditions CeedScalar P_low = 0.1; // Pa CeedScalar rho_low = 0.125; // kg/m^3 // Stabilization parameter CeedScalar c_tau = 0.5; // -, based on Hughes et al (2010) // Tuning parameters for the YZB shock capturing CeedScalar Cyzb = 0.1; // -, used in approximation of (Na),x CeedScalar Byzb = 2.0; // -, 1 for smooth shocks // 2 for sharp shocks PetscReal domain_min[3], domain_max[3], domain_size[3]; PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; // ------------------------------------------------------ // Command line Options // ------------------------------------------------------ PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL); // -- Numerical formulation options PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL)); // -- Warnings if (stab == STAB_SUPG) { PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n")); } if (yzb && implicit) { PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n")); } PetscOptionsEnd(); // ------------------------------------------------------ // Set up the QFunction context // ------------------------------------------------------ // -- Scale variables to desired units CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]); // -- Setup Context setup_context->mid_point = mid_point; setup_context->time = 0.0; setup_context->P_high = P_high; setup_context->rho_high = rho_high; setup_context->P_low = P_low; setup_context->rho_low = rho_low; // -- QFunction Context honee->phys->implicit = implicit; shocktube_ctx->implicit = implicit; shocktube_ctx->stabilization = stab; shocktube_ctx->yzb = yzb; shocktube_ctx->Cyzb = Cyzb; shocktube_ctx->Byzb = Byzb; shocktube_ctx->c_tau = c_tau; PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx)); PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc)); PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &shocktube_qfctx)); PetscCallCeed(ceed, CeedQFunctionContextSetData(shocktube_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx)); PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(shocktube_qfctx, CEED_MEM_HOST, FreeContextPetsc)); // ------------------------------------------------------ // SET UP SHOCKTUBE // ------------------------------------------------------ problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsShockTube, .qf_loc = ICsShockTube_loc, .qfctx = ics_qfctx}; problem->apply_vol_rhs = (HoneeQFSpec){.qf_func_ptr = EulerShockTube, .qf_loc = EulerShockTube_loc, .qfctx = shocktube_qfctx}; problem->num_comps_jac_data = 0; problem->compute_exact_solution_error = PETSC_FALSE; problem->print_info = PRINT_SHOCKTUBE; problem->num_components = 5; PetscCall(PetscMalloc1(problem->num_components, &problem->component_names)); for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i])); PetscFunctionReturn(PETSC_SUCCESS); }