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 196dfcbb05SJames Wright static const char *const RiemannSolverTypes[] = {"HLL", "HLLC", "RiemannSolverTypes", "RIEMANN_", NULL}; 209ed3d70dSJames Wright 21*3e3353edSJames Wright static PetscErrorCode RiemannSolverUnitTests(NewtonianIdealGasContext gas, CeedScalar rtol); 22*3e3353edSJames Wright 23991aef52SJames Wright PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) { 249ed3d70dSJames Wright User user = *(User *)ctx; 259ed3d70dSJames Wright MPI_Comm comm = user->comm; 269ed3d70dSJames Wright Ceed ceed = user->ceed; 279ed3d70dSJames Wright FreestreamContext freestream_ctx; 289ed3d70dSJames Wright CeedQFunctionContext freestream_context; 299ed3d70dSJames Wright RiemannFluxType riemann = RIEMANN_HLLC; 309ed3d70dSJames Wright PetscScalar meter = user->units->meter; 319ed3d70dSJames Wright PetscScalar second = user->units->second; 329ed3d70dSJames Wright PetscScalar Kelvin = user->units->Kelvin; 339ed3d70dSJames Wright PetscScalar Pascal = user->units->Pascal; 349ed3d70dSJames Wright 359ed3d70dSJames Wright PetscFunctionBeginUser; 369ed3d70dSJames Wright // Freestream inherits reference state. We re-dimensionalize so the defaults 379ed3d70dSJames Wright // in -help will be visible in SI units. 389ed3d70dSJames Wright StatePrimitive Y_inf = {.pressure = reference->pressure / Pascal, .velocity = {0}, .temperature = reference->temperature / Kelvin}; 399ed3d70dSJames Wright for (int i = 0; i < 3; i++) Y_inf.velocity[i] = reference->velocity[i] * second / meter; 409ed3d70dSJames Wright 419ed3d70dSJames Wright PetscOptionsBegin(comm, NULL, "Options for Freestream boundary condition", NULL); 429ed3d70dSJames Wright PetscCall(PetscOptionsEnum("-freestream_riemann", "Riemann solver to use in freestream boundary condition", NULL, RiemannSolverTypes, 439ed3d70dSJames Wright (PetscEnum)riemann, (PetscEnum *)&riemann, NULL)); 449ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-freestream_pressure", "Pressure at freestream condition", NULL, Y_inf.pressure, &Y_inf.pressure, NULL)); 459ed3d70dSJames Wright PetscInt narray = 3; 469ed3d70dSJames Wright PetscCall(PetscOptionsScalarArray("-freestream_velocity", "Velocity at freestream condition", NULL, Y_inf.velocity, &narray, NULL)); 479ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-freestream_temperature", "Temperature at freestream condition", NULL, Y_inf.temperature, &Y_inf.temperature, NULL)); 489ed3d70dSJames Wright PetscOptionsEnd(); 499ed3d70dSJames Wright 509ed3d70dSJames Wright switch (user->phys->state_var) { 519ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 529ed3d70dSJames Wright switch (riemann) { 539ed3d70dSJames Wright case RIEMANN_HLL: 549ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Conserv_HLL; 559ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Conserv_HLL_loc; 569ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Conserv_HLL; 579ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLL_loc; 589ed3d70dSJames Wright break; 599ed3d70dSJames Wright case RIEMANN_HLLC: 609ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Conserv_HLLC; 619ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Conserv_HLLC_loc; 629ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Conserv_HLLC; 639ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Conserv_HLLC_loc; 649ed3d70dSJames Wright break; 659ed3d70dSJames Wright } 669ed3d70dSJames Wright break; 679ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 689ed3d70dSJames Wright switch (riemann) { 699ed3d70dSJames Wright case RIEMANN_HLL: 709ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Prim_HLL; 719ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Prim_HLL_loc; 729ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Prim_HLL; 739ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLL_loc; 749ed3d70dSJames Wright break; 759ed3d70dSJames Wright case RIEMANN_HLLC: 769ed3d70dSJames Wright problem->apply_freestream.qfunction = Freestream_Prim_HLLC; 779ed3d70dSJames Wright problem->apply_freestream.qfunction_loc = Freestream_Prim_HLLC_loc; 789ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Prim_HLLC; 799ed3d70dSJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Prim_HLLC_loc; 809ed3d70dSJames Wright break; 819ed3d70dSJames Wright } 829ed3d70dSJames Wright break; 839b103f75SJames Wright case STATEVAR_ENTROPY: 849b103f75SJames Wright switch (riemann) { 859b103f75SJames Wright case RIEMANN_HLL: 869b103f75SJames Wright problem->apply_freestream.qfunction = Freestream_Entropy_HLL; 879b103f75SJames Wright problem->apply_freestream.qfunction_loc = Freestream_Entropy_HLL_loc; 889b103f75SJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Entropy_HLL; 899b103f75SJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Entropy_HLL_loc; 909b103f75SJames Wright break; 919b103f75SJames Wright case RIEMANN_HLLC: 929b103f75SJames Wright problem->apply_freestream.qfunction = Freestream_Entropy_HLLC; 939b103f75SJames Wright problem->apply_freestream.qfunction_loc = Freestream_Entropy_HLLC_loc; 949b103f75SJames Wright problem->apply_freestream_jacobian.qfunction = Freestream_Jacobian_Entropy_HLLC; 959b103f75SJames Wright problem->apply_freestream_jacobian.qfunction_loc = Freestream_Jacobian_Entropy_HLLC_loc; 969b103f75SJames Wright break; 979b103f75SJames Wright } 989b103f75SJames Wright break; 999ed3d70dSJames Wright } 1009ed3d70dSJames Wright 1019ed3d70dSJames Wright Y_inf.pressure *= Pascal; 1029ed3d70dSJames Wright for (int i = 0; i < 3; i++) Y_inf.velocity[i] *= meter / second; 1039ed3d70dSJames Wright Y_inf.temperature *= Kelvin; 1049ed3d70dSJames Wright 1059ed3d70dSJames Wright State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf); 1069ed3d70dSJames Wright 1079ed3d70dSJames Wright // -- Set freestream_ctx struct values 1089ed3d70dSJames Wright PetscCall(PetscCalloc1(1, &freestream_ctx)); 1099ed3d70dSJames Wright freestream_ctx->newtonian_ctx = *newtonian_ig_ctx; 1109ed3d70dSJames Wright freestream_ctx->S_infty = S_infty; 1119ed3d70dSJames Wright 1129ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &freestream_context)); 1139ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(freestream_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*freestream_ctx), freestream_ctx)); 1149ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(freestream_context, CEED_MEM_HOST, FreeContextPetsc)); 1159ed3d70dSJames Wright problem->apply_freestream.qfunction_context = freestream_context; 1169ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(freestream_context, &problem->apply_freestream_jacobian.qfunction_context)); 117*3e3353edSJames Wright 118*3e3353edSJames Wright { 119*3e3353edSJames Wright PetscBool run_unit_tests = PETSC_FALSE; 120*3e3353edSJames Wright 121*3e3353edSJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-riemann_solver_unit_tests", &run_unit_tests, NULL)); 122*3e3353edSJames Wright if (run_unit_tests) PetscCall(RiemannSolverUnitTests(newtonian_ig_ctx, 5e-7)); 123*3e3353edSJames Wright } 1249ed3d70dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1259ed3d70dSJames Wright } 1269ed3d70dSJames Wright 1279ed3d70dSJames Wright static const char *const OutflowTypes[] = {"RIEMANN", "PRESSURE", "OutflowType", "OUTFLOW_", NULL}; 1289ed3d70dSJames Wright typedef enum { 1299ed3d70dSJames Wright OUTFLOW_RIEMANN, 1309ed3d70dSJames Wright OUTFLOW_PRESSURE, 1319ed3d70dSJames Wright } OutflowType; 1329ed3d70dSJames Wright 133991aef52SJames Wright PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) { 1349ed3d70dSJames Wright User user = *(User *)ctx; 1359ed3d70dSJames Wright Ceed ceed = user->ceed; 1369ed3d70dSJames Wright OutflowContext outflow_ctx; 1379ed3d70dSJames Wright OutflowType outflow_type = OUTFLOW_RIEMANN; 1389ed3d70dSJames Wright CeedQFunctionContext outflow_context; 1399ed3d70dSJames Wright const PetscScalar Kelvin = user->units->Kelvin; 1409ed3d70dSJames Wright const PetscScalar Pascal = user->units->Pascal; 1419ed3d70dSJames Wright 1429ed3d70dSJames Wright PetscFunctionBeginUser; 1439ed3d70dSJames Wright CeedScalar pressure = reference->pressure / Pascal; 1449ed3d70dSJames Wright CeedScalar temperature = reference->temperature / Kelvin; 1459ed3d70dSJames Wright CeedScalar recirc = 1, softplus_velocity = 1e-2; 1469ed3d70dSJames Wright PetscOptionsBegin(user->comm, NULL, "Options for Outflow boundary condition", NULL); 1479ed3d70dSJames Wright PetscCall( 1489ed3d70dSJames Wright PetscOptionsEnum("-outflow_type", "Type of outflow condition", NULL, OutflowTypes, (PetscEnum)outflow_type, (PetscEnum *)&outflow_type, NULL)); 1499ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-outflow_pressure", "Pressure at outflow condition", NULL, pressure, &pressure, NULL)); 1509ed3d70dSJames Wright if (outflow_type == OUTFLOW_RIEMANN) { 1519ed3d70dSJames Wright PetscCall(PetscOptionsScalar("-outflow_temperature", "Temperature at outflow condition", NULL, temperature, &temperature, NULL)); 1529ed3d70dSJames Wright PetscCall( 1539ed3d70dSJames Wright PetscOptionsReal("-outflow_recirc", "Fraction of recirculation to allow in exterior velocity state [0,1]", NULL, recirc, &recirc, NULL)); 1549ed3d70dSJames Wright PetscCall(PetscOptionsReal("-outflow_softplus_velocity", "Characteristic velocity of softplus regularization", NULL, softplus_velocity, 1559ed3d70dSJames Wright &softplus_velocity, NULL)); 1569ed3d70dSJames Wright } 1579ed3d70dSJames Wright PetscOptionsEnd(); 1589ed3d70dSJames Wright pressure *= Pascal; 1599ed3d70dSJames Wright temperature *= Kelvin; 1609ed3d70dSJames Wright 1619ed3d70dSJames Wright switch (outflow_type) { 1629ed3d70dSJames Wright case OUTFLOW_RIEMANN: 1639ed3d70dSJames Wright switch (user->phys->state_var) { 1649ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 1659ed3d70dSJames Wright problem->apply_outflow.qfunction = RiemannOutflow_Conserv; 1669ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = RiemannOutflow_Conserv_loc; 1679ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = RiemannOutflow_Jacobian_Conserv; 1689ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Conserv_loc; 1699ed3d70dSJames Wright break; 1709ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 1719ed3d70dSJames Wright problem->apply_outflow.qfunction = RiemannOutflow_Prim; 1729ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = RiemannOutflow_Prim_loc; 1739ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = RiemannOutflow_Jacobian_Prim; 1749ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Prim_loc; 1759ed3d70dSJames Wright break; 1769b103f75SJames Wright case STATEVAR_ENTROPY: 1779b103f75SJames Wright problem->apply_outflow.qfunction = RiemannOutflow_Entropy; 1789b103f75SJames Wright problem->apply_outflow.qfunction_loc = RiemannOutflow_Entropy_loc; 1799b103f75SJames Wright problem->apply_outflow_jacobian.qfunction = RiemannOutflow_Jacobian_Entropy; 1809b103f75SJames Wright problem->apply_outflow_jacobian.qfunction_loc = RiemannOutflow_Jacobian_Entropy_loc; 1819b103f75SJames Wright break; 1829ed3d70dSJames Wright } 1839ed3d70dSJames Wright break; 1849ed3d70dSJames Wright case OUTFLOW_PRESSURE: 1859ed3d70dSJames Wright switch (user->phys->state_var) { 1869ed3d70dSJames Wright case STATEVAR_CONSERVATIVE: 1879ed3d70dSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Conserv; 1889ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Conserv_loc; 1899ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Conserv; 1909ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Conserv_loc; 1919ed3d70dSJames Wright break; 1929ed3d70dSJames Wright case STATEVAR_PRIMITIVE: 1939ed3d70dSJames Wright problem->apply_outflow.qfunction = PressureOutflow_Prim; 1949ed3d70dSJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Prim_loc; 1959ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Prim; 1969ed3d70dSJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Prim_loc; 1979ed3d70dSJames Wright break; 1989b103f75SJames Wright case STATEVAR_ENTROPY: 1999b103f75SJames Wright problem->apply_outflow.qfunction = PressureOutflow_Entropy; 2009b103f75SJames Wright problem->apply_outflow.qfunction_loc = PressureOutflow_Entropy_loc; 2019b103f75SJames Wright problem->apply_outflow_jacobian.qfunction = PressureOutflow_Jacobian_Entropy; 2029b103f75SJames Wright problem->apply_outflow_jacobian.qfunction_loc = PressureOutflow_Jacobian_Entropy_loc; 2039b103f75SJames Wright break; 2049ed3d70dSJames Wright } 2059ed3d70dSJames Wright break; 2069ed3d70dSJames Wright } 2079ed3d70dSJames Wright PetscCall(PetscCalloc1(1, &outflow_ctx)); 2089ed3d70dSJames Wright outflow_ctx->gas = *newtonian_ig_ctx; 2099ed3d70dSJames Wright outflow_ctx->recirc = recirc; 2109ed3d70dSJames Wright outflow_ctx->softplus_velocity = softplus_velocity; 2119ed3d70dSJames Wright outflow_ctx->pressure = pressure; 2129ed3d70dSJames Wright outflow_ctx->temperature = temperature; 2139ed3d70dSJames Wright 2149ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &outflow_context)); 2159ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(outflow_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*outflow_ctx), outflow_ctx)); 2169ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(outflow_context, CEED_MEM_HOST, FreeContextPetsc)); 2179ed3d70dSJames Wright problem->apply_outflow.qfunction_context = outflow_context; 2189ed3d70dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(outflow_context, &problem->apply_outflow_jacobian.qfunction_context)); 2199ed3d70dSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2209ed3d70dSJames Wright } 221*3e3353edSJames Wright 222*3e3353edSJames Wright // @brief Calculate relative error, (A - B) / S 223*3e3353edSJames Wright // If S < threshold, then set S=1 224*3e3353edSJames Wright static inline CeedScalar RelativeError(CeedScalar S, CeedScalar A, CeedScalar B, CeedScalar threshold) { 225*3e3353edSJames Wright return (A - B) / (fabs(S) > threshold ? S : 1); 226*3e3353edSJames Wright } 227*3e3353edSJames Wright 228*3e3353edSJames Wright // @brief Check errors of a State vector and print if above tolerance 229*3e3353edSJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name, 230*3e3353edSJames Wright PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) { 231*3e3353edSJames Wright CeedScalar relative_error[5]; // relative error 232*3e3353edSJames Wright CeedScalar divisor_threshold = 10 * CEED_EPSILON; 233*3e3353edSJames Wright 234*3e3353edSJames Wright PetscFunctionBeginUser; 235*3e3353edSJames Wright relative_error[0] = RelativeError(Q_s[0], Q_a[0], Q_b[0], divisor_threshold); 236*3e3353edSJames Wright relative_error[4] = RelativeError(Q_s[4], Q_a[4], Q_b[4], divisor_threshold); 237*3e3353edSJames Wright 238*3e3353edSJames Wright CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3])); 239*3e3353edSJames Wright for (int i = 1; i < 4; i++) { 240*3e3353edSJames Wright relative_error[i] = RelativeError(u_magnitude, Q_a[i], Q_b[i], divisor_threshold); 241*3e3353edSJames Wright } 242*3e3353edSJames Wright 243*3e3353edSJames Wright if (fabs(relative_error[0]) >= rtol_0) { 244*3e3353edSJames Wright printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]); 245*3e3353edSJames Wright } 246*3e3353edSJames Wright for (int i = 1; i < 4; i++) { 247*3e3353edSJames Wright if (fabs(relative_error[i]) >= rtol_u) { 248*3e3353edSJames Wright printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]); 249*3e3353edSJames Wright } 250*3e3353edSJames Wright } 251*3e3353edSJames Wright if (fabs(relative_error[4]) >= rtol_4) { 252*3e3353edSJames Wright printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]); 253*3e3353edSJames Wright } 254*3e3353edSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 255*3e3353edSJames Wright } 256*3e3353edSJames Wright 257*3e3353edSJames Wright // @brief Verify RiemannFlux_HLL_fwd function against finite-difference approximation 258*3e3353edSJames Wright static PetscErrorCode TestRiemannHLL_fwd(NewtonianIdealGasContext gas, CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 259*3e3353edSJames Wright CeedScalar eps = 4e-7; // Finite difference step 260*3e3353edSJames Wright char buf[128]; 261*3e3353edSJames Wright const CeedScalar T = 200; 262*3e3353edSJames Wright const CeedScalar rho = 1.2; 263*3e3353edSJames Wright const CeedScalar p = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T; 264*3e3353edSJames Wright const CeedScalar u_base = 40; 265*3e3353edSJames Wright const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 266*3e3353edSJames Wright const CeedScalar Y0_left[5] = {p, u[0], u[1], u[2], T}; 267*3e3353edSJames Wright const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T}; 268*3e3353edSJames Wright CeedScalar normal[3] = {1, 2, 3}; 269*3e3353edSJames Wright 270*3e3353edSJames Wright PetscFunctionBeginUser; 271*3e3353edSJames Wright State left0 = StateFromY(gas, Y0_left); 272*3e3353edSJames Wright State right0 = StateFromY(gas, Y0_right); 273*3e3353edSJames Wright ScaleN(normal, 1 / sqrt(Dot3(normal, normal)), 3); 274*3e3353edSJames Wright 275*3e3353edSJames Wright for (int i = 0; i < 10; i++) { 276*3e3353edSJames Wright CeedScalar dFlux[5] = {0.}, dFlux_fd[5] = {0.}; 277*3e3353edSJames Wright { // Calculate dFlux using *_fwd function 278*3e3353edSJames Wright CeedScalar dY_right[5] = {0}; 279*3e3353edSJames Wright CeedScalar dY_left[5] = {0}; 280*3e3353edSJames Wright 281*3e3353edSJames Wright if (i < 5) { 282*3e3353edSJames Wright dY_left[i] = Y0_left[i]; 283*3e3353edSJames Wright } else { 284*3e3353edSJames Wright dY_right[i % 5] = Y0_right[i % 5]; 285*3e3353edSJames Wright } 286*3e3353edSJames Wright State dleft0 = StateFromY_fwd(gas, left0, dY_left); 287*3e3353edSJames Wright State dright0 = StateFromY_fwd(gas, right0, dY_right); 288*3e3353edSJames Wright 289*3e3353edSJames Wright StateConservative dFlux_state = RiemannFlux_HLL_fwd(gas, left0, dleft0, right0, dright0, normal); 290*3e3353edSJames Wright UnpackState_U(dFlux_state, dFlux); 291*3e3353edSJames Wright } 292*3e3353edSJames Wright 293*3e3353edSJames Wright { // Calculate dFlux_fd via finite difference approximation 294*3e3353edSJames Wright CeedScalar Y1_left[5] = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]}; 295*3e3353edSJames Wright CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]}; 296*3e3353edSJames Wright CeedScalar Flux0[5], Flux1[5]; 297*3e3353edSJames Wright 298*3e3353edSJames Wright if (i < 5) { 299*3e3353edSJames Wright Y1_left[i] *= 1 + eps; 300*3e3353edSJames Wright } else { 301*3e3353edSJames Wright Y1_right[i % 5] *= 1 + eps; 302*3e3353edSJames Wright } 303*3e3353edSJames Wright State left1 = StateFromY(gas, Y1_left); 304*3e3353edSJames Wright State right1 = StateFromY(gas, Y1_right); 305*3e3353edSJames Wright 306*3e3353edSJames Wright StateConservative Flux0_state = RiemannFlux_HLL(gas, left0, right0, normal); 307*3e3353edSJames Wright StateConservative Flux1_state = RiemannFlux_HLL(gas, left1, right1, normal); 308*3e3353edSJames Wright UnpackState_U(Flux0_state, Flux0); 309*3e3353edSJames Wright UnpackState_U(Flux1_state, Flux1); 310*3e3353edSJames Wright for (int j = 0; j < 5; j++) dFlux_fd[j] = (Flux1[j] - Flux0[j]) / eps; 311*3e3353edSJames Wright } 312*3e3353edSJames Wright 313*3e3353edSJames Wright snprintf(buf, sizeof buf, "RiemannFlux_HLL i=%d: Flux", i); 314*3e3353edSJames Wright PetscCall(CheckQWithTolerance(dFlux_fd, dFlux, dFlux_fd, buf, rtol_0, rtol_u, rtol_4)); 315*3e3353edSJames Wright } 316*3e3353edSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 317*3e3353edSJames Wright } 318*3e3353edSJames Wright 319*3e3353edSJames Wright // @brief Verify RiemannFlux_HLLC_fwd function against finite-difference approximation 320*3e3353edSJames Wright static PetscErrorCode TestRiemannHLLC_fwd(NewtonianIdealGasContext gas, CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) { 321*3e3353edSJames Wright CeedScalar eps = 4e-7; // Finite difference step 322*3e3353edSJames Wright char buf[128]; 323*3e3353edSJames Wright const CeedScalar T = 200; 324*3e3353edSJames Wright const CeedScalar rho = 1.2; 325*3e3353edSJames Wright const CeedScalar p = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T; 326*3e3353edSJames Wright const CeedScalar u_base = 40; 327*3e3353edSJames Wright const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 328*3e3353edSJames Wright const CeedScalar Y0_left[5] = {p, u[0], u[1], u[2], T}; 329*3e3353edSJames Wright const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T}; 330*3e3353edSJames Wright CeedScalar normal[3] = {1, 2, 3}; 331*3e3353edSJames Wright 332*3e3353edSJames Wright PetscFunctionBeginUser; 333*3e3353edSJames Wright State left0 = StateFromY(gas, Y0_left); 334*3e3353edSJames Wright State right0 = StateFromY(gas, Y0_right); 335*3e3353edSJames Wright ScaleN(normal, 1 / sqrt(Dot3(normal, normal)), 3); 336*3e3353edSJames Wright 337*3e3353edSJames Wright for (int i = 0; i < 10; i++) { 338*3e3353edSJames Wright CeedScalar dFlux[5] = {0.}, dFlux_fd[5] = {0.}; 339*3e3353edSJames Wright { // Calculate dFlux using *_fwd function 340*3e3353edSJames Wright CeedScalar dY_right[5] = {0}; 341*3e3353edSJames Wright CeedScalar dY_left[5] = {0}; 342*3e3353edSJames Wright 343*3e3353edSJames Wright if (i < 5) { 344*3e3353edSJames Wright dY_left[i] = Y0_left[i]; 345*3e3353edSJames Wright } else { 346*3e3353edSJames Wright dY_right[i % 5] = Y0_right[i % 5]; 347*3e3353edSJames Wright } 348*3e3353edSJames Wright State dleft0 = StateFromY_fwd(gas, left0, dY_left); 349*3e3353edSJames Wright State dright0 = StateFromY_fwd(gas, right0, dY_right); 350*3e3353edSJames Wright 351*3e3353edSJames Wright StateConservative dFlux_state = RiemannFlux_HLLC_fwd(gas, left0, dleft0, right0, dright0, normal); 352*3e3353edSJames Wright UnpackState_U(dFlux_state, dFlux); 353*3e3353edSJames Wright } 354*3e3353edSJames Wright 355*3e3353edSJames Wright { // Calculate dFlux_fd via finite difference approximation 356*3e3353edSJames Wright CeedScalar Y1_left[5] = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]}; 357*3e3353edSJames Wright CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]}; 358*3e3353edSJames Wright CeedScalar Flux0[5], Flux1[5]; 359*3e3353edSJames Wright 360*3e3353edSJames Wright if (i < 5) { 361*3e3353edSJames Wright Y1_left[i] *= 1 + eps; 362*3e3353edSJames Wright } else { 363*3e3353edSJames Wright Y1_right[i % 5] *= 1 + eps; 364*3e3353edSJames Wright } 365*3e3353edSJames Wright State left1 = StateFromY(gas, Y1_left); 366*3e3353edSJames Wright State right1 = StateFromY(gas, Y1_right); 367*3e3353edSJames Wright 368*3e3353edSJames Wright StateConservative Flux0_state = RiemannFlux_HLLC(gas, left0, right0, normal); 369*3e3353edSJames Wright StateConservative Flux1_state = RiemannFlux_HLLC(gas, left1, right1, normal); 370*3e3353edSJames Wright UnpackState_U(Flux0_state, Flux0); 371*3e3353edSJames Wright UnpackState_U(Flux1_state, Flux1); 372*3e3353edSJames Wright for (int j = 0; j < 5; j++) dFlux_fd[j] = (Flux1[j] - Flux0[j]) / eps; 373*3e3353edSJames Wright } 374*3e3353edSJames Wright 375*3e3353edSJames Wright snprintf(buf, sizeof buf, "RiemannFlux_HLLC i=%d: Flux", i); 376*3e3353edSJames Wright PetscCall(CheckQWithTolerance(dFlux_fd, dFlux, dFlux_fd, buf, rtol_0, rtol_u, rtol_4)); 377*3e3353edSJames Wright } 378*3e3353edSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 379*3e3353edSJames Wright } 380*3e3353edSJames Wright 381*3e3353edSJames Wright // @brief Verify ComputeHLLSpeeds_Roe_fwd function against finite-difference approximation 382*3e3353edSJames Wright static PetscErrorCode TestComputeHLLSpeeds_Roe_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) { 383*3e3353edSJames Wright CeedScalar eps = 4e-7; // Finite difference step 384*3e3353edSJames Wright char buf[128]; 385*3e3353edSJames Wright const CeedScalar T = 200; 386*3e3353edSJames Wright const CeedScalar rho = 1.2; 387*3e3353edSJames Wright const CeedScalar p = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T; 388*3e3353edSJames Wright const CeedScalar u_base = 40; 389*3e3353edSJames Wright const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 390*3e3353edSJames Wright const CeedScalar Y0_left[5] = {p, u[0], u[1], u[2], T}; 391*3e3353edSJames Wright const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T}; 392*3e3353edSJames Wright CeedScalar normal[3] = {1, 2, 3}; 393*3e3353edSJames Wright 394*3e3353edSJames Wright PetscFunctionBeginUser; 395*3e3353edSJames Wright State left0 = StateFromY(gas, Y0_left); 396*3e3353edSJames Wright State right0 = StateFromY(gas, Y0_right); 397*3e3353edSJames Wright ScaleN(normal, 1 / sqrt(Dot3(normal, normal)), 3); 398*3e3353edSJames Wright CeedScalar u_left0 = Dot3(left0.Y.velocity, normal); 399*3e3353edSJames Wright CeedScalar u_right0 = Dot3(right0.Y.velocity, normal); 400*3e3353edSJames Wright 401*3e3353edSJames Wright for (int i = 0; i < 10; i++) { 402*3e3353edSJames Wright CeedScalar ds_left, ds_right, ds_left_fd, ds_right_fd; 403*3e3353edSJames Wright { // Calculate ds_{left,right} using *_fwd function 404*3e3353edSJames Wright CeedScalar dY_right[5] = {0}; 405*3e3353edSJames Wright CeedScalar dY_left[5] = {0}; 406*3e3353edSJames Wright 407*3e3353edSJames Wright if (i < 5) { 408*3e3353edSJames Wright dY_left[i] = Y0_left[i]; 409*3e3353edSJames Wright } else { 410*3e3353edSJames Wright dY_right[i % 5] = Y0_right[i % 5]; 411*3e3353edSJames Wright } 412*3e3353edSJames Wright State dleft0 = StateFromY_fwd(gas, left0, dY_left); 413*3e3353edSJames Wright State dright0 = StateFromY_fwd(gas, right0, dY_right); 414*3e3353edSJames Wright CeedScalar du_left = Dot3(dleft0.Y.velocity, normal); 415*3e3353edSJames Wright CeedScalar du_right = Dot3(dright0.Y.velocity, normal); 416*3e3353edSJames Wright 417*3e3353edSJames Wright CeedScalar s_left, s_right; // Throw away 418*3e3353edSJames Wright ComputeHLLSpeeds_Roe_fwd(gas, left0, dleft0, u_left0, du_left, right0, dright0, u_right0, du_right, &s_left, &ds_left, &s_right, &ds_right); 419*3e3353edSJames Wright } 420*3e3353edSJames Wright 421*3e3353edSJames Wright { // Calculate ds_{left,right}_fd via finite difference approximation 422*3e3353edSJames Wright CeedScalar Y1_left[5] = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]}; 423*3e3353edSJames Wright CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]}; 424*3e3353edSJames Wright 425*3e3353edSJames Wright if (i < 5) { 426*3e3353edSJames Wright Y1_left[i] *= 1 + eps; 427*3e3353edSJames Wright } else { 428*3e3353edSJames Wright Y1_right[i % 5] *= 1 + eps; 429*3e3353edSJames Wright } 430*3e3353edSJames Wright State left1 = StateFromY(gas, Y1_left); 431*3e3353edSJames Wright State right1 = StateFromY(gas, Y1_right); 432*3e3353edSJames Wright CeedScalar u_left1 = Dot3(left1.Y.velocity, normal); 433*3e3353edSJames Wright CeedScalar u_right1 = Dot3(right1.Y.velocity, normal); 434*3e3353edSJames Wright 435*3e3353edSJames Wright CeedScalar s_left0, s_right0, s_left1, s_right1; 436*3e3353edSJames Wright ComputeHLLSpeeds_Roe(gas, left0, u_left0, right0, u_right0, &s_left0, &s_right0); 437*3e3353edSJames Wright ComputeHLLSpeeds_Roe(gas, left1, u_left1, right1, u_right1, &s_left1, &s_right1); 438*3e3353edSJames Wright ds_left_fd = (s_left1 - s_left0) / eps; 439*3e3353edSJames Wright ds_right_fd = (s_right1 - s_right0) / eps; 440*3e3353edSJames Wright } 441*3e3353edSJames Wright 442*3e3353edSJames Wright snprintf(buf, sizeof buf, "ComputeHLLSpeeds_Roe i=%d:", i); 443*3e3353edSJames Wright { 444*3e3353edSJames Wright CeedScalar divisor_threshold = 10 * CEED_EPSILON; 445*3e3353edSJames Wright CeedScalar ds_left_err, ds_right_err; 446*3e3353edSJames Wright 447*3e3353edSJames Wright ds_left_err = RelativeError(ds_left_fd, ds_left, ds_left_fd, divisor_threshold); 448*3e3353edSJames Wright ds_right_err = RelativeError(ds_right_fd, ds_right, ds_right_fd, divisor_threshold); 449*3e3353edSJames Wright if (fabs(ds_left_err) >= rtol) printf("%s ds_left error %g (expected %.10e, got %.10e)\n", buf, ds_left_err, ds_left_fd, ds_left); 450*3e3353edSJames Wright if (fabs(ds_right_err) >= rtol) printf("%s ds_right error %g (expected %.10e, got %.10e)\n", buf, ds_right_err, ds_right_fd, ds_right); 451*3e3353edSJames Wright } 452*3e3353edSJames Wright } 453*3e3353edSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 454*3e3353edSJames Wright } 455*3e3353edSJames Wright 456*3e3353edSJames Wright // @brief Verify TotalSpecificEnthalpy_fwd function against finite-difference approximation 457*3e3353edSJames Wright static PetscErrorCode TestTotalSpecificEnthalpy_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) { 458*3e3353edSJames Wright CeedScalar eps = 4e-7; // Finite difference step 459*3e3353edSJames Wright char buf[128]; 460*3e3353edSJames Wright const CeedScalar T = 200; 461*3e3353edSJames Wright const CeedScalar rho = 1.2; 462*3e3353edSJames Wright const CeedScalar p = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T; 463*3e3353edSJames Wright const CeedScalar u_base = 40; 464*3e3353edSJames Wright const CeedScalar u[3] = {u_base, u_base * 1.1, u_base * 1.2}; 465*3e3353edSJames Wright const CeedScalar Y0[5] = {p, u[0], u[1], u[2], T}; 466*3e3353edSJames Wright 467*3e3353edSJames Wright PetscFunctionBeginUser; 468*3e3353edSJames Wright State state0 = StateFromY(gas, Y0); 469*3e3353edSJames Wright 470*3e3353edSJames Wright for (int i = 0; i < 5; i++) { 471*3e3353edSJames Wright CeedScalar dH, dH_fd; 472*3e3353edSJames Wright { // Calculate dH using *_fwd function 473*3e3353edSJames Wright CeedScalar dY[5] = {0}; 474*3e3353edSJames Wright 475*3e3353edSJames Wright dY[i] = Y0[i]; 476*3e3353edSJames Wright State dstate0 = StateFromY_fwd(gas, state0, dY); 477*3e3353edSJames Wright dH = TotalSpecificEnthalpy_fwd(gas, state0, dstate0); 478*3e3353edSJames Wright } 479*3e3353edSJames Wright 480*3e3353edSJames Wright { // Calculate dH_fd via finite difference approximation 481*3e3353edSJames Wright CeedScalar H0, H1; 482*3e3353edSJames Wright CeedScalar Y1[5] = {Y0[0], Y0[1], Y0[2], Y0[3], Y0[4]}; 483*3e3353edSJames Wright Y1[i] *= 1 + eps; 484*3e3353edSJames Wright State state1 = StateFromY(gas, Y1); 485*3e3353edSJames Wright 486*3e3353edSJames Wright H0 = TotalSpecificEnthalpy(gas, state0); 487*3e3353edSJames Wright H1 = TotalSpecificEnthalpy(gas, state1); 488*3e3353edSJames Wright dH_fd = (H1 - H0) / eps; 489*3e3353edSJames Wright } 490*3e3353edSJames Wright 491*3e3353edSJames Wright snprintf(buf, sizeof buf, "TotalSpecificEnthalpy i=%d:", i); 492*3e3353edSJames Wright { 493*3e3353edSJames Wright CeedScalar divisor_threshold = 10 * CEED_EPSILON; 494*3e3353edSJames Wright CeedScalar dH_err; 495*3e3353edSJames Wright 496*3e3353edSJames Wright dH_err = RelativeError(dH_fd, dH, dH_fd, divisor_threshold); 497*3e3353edSJames Wright if (fabs(dH_err) >= rtol) printf("%s dH error %g (expected %.10e, got %.10e)\n", buf, dH_err, dH_fd, dH); 498*3e3353edSJames Wright } 499*3e3353edSJames Wright } 500*3e3353edSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 501*3e3353edSJames Wright } 502*3e3353edSJames Wright 503*3e3353edSJames Wright // @brief Verify RoeSetup_fwd function against finite-difference approximation 504*3e3353edSJames Wright static PetscErrorCode TestRowSetup_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) { 505*3e3353edSJames Wright CeedScalar eps = 4e-7; // Finite difference step 506*3e3353edSJames Wright char buf[128]; 507*3e3353edSJames Wright const CeedScalar rho0[2] = {1.2, 1.4}; 508*3e3353edSJames Wright 509*3e3353edSJames Wright PetscFunctionBeginUser; 510*3e3353edSJames Wright for (int i = 0; i < 2; i++) { 511*3e3353edSJames Wright RoeWeights dR, dR_fd; 512*3e3353edSJames Wright { // Calculate using *_fwd function 513*3e3353edSJames Wright CeedScalar drho[5] = {0}; 514*3e3353edSJames Wright 515*3e3353edSJames Wright drho[i] = rho0[i]; 516*3e3353edSJames Wright dR = RoeSetup_fwd(rho0[0], rho0[1], drho[0], drho[1]); 517*3e3353edSJames Wright } 518*3e3353edSJames Wright 519*3e3353edSJames Wright { // Calculate via finite difference approximation 520*3e3353edSJames Wright RoeWeights R0, R1; 521*3e3353edSJames Wright CeedScalar rho1[5] = {rho0[0], rho0[1]}; 522*3e3353edSJames Wright rho1[i] *= 1 + eps; 523*3e3353edSJames Wright 524*3e3353edSJames Wright R0 = RoeSetup(rho0[0], rho0[1]); 525*3e3353edSJames Wright R1 = RoeSetup(rho1[0], rho1[1]); 526*3e3353edSJames Wright dR_fd.left = (R1.left - R0.left) / eps; 527*3e3353edSJames Wright dR_fd.right = (R1.right - R0.right) / eps; 528*3e3353edSJames Wright } 529*3e3353edSJames Wright 530*3e3353edSJames Wright snprintf(buf, sizeof buf, "RoeSetup i=%d:", i); 531*3e3353edSJames Wright { 532*3e3353edSJames Wright CeedScalar divisor_threshold = 10 * CEED_EPSILON; 533*3e3353edSJames Wright RoeWeights dR_err; 534*3e3353edSJames Wright 535*3e3353edSJames Wright dR_err.left = RelativeError(dR_fd.left, dR.left, dR_fd.left, divisor_threshold); 536*3e3353edSJames Wright dR_err.right = RelativeError(dR_fd.right, dR.right, dR_fd.right, divisor_threshold); 537*3e3353edSJames Wright if (fabs(dR_err.left) >= rtol) printf("%s dR.left error %g (expected %.10e, got %.10e)\n", buf, dR_err.left, dR_fd.left, dR.left); 538*3e3353edSJames Wright if (fabs(dR_err.right) >= rtol) printf("%s dR.right error %g (expected %.10e, got %.10e)\n", buf, dR_err.right, dR_fd.right, dR.right); 539*3e3353edSJames Wright } 540*3e3353edSJames Wright } 541*3e3353edSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 542*3e3353edSJames Wright } 543*3e3353edSJames Wright 544*3e3353edSJames Wright // @brief Test Riemann solver related `*_fwd` functions via finite-difference approximation 545*3e3353edSJames Wright static PetscErrorCode RiemannSolverUnitTests(NewtonianIdealGasContext gas, CeedScalar rtol) { 546*3e3353edSJames Wright PetscFunctionBeginUser; 547*3e3353edSJames Wright PetscCall(TestRiemannHLL_fwd(gas, rtol, rtol, rtol)); 548*3e3353edSJames Wright PetscCall(TestRiemannHLLC_fwd(gas, rtol, rtol, rtol)); 549*3e3353edSJames Wright PetscCall(TestComputeHLLSpeeds_Roe_fwd(gas, rtol)); 550*3e3353edSJames Wright PetscCall(TestTotalSpecificEnthalpy_fwd(gas, rtol)); 551*3e3353edSJames Wright PetscCall(TestRowSetup_fwd(gas, rtol)); 552*3e3353edSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 553*3e3353edSJames Wright } 554