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