1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 204e40bb6SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3af8870a9STimothy Aiken // 404e40bb6SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5af8870a9STimothy Aiken // 604e40bb6SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7af8870a9STimothy Aiken 8af8870a9STimothy Aiken /// @file 9af8870a9STimothy Aiken /// Utility functions for setting up SHOCKTUBE 10af8870a9STimothy Aiken 11af8870a9STimothy Aiken #include "../qfunctions/shocktube.h" 12af8870a9STimothy Aiken 13e419654dSJeremy L Thompson #include <ceed.h> 14e419654dSJeremy L Thompson #include <petscdm.h> 15e419654dSJeremy L Thompson 16*149fb536SJames Wright #include <navierstokes.h> 176dd99bedSLeila Ghaffari 18991aef52SJames Wright PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 193636f6a4SJames Wright SetupContextShock setup_context; 20af8870a9STimothy Aiken User user = *(User *)ctx; 21b4c37c5cSJames Wright MPI_Comm comm = user->comm; 22b4c37c5cSJames Wright Ceed ceed = user->ceed; 23af8870a9STimothy Aiken PetscBool implicit; 24af8870a9STimothy Aiken PetscBool yzb; 25af8870a9STimothy Aiken PetscInt stab; 2615a3537eSJed Brown ShockTubeContext shocktube_ctx; 2715a3537eSJed Brown CeedQFunctionContext shocktube_context; 28af8870a9STimothy Aiken 2915a3537eSJed Brown PetscFunctionBeginUser; 302b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 312b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &shocktube_ctx)); 32af8870a9STimothy Aiken 33af8870a9STimothy Aiken // ------------------------------------------------------ 34af8870a9STimothy Aiken // SET UP SHOCKTUBE 35af8870a9STimothy Aiken // ------------------------------------------------------ 36af8870a9STimothy Aiken problem->dim = 3; 379785fe93SJed Brown problem->ics.qfunction = ICsShockTube; 389785fe93SJed Brown problem->ics.qfunction_loc = ICsShockTube_loc; 399785fe93SJed Brown problem->apply_vol_rhs.qfunction = EulerShockTube; 409785fe93SJed Brown problem->apply_vol_rhs.qfunction_loc = EulerShockTube_loc; 419785fe93SJed Brown problem->apply_vol_ifunction.qfunction = NULL; 429785fe93SJed Brown problem->apply_vol_ifunction.qfunction_loc = NULL; 4358ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 44af8870a9STimothy Aiken problem->print_info = PRINT_SHOCKTUBE; 45af8870a9STimothy Aiken 46af8870a9STimothy Aiken // ------------------------------------------------------ 47af8870a9STimothy Aiken // Create the libCEED context 48af8870a9STimothy Aiken // ------------------------------------------------------ 49af8870a9STimothy Aiken // Driver section initial conditions 50af8870a9STimothy Aiken CeedScalar P_high = 1.0; // Pa 51af8870a9STimothy Aiken CeedScalar rho_high = 1.0; // kg/m^3 52af8870a9STimothy Aiken // Driven section initial conditions 53af8870a9STimothy Aiken CeedScalar P_low = 0.1; // Pa 54af8870a9STimothy Aiken CeedScalar rho_low = 0.125; // kg/m^3 55af8870a9STimothy Aiken // Stabilization parameter 56af8870a9STimothy Aiken CeedScalar c_tau = 0.5; // -, based on Hughes et al (2010) 57af8870a9STimothy Aiken // Tuning parameters for the YZB shock capturing 58af8870a9STimothy Aiken CeedScalar Cyzb = 0.1; // -, used in approximation of (Na),x 59af8870a9STimothy Aiken CeedScalar Byzb = 2.0; // -, 1 for smooth shocks 60af8870a9STimothy Aiken // 2 for sharp shocks 61af8870a9STimothy Aiken PetscReal domain_min[3], domain_max[3], domain_size[3]; 622b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 63493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 64af8870a9STimothy Aiken 65af8870a9STimothy Aiken // ------------------------------------------------------ 66af8870a9STimothy Aiken // Create the PETSc context 67af8870a9STimothy Aiken // ------------------------------------------------------ 68af8870a9STimothy Aiken PetscScalar meter = 1e-2; // 1 meter in scaled length units 69af8870a9STimothy Aiken PetscScalar second = 1e-2; // 1 second in scaled time units 70af8870a9STimothy Aiken 71af8870a9STimothy Aiken // ------------------------------------------------------ 72af8870a9STimothy Aiken // Command line Options 73af8870a9STimothy Aiken // ------------------------------------------------------ 74af8870a9STimothy Aiken PetscOptionsBegin(comm, NULL, "Options for SHOCKTUBE problem", NULL); 75af8870a9STimothy Aiken 76af8870a9STimothy Aiken // -- Numerical formulation options 772b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 782b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 792b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-c_tau", "Stabilization constant", NULL, c_tau, &c_tau, NULL)); 802b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-yzb", "Use YZB discontinuity capturing", NULL, yzb = PETSC_FALSE, &yzb, NULL)); 81af8870a9STimothy Aiken 82af8870a9STimothy Aiken // -- Units 832b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 84af8870a9STimothy Aiken meter = fabs(meter); 852b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 86af8870a9STimothy Aiken second = fabs(second); 87af8870a9STimothy Aiken 88af8870a9STimothy Aiken // -- Warnings 89af8870a9STimothy Aiken if (stab == STAB_SUPG) { 902b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! -stab supg not implemented for the shocktube problem. \n")); 91af8870a9STimothy Aiken } 92af8870a9STimothy Aiken if (yzb && implicit) { 932b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! -yzb only implemented for explicit timestepping. \n")); 94af8870a9STimothy Aiken } 95af8870a9STimothy Aiken 96af8870a9STimothy Aiken PetscOptionsEnd(); 97af8870a9STimothy Aiken 98af8870a9STimothy Aiken // ------------------------------------------------------ 99af8870a9STimothy Aiken // Set up the PETSc context 100af8870a9STimothy Aiken // ------------------------------------------------------ 101af8870a9STimothy Aiken user->units->meter = meter; 102af8870a9STimothy Aiken user->units->second = second; 103af8870a9STimothy Aiken 104af8870a9STimothy Aiken // ------------------------------------------------------ 105af8870a9STimothy Aiken // Set up the libCEED context 106af8870a9STimothy Aiken // ------------------------------------------------------ 107af8870a9STimothy Aiken // -- Scale variables to desired units 108493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) { 109af8870a9STimothy Aiken domain_size[i] *= meter; 110af8870a9STimothy Aiken domain_min[i] *= meter; 111af8870a9STimothy Aiken } 112af8870a9STimothy Aiken problem->dm_scale = meter; 113af8870a9STimothy Aiken CeedScalar mid_point = 0.5 * (domain_size[0] + domain_min[0]); 114af8870a9STimothy Aiken 115af8870a9STimothy Aiken // -- Setup Context 116af8870a9STimothy Aiken setup_context->mid_point = mid_point; 117af8870a9STimothy Aiken setup_context->time = 0.0; 118af8870a9STimothy Aiken setup_context->P_high = P_high; 119af8870a9STimothy Aiken setup_context->rho_high = rho_high; 120af8870a9STimothy Aiken setup_context->P_low = P_low; 121af8870a9STimothy Aiken setup_context->rho_low = rho_low; 122af8870a9STimothy Aiken 123af8870a9STimothy Aiken // -- QFunction Context 124af8870a9STimothy Aiken user->phys->implicit = implicit; 12515a3537eSJed Brown shocktube_ctx->implicit = implicit; 12615a3537eSJed Brown shocktube_ctx->stabilization = stab; 12715a3537eSJed Brown shocktube_ctx->yzb = yzb; 12815a3537eSJed Brown shocktube_ctx->Cyzb = Cyzb; 12915a3537eSJed Brown shocktube_ctx->Byzb = Byzb; 13015a3537eSJed Brown shocktube_ctx->c_tau = c_tau; 131af8870a9STimothy Aiken 132b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 133b4c37c5cSJames Wright PetscCallCeed(ceed, 134b4c37c5cSJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 135b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 136af8870a9STimothy Aiken 137b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &shocktube_context)); 138b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(shocktube_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*shocktube_ctx), shocktube_ctx)); 139b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(shocktube_context, CEED_MEM_HOST, FreeContextPetsc)); 14015a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = shocktube_context; 141d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 142af8870a9STimothy Aiken } 143af8870a9STimothy Aiken 144991aef52SJames Wright PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx) { 1452d49c0afSJames Wright MPI_Comm comm = user->comm; 146af8870a9STimothy Aiken 14706f41313SJames Wright PetscFunctionBeginUser; 1482b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 149af8870a9STimothy Aiken " Problem:\n" 150af8870a9STimothy Aiken " Problem Name : %s\n", 1512b916ea7SJeremy L Thompson app_ctx->problem_name)); 152d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 153af8870a9STimothy Aiken } 154