1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2*ae2b091fSJames 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 14991aef52SJames Wright PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 153636f6a4SJames Wright SetupContextShock setup_context; 16af8870a9STimothy Aiken User user = *(User *)ctx; 17b4c37c5cSJames Wright MPI_Comm comm = user->comm; 18b4c37c5cSJames Wright Ceed ceed = user->ceed; 19af8870a9STimothy Aiken PetscBool implicit; 20af8870a9STimothy Aiken PetscBool yzb; 21af8870a9STimothy Aiken PetscInt stab; 2215a3537eSJed Brown ShockTubeContext shocktube_ctx; 2315a3537eSJed Brown CeedQFunctionContext shocktube_context; 24af8870a9STimothy Aiken 2515a3537eSJed Brown PetscFunctionBeginUser; 262b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 272b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &shocktube_ctx)); 28af8870a9STimothy Aiken 29af8870a9STimothy Aiken // ------------------------------------------------------ 30af8870a9STimothy Aiken // SET UP SHOCKTUBE 31af8870a9STimothy Aiken // ------------------------------------------------------ 32af8870a9STimothy Aiken problem->dim = 3; 339785fe93SJed Brown problem->ics.qfunction = ICsShockTube; 349785fe93SJed Brown problem->ics.qfunction_loc = ICsShockTube_loc; 359785fe93SJed Brown problem->apply_vol_rhs.qfunction = EulerShockTube; 369785fe93SJed Brown problem->apply_vol_rhs.qfunction_loc = EulerShockTube_loc; 379785fe93SJed Brown problem->apply_vol_ifunction.qfunction = NULL; 389785fe93SJed Brown problem->apply_vol_ifunction.qfunction_loc = NULL; 3958ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 40af8870a9STimothy Aiken problem->print_info = PRINT_SHOCKTUBE; 41af8870a9STimothy Aiken 42af8870a9STimothy Aiken // ------------------------------------------------------ 43af8870a9STimothy Aiken // Create the libCEED 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 // Create the PETSc context 63af8870a9STimothy Aiken // ------------------------------------------------------ 64af8870a9STimothy Aiken PetscScalar meter = 1e-2; // 1 meter in scaled length units 65af8870a9STimothy Aiken PetscScalar second = 1e-2; // 1 second in scaled time units 66af8870a9STimothy Aiken 67af8870a9STimothy Aiken // ------------------------------------------------------ 68af8870a9STimothy Aiken // Command line Options 69af8870a9STimothy Aiken // ------------------------------------------------------ 70af8870a9STimothy Aiken PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL); 71af8870a9STimothy Aiken 72af8870a9STimothy Aiken // -- Numerical formulation options 732b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 742b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 752b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 762b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL)); 77af8870a9STimothy Aiken 78af8870a9STimothy Aiken // -- Units 792b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 80af8870a9STimothy Aiken meter = fabs(meter); 812b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 82af8870a9STimothy Aiken second = fabs(second); 83af8870a9STimothy Aiken 84af8870a9STimothy Aiken // -- Warnings 85af8870a9STimothy Aiken if (stab == STAB_SUPG) { 862b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n")); 87af8870a9STimothy Aiken } 88af8870a9STimothy Aiken if (yzb && implicit) { 892b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n")); 90af8870a9STimothy Aiken } 91af8870a9STimothy Aiken 92af8870a9STimothy Aiken PetscOptionsEnd(); 93af8870a9STimothy Aiken 94af8870a9STimothy Aiken // ------------------------------------------------------ 95af8870a9STimothy Aiken // Set up the PETSc context 96af8870a9STimothy Aiken // ------------------------------------------------------ 97af8870a9STimothy Aiken user->units->meter = meter; 98af8870a9STimothy Aiken user->units->second = second; 99af8870a9STimothy Aiken 100af8870a9STimothy Aiken // ------------------------------------------------------ 101af8870a9STimothy Aiken // Set up the libCEED context 102af8870a9STimothy Aiken // ------------------------------------------------------ 103af8870a9STimothy Aiken // -- Scale variables to desired units 104493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) { 105af8870a9STimothy Aiken domain_size[i] *= meter; 106af8870a9STimothy Aiken domain_min[i] *= meter; 107af8870a9STimothy Aiken } 108af8870a9STimothy Aiken problem->dm_scale = meter; 109af8870a9STimothy Aiken CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]); 110af8870a9STimothy Aiken 111af8870a9STimothy Aiken // -- Setup Context 112af8870a9STimothy Aiken setup_context->mid_point = mid_point; 113af8870a9STimothy Aiken setup_context->time = 0.0; 114af8870a9STimothy Aiken setup_context->P_high = P_high; 115af8870a9STimothy Aiken setup_context->rho_high = rho_high; 116af8870a9STimothy Aiken setup_context->P_low = P_low; 117af8870a9STimothy Aiken setup_context->rho_low = rho_low; 118af8870a9STimothy Aiken 119af8870a9STimothy Aiken // -- QFunction Context 120af8870a9STimothy Aiken user->phys->implicit = implicit; 12115a3537eSJed Brown shocktube_ctx->implicit = implicit; 12215a3537eSJed Brown shocktube_ctx->stabilization = stab; 12315a3537eSJed Brown shocktube_ctx->yzb = yzb; 12415a3537eSJed Brown shocktube_ctx->Cyzb = Cyzb; 12515a3537eSJed Brown shocktube_ctx->Byzb = Byzb; 12615a3537eSJed Brown shocktube_ctx->c_tau = c_tau; 127af8870a9STimothy Aiken 128b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 129b4c37c5cSJames Wright PetscCallCeed(ceed, 130b4c37c5cSJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 131b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 132af8870a9STimothy Aiken 133b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &shocktube_context)); 134b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(shocktube_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx)); 135b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(shocktube_context, CEED_MEM_HOST, FreeContextPetsc)); 13615a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = shocktube_context; 137d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 138af8870a9STimothy Aiken } 139af8870a9STimothy Aiken 140991aef52SJames Wright PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx) { 1412d49c0afSJames Wright MPI_Comm comm = user->comm; 142af8870a9STimothy Aiken 14306f41313SJames Wright PetscFunctionBeginUser; 1442b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 145af8870a9STimothy Aiken " Problem:\n" 146af8870a9STimothy Aiken " Problem Name : %s\n", 1472b916ea7SJeremy L Thompson app_ctx->problem_name)); 148d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 149af8870a9STimothy Aiken } 150