xref: /libCEED/examples/fluids/problems/bc_freestream.c (revision a2d72b6f1ed489cbeb0eb5f72cf8bf977e7ff50a)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
29f844368SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
39f844368SJames Wright //
49f844368SJames Wright // SPDX-License-Identifier: BSD-2-Clause
59f844368SJames Wright //
69f844368SJames Wright // This file is part of CEED:  http://github.com/ceed
79f844368SJames Wright 
89f844368SJames Wright /// @file
99f844368SJames Wright /// Utility functions for setting up Freestream boundary condition
109f844368SJames Wright 
119f844368SJames Wright #include "../qfunctions/bc_freestream.h"
129f844368SJames Wright 
139f844368SJames Wright #include <ceed.h>
149f844368SJames Wright #include <petscdm.h>
159f844368SJames Wright 
169f844368SJames Wright #include "../navierstokes.h"
179f844368SJames Wright #include "../qfunctions/newtonian_types.h"
189f844368SJames Wright 
199f844368SJames Wright static const char *const RiemannSolverTypes[] = {"hll", "hllc", "RiemannSolverTypes", "RIEMANN_", NULL};
209f844368SJames Wright 
21731c13d7SJames Wright PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) {
229f844368SJames Wright   User                 user = *(User *)ctx;
239f844368SJames Wright   MPI_Comm             comm = user->comm;
249f844368SJames Wright   Ceed                 ceed = user->ceed;
259f844368SJames Wright   FreestreamContext    freestream_ctx;
269f844368SJames Wright   CeedQFunctionContext freestream_context;
279f844368SJames Wright   RiemannFluxType      riemann = RIEMANN_HLLC;
289f844368SJames Wright   PetscScalar          meter   = user->units->meter;
299f844368SJames Wright   PetscScalar          second  = user->units->second;
309f844368SJames Wright   PetscScalar          Kelvin  = user->units->Kelvin;
319f844368SJames Wright   PetscScalar          Pascal  = user->units->Pascal;
329f844368SJames Wright 
339f844368SJames Wright   PetscFunctionBeginUser;
349f844368SJames Wright   // Freestream inherits reference state. We re-dimensionalize so the defaults
359f844368SJames Wright   // in -help will be visible in SI units.
369f844368SJames Wright   StatePrimitive Y_inf = {.pressure = reference->pressure / Pascal, .velocity = {0}, .temperature = reference->temperature / Kelvin};
379f844368SJames Wright   for (int i = 0; i < 3; i++) Y_inf.velocity[i] = reference->velocity[i] * second / meter;
389f844368SJames Wright 
399f844368SJames Wright   PetscOptionsBegin(comm, NULL, "Options for Freestream boundary condition", NULL);
409f844368SJames Wright   PetscCall(PetscOptionsEnum("-freestream_riemann", "Riemann solver to use in freestream boundary condition", NULL, RiemannSolverTypes,
419f844368SJames Wright                              (PetscEnum)riemann, (PetscEnum *)&riemann, NULL));
429f844368SJames Wright   PetscCall(PetscOptionsScalar("-freestream_pressure", "Pressure at freestream condition", NULL, Y_inf.pressure, &Y_inf.pressure, NULL));
439f844368SJames Wright   PetscInt narray = 3;
449f844368SJames Wright   PetscCall(PetscOptionsScalarArray("-freestream_velocity", "Velocity at freestream condition", NULL, Y_inf.velocity, &narray, NULL));
459f844368SJames Wright   PetscCall(PetscOptionsScalar("-freestream_temperature", "Temperature at freestream condition", NULL, Y_inf.temperature, &Y_inf.temperature, NULL));
469f844368SJames Wright   PetscOptionsEnd();
479f844368SJames Wright 
489f844368SJames Wright   switch (user->phys->state_var) {
499f844368SJames Wright     case STATEVAR_CONSERVATIVE:
509f844368SJames Wright       switch (riemann) {
519f844368SJames Wright         case RIEMANN_HLL:
529f844368SJames Wright           problem->apply_freestream.qfunction              = Freestream_Conserv_HLL;
539f844368SJames Wright           problem->apply_freestream.qfunction_loc          = Freestream_Conserv_HLL_loc;
549f844368SJames Wright           problem->apply_freestream_jacobian.qfunction     = Freestream_Jacobian_Conserv_HLL;
559f844368SJames Wright           problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLL_loc;
569f844368SJames Wright           break;
579f844368SJames Wright         case RIEMANN_HLLC:
589f844368SJames Wright           problem->apply_freestream.qfunction              = Freestream_Conserv_HLLC;
599f844368SJames Wright           problem->apply_freestream.qfunction_loc          = Freestream_Conserv_HLLC_loc;
609f844368SJames Wright           problem->apply_freestream_jacobian.qfunction     = Freestream_Jacobian_Conserv_HLLC;
619f844368SJames Wright           problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLLC_loc;
629f844368SJames Wright           break;
639f844368SJames Wright       }
649f844368SJames Wright       break;
659f844368SJames Wright     case STATEVAR_PRIMITIVE:
669f844368SJames Wright       switch (riemann) {
679f844368SJames Wright         case RIEMANN_HLL:
689f844368SJames Wright           problem->apply_freestream.qfunction              = Freestream_Prim_HLL;
699f844368SJames Wright           problem->apply_freestream.qfunction_loc          = Freestream_Prim_HLL_loc;
709f844368SJames Wright           problem->apply_freestream_jacobian.qfunction     = Freestream_Jacobian_Prim_HLL;
719f844368SJames Wright           problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLL_loc;
729f844368SJames Wright           break;
739f844368SJames Wright         case RIEMANN_HLLC:
749f844368SJames Wright           problem->apply_freestream.qfunction              = Freestream_Prim_HLLC;
759f844368SJames Wright           problem->apply_freestream.qfunction_loc          = Freestream_Prim_HLLC_loc;
769f844368SJames Wright           problem->apply_freestream_jacobian.qfunction     = Freestream_Jacobian_Prim_HLLC;
779f844368SJames Wright           problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLLC_loc;
789f844368SJames Wright           break;
799f844368SJames Wright       }
809f844368SJames Wright       break;
81*a2d72b6fSJames Wright     case STATEVAR_ENTROPY:
82*a2d72b6fSJames Wright       switch (riemann) {
83*a2d72b6fSJames Wright         case RIEMANN_HLL:
84*a2d72b6fSJames Wright           problem->apply_freestream.qfunction              = Freestream_Entropy_HLL;
85*a2d72b6fSJames Wright           problem->apply_freestream.qfunction_loc          = Freestream_Entropy_HLL_loc;
86*a2d72b6fSJames Wright           problem->apply_freestream_jacobian.qfunction     = Freestream_Jacobian_Entropy_HLL;
87*a2d72b6fSJames Wright           problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Entropy_HLL_loc;
88*a2d72b6fSJames Wright           break;
89*a2d72b6fSJames Wright         case RIEMANN_HLLC:
90*a2d72b6fSJames Wright           problem->apply_freestream.qfunction              = Freestream_Entropy_HLLC;
91*a2d72b6fSJames Wright           problem->apply_freestream.qfunction_loc          = Freestream_Entropy_HLLC_loc;
92*a2d72b6fSJames Wright           problem->apply_freestream_jacobian.qfunction     = Freestream_Jacobian_Entropy_HLLC;
93*a2d72b6fSJames Wright           problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Entropy_HLLC_loc;
94*a2d72b6fSJames Wright           break;
95*a2d72b6fSJames Wright       }
96*a2d72b6fSJames Wright       break;
979f844368SJames Wright   }
989f844368SJames Wright 
999f844368SJames Wright   Y_inf.pressure *= Pascal;
1009f844368SJames Wright   for (int i = 0; i < 3; i++) Y_inf.velocity[i] *= meter / second;
1019f844368SJames Wright   Y_inf.temperature *= Kelvin;
1029f844368SJames Wright 
1039f844368SJames Wright   State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf);
1049f844368SJames Wright 
1059f844368SJames Wright   // -- Set freestream_ctx struct values
1069f844368SJames Wright   PetscCall(PetscCalloc1(1, &freestream_ctx));
1079f844368SJames Wright   freestream_ctx->newtonian_ctx = *newtonian_ig_ctx;
1089f844368SJames Wright   freestream_ctx->S_infty       = S_infty;
1099f844368SJames Wright 
1109f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &freestream_context));
1119f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(freestream_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*freestream_ctx), freestream_ctx));
1129f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(freestream_context, CEED_MEM_HOST, FreeContextPetsc));
1139f844368SJames Wright   problem->apply_freestream.qfunction_context = freestream_context;
1149f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(freestream_context, &problem->apply_freestream_jacobian.qfunction_context));
1159f844368SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1169f844368SJames Wright }
1179f844368SJames Wright 
1189f844368SJames Wright static const char *const OutflowTypes[] = {"RIEMANN", "PRESSURE", "OutflowType", "OUTFLOW_", NULL};
1199f844368SJames Wright typedef enum {
1209f844368SJames Wright   OUTFLOW_RIEMANN,
1219f844368SJames Wright   OUTFLOW_PRESSURE,
1229f844368SJames Wright } OutflowType;
1239f844368SJames Wright 
124731c13d7SJames Wright PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) {
1259f844368SJames Wright   User                 user = *(User *)ctx;
1269f844368SJames Wright   Ceed                 ceed = user->ceed;
1279f844368SJames Wright   OutflowContext       outflow_ctx;
1289f844368SJames Wright   OutflowType          outflow_type = OUTFLOW_RIEMANN;
1299f844368SJames Wright   CeedQFunctionContext outflow_context;
1309f844368SJames Wright   const PetscScalar    Kelvin = user->units->Kelvin;
1319f844368SJames Wright   const PetscScalar    Pascal = user->units->Pascal;
1329f844368SJames Wright 
1339f844368SJames Wright   PetscFunctionBeginUser;
1349f844368SJames Wright   CeedScalar pressure    = reference->pressure / Pascal;
1359f844368SJames Wright   CeedScalar temperature = reference->temperature / Kelvin;
1369f844368SJames Wright   CeedScalar recirc = 1, softplus_velocity = 1e-2;
1379f844368SJames Wright   PetscOptionsBegin(user->comm, NULL, "Options for Outflow boundary condition", NULL);
1389f844368SJames Wright   PetscCall(
1399f844368SJames Wright       PetscOptionsEnum("-outflow_type", "Type of outflow condition", NULL, OutflowTypes, (PetscEnum)outflow_type, (PetscEnum *)&outflow_type, NULL));
1409f844368SJames Wright   PetscCall(PetscOptionsScalar("-outflow_pressure", "Pressure at outflow condition", NULL, pressure, &pressure, NULL));
1419f844368SJames Wright   if (outflow_type == OUTFLOW_RIEMANN) {
1429f844368SJames Wright     PetscCall(PetscOptionsScalar("-outflow_temperature", "Temperature at outflow condition", NULL, temperature, &temperature, NULL));
1439f844368SJames Wright     PetscCall(
1449f844368SJames Wright         PetscOptionsReal("-outflow_recirc", "Fraction of recirculation to allow in exterior velocity state [0,1]", NULL, recirc, &recirc, NULL));
1459f844368SJames Wright     PetscCall(PetscOptionsReal("-outflow_softplus_velocity", "Characteristic velocity of softplus regularization", NULL, softplus_velocity,
1469f844368SJames Wright                                &softplus_velocity, NULL));
1479f844368SJames Wright   }
1489f844368SJames Wright   PetscOptionsEnd();
1499f844368SJames Wright   pressure *= Pascal;
1509f844368SJames Wright   temperature *= Kelvin;
1519f844368SJames Wright 
1529f844368SJames Wright   switch (outflow_type) {
1539f844368SJames Wright     case OUTFLOW_RIEMANN:
1549f844368SJames Wright       switch (user->phys->state_var) {
1559f844368SJames Wright         case STATEVAR_CONSERVATIVE:
1569f844368SJames Wright           problem->apply_outflow.qfunction              = RiemannOutflow_Conserv;
1579f844368SJames Wright           problem->apply_outflow.qfunction_loc          = RiemannOutflow_Conserv_loc;
1589f844368SJames Wright           problem->apply_outflow_jacobian.qfunction     = RiemannOutflow_Jacobian_Conserv;
1599f844368SJames Wright           problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Conserv_loc;
1609f844368SJames Wright           break;
1619f844368SJames Wright         case STATEVAR_PRIMITIVE:
1629f844368SJames Wright           problem->apply_outflow.qfunction              = RiemannOutflow_Prim;
1639f844368SJames Wright           problem->apply_outflow.qfunction_loc          = RiemannOutflow_Prim_loc;
1649f844368SJames Wright           problem->apply_outflow_jacobian.qfunction     = RiemannOutflow_Jacobian_Prim;
1659f844368SJames Wright           problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Prim_loc;
1669f844368SJames Wright           break;
167*a2d72b6fSJames Wright         case STATEVAR_ENTROPY:
168*a2d72b6fSJames Wright           problem->apply_outflow.qfunction              = RiemannOutflow_Entropy;
169*a2d72b6fSJames Wright           problem->apply_outflow.qfunction_loc          = RiemannOutflow_Entropy_loc;
170*a2d72b6fSJames Wright           problem->apply_outflow_jacobian.qfunction     = RiemannOutflow_Jacobian_Entropy;
171*a2d72b6fSJames Wright           problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Entropy_loc;
172*a2d72b6fSJames Wright           break;
1739f844368SJames Wright       }
1749f844368SJames Wright       break;
1759f844368SJames Wright     case OUTFLOW_PRESSURE:
1769f844368SJames Wright       switch (user->phys->state_var) {
1779f844368SJames Wright         case STATEVAR_CONSERVATIVE:
1789f844368SJames Wright           problem->apply_outflow.qfunction              = PressureOutflow_Conserv;
1799f844368SJames Wright           problem->apply_outflow.qfunction_loc          = PressureOutflow_Conserv_loc;
1809f844368SJames Wright           problem->apply_outflow_jacobian.qfunction     = PressureOutflow_Jacobian_Conserv;
1819f844368SJames Wright           problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc;
1829f844368SJames Wright           break;
1839f844368SJames Wright         case STATEVAR_PRIMITIVE:
1849f844368SJames Wright           problem->apply_outflow.qfunction              = PressureOutflow_Prim;
1859f844368SJames Wright           problem->apply_outflow.qfunction_loc          = PressureOutflow_Prim_loc;
1869f844368SJames Wright           problem->apply_outflow_jacobian.qfunction     = PressureOutflow_Jacobian_Prim;
1879f844368SJames Wright           problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc;
1889f844368SJames Wright           break;
189*a2d72b6fSJames Wright         case STATEVAR_ENTROPY:
190*a2d72b6fSJames Wright           problem->apply_outflow.qfunction              = PressureOutflow_Entropy;
191*a2d72b6fSJames Wright           problem->apply_outflow.qfunction_loc          = PressureOutflow_Entropy_loc;
192*a2d72b6fSJames Wright           problem->apply_outflow_jacobian.qfunction     = PressureOutflow_Jacobian_Entropy;
193*a2d72b6fSJames Wright           problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Entropy_loc;
194*a2d72b6fSJames Wright           break;
1959f844368SJames Wright       }
1969f844368SJames Wright       break;
1979f844368SJames Wright   }
1989f844368SJames Wright   PetscCall(PetscCalloc1(1, &outflow_ctx));
1999f844368SJames Wright   outflow_ctx->gas               = *newtonian_ig_ctx;
2009f844368SJames Wright   outflow_ctx->recirc            = recirc;
2019f844368SJames Wright   outflow_ctx->softplus_velocity = softplus_velocity;
2029f844368SJames Wright   outflow_ctx->pressure          = pressure;
2039f844368SJames Wright   outflow_ctx->temperature       = temperature;
2049f844368SJames Wright 
2059f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &outflow_context));
2069f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(outflow_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*outflow_ctx), outflow_ctx));
2079f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(outflow_context, CEED_MEM_HOST, FreeContextPetsc));
2089f844368SJames Wright   problem->apply_outflow.qfunction_context = outflow_context;
2099f844368SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(outflow_context, &problem->apply_outflow_jacobian.qfunction_context));
2109f844368SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2119f844368SJames Wright }
212