1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 29ed3d70dSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 39ed3d70dSJames Wright // 49ed3d70dSJames Wright // SPDX-License-Identifier: BSD-2-Clause 59ed3d70dSJames Wright // 69ed3d70dSJames Wright // This file is part of CEED: http://github.com/ceed 79ed3d70dSJames Wright 89ed3d70dSJames Wright /// @file 99ed3d70dSJames Wright /// Utility functions for setting up Freestream boundary condition 109ed3d70dSJames Wright 119ed3d70dSJames Wright #include "../qfunctions/bc_freestream.h" 129ed3d70dSJames Wright 139ed3d70dSJames Wright #include <ceed.h> 149ed3d70dSJames Wright #include <petscdm.h> 159ed3d70dSJames Wright 169ed3d70dSJames Wright #include "../navierstokes.h" 179ed3d70dSJames Wright #include "../qfunctions/newtonian_types.h" 189ed3d70dSJames Wright 199ed3d70dSJames Wright static const char *const RiemannSolverTypes[] = {"hll", "hllc", "RiemannSolverTypes", "RIEMANN_", NULL}; 209ed3d70dSJames Wright 21*991aef52SJames Wright PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) { 229ed3d70dSJames Wright User user = *(User *)ctx; 239ed3d70dSJames Wright MPI_Comm comm = user->comm; 249ed3d70dSJames Wright Ceed ceed = user->ceed; 259ed3d70dSJames Wright FreestreamContext freestream_ctx; 269ed3d70dSJames Wright CeedQFunctionContext freestream_context; 279ed3d70dSJames Wright RiemannFluxType riemann = RIEMANN_HLLC; 289ed3d70dSJames Wright PetscScalar meter = user->units->meter; 299ed3d70dSJames Wright PetscScalar second = user->units->second; 309ed3d70dSJames Wright PetscScalar Kelvin = user->units->Kelvin; 319ed3d70dSJames Wright PetscScalar Pascal = user->units->Pascal; 329ed3d70dSJames Wright 339ed3d70dSJames Wright PetscFunctionBeginUser; 349ed3d70dSJames Wright // Freestream inherits reference state. We re-dimensionalize so the defaults 359ed3d70dSJames Wright // in -help will be visible in SI units. 369ed3d70dSJames Wright StatePrimitive Y_inf = {.pressure = reference->pressure / Pascal, .velocity = {0}, .temperature = reference->temperature / Kelvin}; 379ed3d70dSJames Wright for (int i = 0; i < 3; i++) Y_inf.velocity[i] = reference->velocity[i] * second / meter; 389ed3d70dSJames Wright 399ed3d70dSJames Wright PetscOptionsBegin(comm, NULL, "Options for Freestream boundary condition", NULL); 409ed3d70dSJames Wright PetscCall(PetscOptionsEnum("-freestream_riemann", "Riemann solver to use in freestream boundary condition", NULL, RiemannSolverTypes, 419ed3d70dSJames Wright (PetscEnum)riemann, (PetscEnum *)&riemann, NULL)); 429ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-freestream_pressure", "Pressure at freestream condition", NULL, Y_inf.pressure, &Y_inf.pressure, NULL)); 439ed3d70dSJames Wright PetscInt narray = 3; 449ed3d70dSJames Wright PetscCall(PetscOptionsScalarArray("-freestream_velocity", "Velocity at freestream condition", NULL, Y_inf.velocity, &narray, NULL)); 459ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-freestream_temperature", "Temperature at freestream condition", NULL, Y_inf.temperature, &Y_inf.temperature, NULL)); 469ed3d70dSJames Wright PetscOptionsEnd(); 479ed3d70dSJames Wright 489ed3d70dSJames Wright switch (user->phys->state_var) { 499ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 509ed3d70dSJames Wright switch (riemann) { 519ed3d70dSJames Wright case RIEMANN_HLL: 529ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Conserv_HLL; 539ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Conserv_HLL_loc; 549ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Conserv_HLL; 559ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLL_loc; 569ed3d70dSJames Wright break; 579ed3d70dSJames Wright case RIEMANN_HLLC: 589ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Conserv_HLLC; 599ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Conserv_HLLC_loc; 609ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Conserv_HLLC; 619ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLLC_loc; 629ed3d70dSJames Wright break; 639ed3d70dSJames Wright } 649ed3d70dSJames Wright break; 659ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 669ed3d70dSJames Wright switch (riemann) { 679ed3d70dSJames Wright case RIEMANN_HLL: 689ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Prim_HLL; 699ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Prim_HLL_loc; 709ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Prim_HLL; 719ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLL_loc; 729ed3d70dSJames Wright break; 739ed3d70dSJames Wright case RIEMANN_HLLC: 749ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Prim_HLLC; 759ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Prim_HLLC_loc; 769ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Prim_HLLC; 779ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLLC_loc; 789ed3d70dSJames Wright break; 799ed3d70dSJames Wright } 809ed3d70dSJames Wright break; 819ed3d70dSJames Wright } 829ed3d70dSJames Wright 839ed3d70dSJames Wright Y_inf.pressure *= Pascal; 849ed3d70dSJames Wright for (int i = 0; i < 3; i++) Y_inf.velocity[i] *= meter / second; 859ed3d70dSJames Wright Y_inf.temperature *= Kelvin; 869ed3d70dSJames Wright 879ed3d70dSJames Wright State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf); 889ed3d70dSJames Wright 899ed3d70dSJames Wright // -- Set freestream_ctx struct values 909ed3d70dSJames Wright PetscCall(PetscCalloc1(1, &freestream_ctx)); 919ed3d70dSJames Wright freestream_ctx->newtonian_ctx = *newtonian_ig_ctx; 929ed3d70dSJames Wright freestream_ctx->S_infty = S_infty; 939ed3d70dSJames Wright 949ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &freestream_context)); 959ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(freestream_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*freestream_ctx), freestream_ctx)); 969ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(freestream_context, CEED_MEM_HOST, FreeContextPetsc)); 979ed3d70dSJames Wright problem->apply_freestream.qfunction_context = freestream_context; 989ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(freestream_context, &problem->apply_freestream_jacobian.qfunction_context)); 999ed3d70dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1009ed3d70dSJames Wright } 1019ed3d70dSJames Wright 1029ed3d70dSJames Wright static const char *const OutflowTypes[] = {"RIEMANN", "PRESSURE", "OutflowType", "OUTFLOW_", NULL}; 1039ed3d70dSJames Wright typedef enum { 1049ed3d70dSJames Wright OUTFLOW_RIEMANN, 1059ed3d70dSJames Wright OUTFLOW_PRESSURE, 1069ed3d70dSJames Wright } OutflowType; 1079ed3d70dSJames Wright 108*991aef52SJames Wright PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) { 1099ed3d70dSJames Wright User user = *(User *)ctx; 1109ed3d70dSJames Wright Ceed ceed = user->ceed; 1119ed3d70dSJames Wright OutflowContext outflow_ctx; 1129ed3d70dSJames Wright OutflowType outflow_type = OUTFLOW_RIEMANN; 1139ed3d70dSJames Wright CeedQFunctionContext outflow_context; 1149ed3d70dSJames Wright const PetscScalar Kelvin = user->units->Kelvin; 1159ed3d70dSJames Wright const PetscScalar Pascal = user->units->Pascal; 1169ed3d70dSJames Wright 1179ed3d70dSJames Wright PetscFunctionBeginUser; 1189ed3d70dSJames Wright CeedScalar pressure = reference->pressure / Pascal; 1199ed3d70dSJames Wright CeedScalar temperature = reference->temperature / Kelvin; 1209ed3d70dSJames Wright CeedScalar recirc = 1, softplus_velocity = 1e-2; 1219ed3d70dSJames Wright PetscOptionsBegin(user->comm, NULL, "Options for Outflow boundary condition", NULL); 1229ed3d70dSJames Wright PetscCall( 1239ed3d70dSJames Wright PetscOptionsEnum("-outflow_type", "Type of outflow condition", NULL, OutflowTypes, (PetscEnum)outflow_type, (PetscEnum *)&outflow_type, NULL)); 1249ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-outflow_pressure", "Pressure at outflow condition", NULL, pressure, &pressure, NULL)); 1259ed3d70dSJames Wright if (outflow_type == OUTFLOW_RIEMANN) { 1269ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-outflow_temperature", "Temperature at outflow condition", NULL, temperature, &temperature, NULL)); 1279ed3d70dSJames Wright PetscCall( 1289ed3d70dSJames Wright PetscOptionsReal("-outflow_recirc", "Fraction of recirculation to allow in exterior velocity state [0,1]", NULL, recirc, &recirc, NULL)); 1299ed3d70dSJames Wright PetscCall(PetscOptionsReal("-outflow_softplus_velocity", "Characteristic velocity of softplus regularization", NULL, softplus_velocity, 1309ed3d70dSJames Wright &softplus_velocity, NULL)); 1319ed3d70dSJames Wright } 1329ed3d70dSJames Wright PetscOptionsEnd(); 1339ed3d70dSJames Wright pressure *= Pascal; 1349ed3d70dSJames Wright temperature *= Kelvin; 1359ed3d70dSJames Wright 1369ed3d70dSJames Wright switch (outflow_type) { 1379ed3d70dSJames Wright case OUTFLOW_RIEMANN: 1389ed3d70dSJames Wright switch (user->phys->state_var) { 1399ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 1409ed3d70dSJames Wright problem->apply_outflow.qfunction = RiemannOutflow_Conserv; 1419ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = RiemannOutflow_Conserv_loc; 1429ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = RiemannOutflow_Jacobian_Conserv; 1439ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Conserv_loc; 1449ed3d70dSJames Wright break; 1459ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 1469ed3d70dSJames Wright problem->apply_outflow.qfunction = RiemannOutflow_Prim; 1479ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = RiemannOutflow_Prim_loc; 1489ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = RiemannOutflow_Jacobian_Prim; 1499ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Prim_loc; 1509ed3d70dSJames Wright break; 1519ed3d70dSJames Wright } 1529ed3d70dSJames Wright break; 1539ed3d70dSJames Wright case OUTFLOW_PRESSURE: 1549ed3d70dSJames Wright switch (user->phys->state_var) { 1559ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 1569ed3d70dSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Conserv; 1579ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 1589ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 1599ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 1609ed3d70dSJames Wright break; 1619ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 1629ed3d70dSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Prim; 1639ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 1649ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 1659ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 1669ed3d70dSJames Wright break; 1679ed3d70dSJames Wright } 1689ed3d70dSJames Wright break; 1699ed3d70dSJames Wright } 1709ed3d70dSJames Wright PetscCall(PetscCalloc1(1, &outflow_ctx)); 1719ed3d70dSJames Wright outflow_ctx->gas = *newtonian_ig_ctx; 1729ed3d70dSJames Wright outflow_ctx->recirc = recirc; 1739ed3d70dSJames Wright outflow_ctx->softplus_velocity = softplus_velocity; 1749ed3d70dSJames Wright outflow_ctx->pressure = pressure; 1759ed3d70dSJames Wright outflow_ctx->temperature = temperature; 1769ed3d70dSJames Wright 1779ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &outflow_context)); 1789ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(outflow_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*outflow_ctx), outflow_ctx)); 1799ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(outflow_context, CEED_MEM_HOST, FreeContextPetsc)); 1809ed3d70dSJames Wright problem->apply_outflow.qfunction_context = outflow_context; 1819ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(outflow_context, &problem->apply_outflow_jacobian.qfunction_context)); 1829ed3d70dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1839ed3d70dSJames Wright } 184