1*9ed3d70dSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2*9ed3d70dSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3*9ed3d70dSJames Wright // 4*9ed3d70dSJames Wright // SPDX-License-Identifier: BSD-2-Clause 5*9ed3d70dSJames Wright // 6*9ed3d70dSJames Wright // This file is part of CEED: http://github.com/ceed 7*9ed3d70dSJames Wright 8*9ed3d70dSJames Wright /// @file 9*9ed3d70dSJames Wright /// Utility functions for setting up Freestream boundary condition 10*9ed3d70dSJames Wright 11*9ed3d70dSJames Wright #include "../qfunctions/bc_freestream.h" 12*9ed3d70dSJames Wright 13*9ed3d70dSJames Wright #include <ceed.h> 14*9ed3d70dSJames Wright #include <petscdm.h> 15*9ed3d70dSJames Wright 16*9ed3d70dSJames Wright #include "../navierstokes.h" 17*9ed3d70dSJames Wright #include "../qfunctions/newtonian_types.h" 18*9ed3d70dSJames Wright 19*9ed3d70dSJames Wright static const char *const RiemannSolverTypes[] = {"hll", "hllc", "RiemannSolverTypes", "RIEMANN_", NULL}; 20*9ed3d70dSJames Wright 21*9ed3d70dSJames Wright PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) { 22*9ed3d70dSJames Wright User user = *(User *)ctx; 23*9ed3d70dSJames Wright MPI_Comm comm = user->comm; 24*9ed3d70dSJames Wright Ceed ceed = user->ceed; 25*9ed3d70dSJames Wright FreestreamContext freestream_ctx; 26*9ed3d70dSJames Wright CeedQFunctionContext freestream_context; 27*9ed3d70dSJames Wright RiemannFluxType riemann = RIEMANN_HLLC; 28*9ed3d70dSJames Wright PetscScalar meter = user->units->meter; 29*9ed3d70dSJames Wright PetscScalar second = user->units->second; 30*9ed3d70dSJames Wright PetscScalar Kelvin = user->units->Kelvin; 31*9ed3d70dSJames Wright PetscScalar Pascal = user->units->Pascal; 32*9ed3d70dSJames Wright 33*9ed3d70dSJames Wright PetscFunctionBeginUser; 34*9ed3d70dSJames Wright // Freestream inherits reference state. We re-dimensionalize so the defaults 35*9ed3d70dSJames Wright // in -help will be visible in SI units. 36*9ed3d70dSJames Wright StatePrimitive Y_inf = {.pressure = reference->pressure / Pascal, .velocity = {0}, .temperature = reference->temperature / Kelvin}; 37*9ed3d70dSJames Wright for (int i = 0; i < 3; i++) Y_inf.velocity[i] = reference->velocity[i] * second / meter; 38*9ed3d70dSJames Wright 39*9ed3d70dSJames Wright PetscOptionsBegin(comm, NULL, "Options for Freestream boundary condition", NULL); 40*9ed3d70dSJames Wright PetscCall(PetscOptionsEnum("-freestream_riemann", "Riemann solver to use in freestream boundary condition", NULL, RiemannSolverTypes, 41*9ed3d70dSJames Wright (PetscEnum)riemann, (PetscEnum *)&riemann, NULL)); 42*9ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-freestream_pressure", "Pressure at freestream condition", NULL, Y_inf.pressure, &Y_inf.pressure, NULL)); 43*9ed3d70dSJames Wright PetscInt narray = 3; 44*9ed3d70dSJames Wright PetscCall(PetscOptionsScalarArray("-freestream_velocity", "Velocity at freestream condition", NULL, Y_inf.velocity, &narray, NULL)); 45*9ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-freestream_temperature", "Temperature at freestream condition", NULL, Y_inf.temperature, &Y_inf.temperature, NULL)); 46*9ed3d70dSJames Wright PetscOptionsEnd(); 47*9ed3d70dSJames Wright 48*9ed3d70dSJames Wright switch (user->phys->state_var) { 49*9ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 50*9ed3d70dSJames Wright switch (riemann) { 51*9ed3d70dSJames Wright case RIEMANN_HLL: 52*9ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Conserv_HLL; 53*9ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Conserv_HLL_loc; 54*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Conserv_HLL; 55*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLL_loc; 56*9ed3d70dSJames Wright break; 57*9ed3d70dSJames Wright case RIEMANN_HLLC: 58*9ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Conserv_HLLC; 59*9ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Conserv_HLLC_loc; 60*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Conserv_HLLC; 61*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLLC_loc; 62*9ed3d70dSJames Wright break; 63*9ed3d70dSJames Wright } 64*9ed3d70dSJames Wright break; 65*9ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 66*9ed3d70dSJames Wright switch (riemann) { 67*9ed3d70dSJames Wright case RIEMANN_HLL: 68*9ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Prim_HLL; 69*9ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Prim_HLL_loc; 70*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Prim_HLL; 71*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLL_loc; 72*9ed3d70dSJames Wright break; 73*9ed3d70dSJames Wright case RIEMANN_HLLC: 74*9ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Prim_HLLC; 75*9ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Prim_HLLC_loc; 76*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Prim_HLLC; 77*9ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLLC_loc; 78*9ed3d70dSJames Wright break; 79*9ed3d70dSJames Wright } 80*9ed3d70dSJames Wright break; 81*9ed3d70dSJames Wright } 82*9ed3d70dSJames Wright 83*9ed3d70dSJames Wright Y_inf.pressure *= Pascal; 84*9ed3d70dSJames Wright for (int i = 0; i < 3; i++) Y_inf.velocity[i] *= meter / second; 85*9ed3d70dSJames Wright Y_inf.temperature *= Kelvin; 86*9ed3d70dSJames Wright 87*9ed3d70dSJames Wright State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf); 88*9ed3d70dSJames Wright 89*9ed3d70dSJames Wright // -- Set freestream_ctx struct values 90*9ed3d70dSJames Wright PetscCall(PetscCalloc1(1, &freestream_ctx)); 91*9ed3d70dSJames Wright freestream_ctx->newtonian_ctx = *newtonian_ig_ctx; 92*9ed3d70dSJames Wright freestream_ctx->S_infty = S_infty; 93*9ed3d70dSJames Wright 94*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &freestream_context)); 95*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(freestream_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*freestream_ctx), freestream_ctx)); 96*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(freestream_context, CEED_MEM_HOST, FreeContextPetsc)); 97*9ed3d70dSJames Wright problem->apply_freestream.qfunction_context = freestream_context; 98*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(freestream_context, &problem->apply_freestream_jacobian.qfunction_context)); 99*9ed3d70dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 100*9ed3d70dSJames Wright } 101*9ed3d70dSJames Wright 102*9ed3d70dSJames Wright static const char *const OutflowTypes[] = {"RIEMANN", "PRESSURE", "OutflowType", "OUTFLOW_", NULL}; 103*9ed3d70dSJames Wright typedef enum { 104*9ed3d70dSJames Wright OUTFLOW_RIEMANN, 105*9ed3d70dSJames Wright OUTFLOW_PRESSURE, 106*9ed3d70dSJames Wright } OutflowType; 107*9ed3d70dSJames Wright 108*9ed3d70dSJames Wright PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) { 109*9ed3d70dSJames Wright User user = *(User *)ctx; 110*9ed3d70dSJames Wright Ceed ceed = user->ceed; 111*9ed3d70dSJames Wright OutflowContext outflow_ctx; 112*9ed3d70dSJames Wright OutflowType outflow_type = OUTFLOW_RIEMANN; 113*9ed3d70dSJames Wright CeedQFunctionContext outflow_context; 114*9ed3d70dSJames Wright const PetscScalar Kelvin = user->units->Kelvin; 115*9ed3d70dSJames Wright const PetscScalar Pascal = user->units->Pascal; 116*9ed3d70dSJames Wright 117*9ed3d70dSJames Wright PetscFunctionBeginUser; 118*9ed3d70dSJames Wright CeedScalar pressure = reference->pressure / Pascal; 119*9ed3d70dSJames Wright CeedScalar temperature = reference->temperature / Kelvin; 120*9ed3d70dSJames Wright CeedScalar recirc = 1, softplus_velocity = 1e-2; 121*9ed3d70dSJames Wright PetscOptionsBegin(user->comm, NULL, "Options for Outflow boundary condition", NULL); 122*9ed3d70dSJames Wright PetscCall( 123*9ed3d70dSJames Wright PetscOptionsEnum("-outflow_type", "Type of outflow condition", NULL, OutflowTypes, (PetscEnum)outflow_type, (PetscEnum *)&outflow_type, NULL)); 124*9ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-outflow_pressure", "Pressure at outflow condition", NULL, pressure, &pressure, NULL)); 125*9ed3d70dSJames Wright if (outflow_type == OUTFLOW_RIEMANN) { 126*9ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-outflow_temperature", "Temperature at outflow condition", NULL, temperature, &temperature, NULL)); 127*9ed3d70dSJames Wright PetscCall( 128*9ed3d70dSJames Wright PetscOptionsReal("-outflow_recirc", "Fraction of recirculation to allow in exterior velocity state [0,1]", NULL, recirc, &recirc, NULL)); 129*9ed3d70dSJames Wright PetscCall(PetscOptionsReal("-outflow_softplus_velocity", "Characteristic velocity of softplus regularization", NULL, softplus_velocity, 130*9ed3d70dSJames Wright &softplus_velocity, NULL)); 131*9ed3d70dSJames Wright } 132*9ed3d70dSJames Wright PetscOptionsEnd(); 133*9ed3d70dSJames Wright pressure *= Pascal; 134*9ed3d70dSJames Wright temperature *= Kelvin; 135*9ed3d70dSJames Wright 136*9ed3d70dSJames Wright switch (outflow_type) { 137*9ed3d70dSJames Wright case OUTFLOW_RIEMANN: 138*9ed3d70dSJames Wright switch (user->phys->state_var) { 139*9ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 140*9ed3d70dSJames Wright problem->apply_outflow.qfunction = RiemannOutflow_Conserv; 141*9ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = RiemannOutflow_Conserv_loc; 142*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = RiemannOutflow_Jacobian_Conserv; 143*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Conserv_loc; 144*9ed3d70dSJames Wright break; 145*9ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 146*9ed3d70dSJames Wright problem->apply_outflow.qfunction = RiemannOutflow_Prim; 147*9ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = RiemannOutflow_Prim_loc; 148*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = RiemannOutflow_Jacobian_Prim; 149*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Prim_loc; 150*9ed3d70dSJames Wright break; 151*9ed3d70dSJames Wright } 152*9ed3d70dSJames Wright break; 153*9ed3d70dSJames Wright case OUTFLOW_PRESSURE: 154*9ed3d70dSJames Wright switch (user->phys->state_var) { 155*9ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 156*9ed3d70dSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Conserv; 157*9ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 158*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 159*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 160*9ed3d70dSJames Wright break; 161*9ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 162*9ed3d70dSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Prim; 163*9ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 164*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 165*9ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 166*9ed3d70dSJames Wright break; 167*9ed3d70dSJames Wright } 168*9ed3d70dSJames Wright break; 169*9ed3d70dSJames Wright } 170*9ed3d70dSJames Wright PetscCall(PetscCalloc1(1, &outflow_ctx)); 171*9ed3d70dSJames Wright outflow_ctx->gas = *newtonian_ig_ctx; 172*9ed3d70dSJames Wright outflow_ctx->recirc = recirc; 173*9ed3d70dSJames Wright outflow_ctx->softplus_velocity = softplus_velocity; 174*9ed3d70dSJames Wright outflow_ctx->pressure = pressure; 175*9ed3d70dSJames Wright outflow_ctx->temperature = temperature; 176*9ed3d70dSJames Wright 177*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &outflow_context)); 178*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(outflow_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*outflow_ctx), outflow_ctx)); 179*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(outflow_context, CEED_MEM_HOST, FreeContextPetsc)); 180*9ed3d70dSJames Wright problem->apply_outflow.qfunction_context = outflow_context; 181*9ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(outflow_context, &problem->apply_outflow_jacobian.qfunction_context)); 182*9ed3d70dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 183*9ed3d70dSJames Wright } 184