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