xref: /honee/problems/bc_freestream.c (revision 2d898fa6800fcddb9cd7cb214e287e8576fd592f)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
39ed3d70dSJames Wright 
49ed3d70dSJames Wright /// @file
59ed3d70dSJames Wright /// Utility functions for setting up Freestream boundary condition
69ed3d70dSJames Wright 
79ed3d70dSJames Wright #include "../qfunctions/bc_freestream.h"
89ed3d70dSJames Wright 
99ed3d70dSJames Wright #include <ceed.h>
109ed3d70dSJames Wright #include <petscdm.h>
119ed3d70dSJames Wright 
12149fb536SJames Wright #include <navierstokes.h>
139ed3d70dSJames Wright #include "../qfunctions/newtonian_types.h"
149ed3d70dSJames Wright 
156dfcbb05SJames Wright static const char *const RiemannSolverTypes[] = {"HLL", "HLLC", "RiemannSolverTypes", "RIEMANN_", NULL};
169ed3d70dSJames Wright 
173e3353edSJames Wright static PetscErrorCode RiemannSolverUnitTests(NewtonianIdealGasContext gas, CeedScalar rtol);
183e3353edSJames Wright 
19d4e423e7SJames Wright typedef struct {
20d4e423e7SJames Wright   RiemannFluxType flux_type;
21d4e423e7SJames Wright } *FreestreamHoneeBCCtx;
22d4e423e7SJames Wright 
23d4e423e7SJames Wright static PetscErrorCode FreestreamBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
24d4e423e7SJames Wright   Honee         honee;
25d4e423e7SJames Wright   HoneeBCStruct honee_bc;
26d4e423e7SJames Wright 
27d4e423e7SJames Wright   PetscFunctionBeginUser;
28d4e423e7SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
29d4e423e7SJames Wright   honee                                    = honee_bc->honee;
30d4e423e7SJames Wright   FreestreamHoneeBCCtx freestream_honee_bc = (FreestreamHoneeBCCtx)honee_bc->ctx;
31d4e423e7SJames Wright 
32d4e423e7SJames Wright   switch (honee->phys->state_var) {
33d4e423e7SJames Wright     case STATEVAR_CONSERVATIVE:
34d4e423e7SJames Wright       switch (freestream_honee_bc->flux_type) {
35d4e423e7SJames Wright         case RIEMANN_HLL:
36d4e423e7SJames Wright           PetscCall(HoneeBCCreateIFunctionQF(bc_def, Freestream_Conserv_HLL, Freestream_Conserv_HLL_loc, honee_bc->qfctx, qf));
37d4e423e7SJames Wright           break;
38d4e423e7SJames Wright         case RIEMANN_HLLC:
39d4e423e7SJames Wright           PetscCall(HoneeBCCreateIFunctionQF(bc_def, Freestream_Conserv_HLLC, Freestream_Conserv_HLLC_loc, honee_bc->qfctx, qf));
40d4e423e7SJames Wright           break;
41d4e423e7SJames Wright       }
42d4e423e7SJames Wright       break;
43d4e423e7SJames Wright     case STATEVAR_PRIMITIVE:
44d4e423e7SJames Wright       switch (freestream_honee_bc->flux_type) {
45d4e423e7SJames Wright         case RIEMANN_HLL:
46d4e423e7SJames Wright           PetscCall(HoneeBCCreateIFunctionQF(bc_def, Freestream_Prim_HLL, Freestream_Prim_HLL_loc, honee_bc->qfctx, qf));
47d4e423e7SJames Wright           break;
48d4e423e7SJames Wright         case RIEMANN_HLLC:
49d4e423e7SJames Wright           PetscCall(HoneeBCCreateIFunctionQF(bc_def, Freestream_Prim_HLLC, Freestream_Prim_HLLC_loc, honee_bc->qfctx, qf));
50d4e423e7SJames Wright           break;
51d4e423e7SJames Wright       }
52d4e423e7SJames Wright       break;
53d4e423e7SJames Wright     case STATEVAR_ENTROPY:
54d4e423e7SJames Wright       switch (freestream_honee_bc->flux_type) {
55d4e423e7SJames Wright         case RIEMANN_HLL:
56d4e423e7SJames Wright           PetscCall(HoneeBCCreateIFunctionQF(bc_def, Freestream_Entropy_HLL, Freestream_Entropy_HLL_loc, honee_bc->qfctx, qf));
57d4e423e7SJames Wright           break;
58d4e423e7SJames Wright         case RIEMANN_HLLC:
59d4e423e7SJames Wright           PetscCall(HoneeBCCreateIFunctionQF(bc_def, Freestream_Entropy_HLLC, Freestream_Entropy_HLLC_loc, honee_bc->qfctx, qf));
60d4e423e7SJames Wright           break;
61d4e423e7SJames Wright       }
62d4e423e7SJames Wright       break;
63d4e423e7SJames Wright   }
64d4e423e7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
65d4e423e7SJames Wright }
66d4e423e7SJames Wright 
67d4e423e7SJames Wright static PetscErrorCode FreestreamBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) {
68d4e423e7SJames Wright   Honee         honee;
69d4e423e7SJames Wright   HoneeBCStruct honee_bc;
70d4e423e7SJames Wright 
71d4e423e7SJames Wright   PetscFunctionBeginUser;
72d4e423e7SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
73d4e423e7SJames Wright   honee                                    = honee_bc->honee;
74d4e423e7SJames Wright   FreestreamHoneeBCCtx freestream_honee_bc = (FreestreamHoneeBCCtx)honee_bc->ctx;
75d4e423e7SJames Wright 
76d4e423e7SJames Wright   switch (honee->phys->state_var) {
77d4e423e7SJames Wright     case STATEVAR_CONSERVATIVE:
78d4e423e7SJames Wright       switch (freestream_honee_bc->flux_type) {
79d4e423e7SJames Wright         case RIEMANN_HLL:
80d4e423e7SJames Wright           PetscCall(HoneeBCCreateIJacobianQF(bc_def, Freestream_Jacobian_Conserv_HLL, Freestream_Jacobian_Conserv_HLL_loc, honee_bc->qfctx, qf));
81d4e423e7SJames Wright           break;
82d4e423e7SJames Wright         case RIEMANN_HLLC:
83d4e423e7SJames Wright           PetscCall(HoneeBCCreateIJacobianQF(bc_def, Freestream_Jacobian_Conserv_HLLC, Freestream_Jacobian_Conserv_HLLC_loc, honee_bc->qfctx, qf));
84d4e423e7SJames Wright           break;
85d4e423e7SJames Wright       }
86d4e423e7SJames Wright       break;
87d4e423e7SJames Wright     case STATEVAR_PRIMITIVE:
88d4e423e7SJames Wright       switch (freestream_honee_bc->flux_type) {
89d4e423e7SJames Wright         case RIEMANN_HLL:
90d4e423e7SJames Wright           PetscCall(HoneeBCCreateIJacobianQF(bc_def, Freestream_Jacobian_Prim_HLL, Freestream_Jacobian_Prim_HLL_loc, honee_bc->qfctx, qf));
91d4e423e7SJames Wright           break;
92d4e423e7SJames Wright         case RIEMANN_HLLC:
93d4e423e7SJames Wright           PetscCall(HoneeBCCreateIJacobianQF(bc_def, Freestream_Jacobian_Prim_HLLC, Freestream_Jacobian_Prim_HLLC_loc, honee_bc->qfctx, qf));
94d4e423e7SJames Wright           break;
95d4e423e7SJames Wright       }
96d4e423e7SJames Wright       break;
97d4e423e7SJames Wright     case STATEVAR_ENTROPY:
98d4e423e7SJames Wright       switch (freestream_honee_bc->flux_type) {
99d4e423e7SJames Wright         case RIEMANN_HLL:
100d4e423e7SJames Wright           PetscCall(HoneeBCCreateIJacobianQF(bc_def, Freestream_Jacobian_Entropy_HLL, Freestream_Jacobian_Entropy_HLL_loc, honee_bc->qfctx, qf));
101d4e423e7SJames Wright           break;
102d4e423e7SJames Wright         case RIEMANN_HLLC:
103d4e423e7SJames Wright           PetscCall(HoneeBCCreateIJacobianQF(bc_def, Freestream_Jacobian_Entropy_HLLC, Freestream_Jacobian_Entropy_HLLC_loc, honee_bc->qfctx, qf));
104d4e423e7SJames Wright           break;
105d4e423e7SJames Wright       }
106d4e423e7SJames Wright       break;
107d4e423e7SJames Wright   }
108d4e423e7SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
109d4e423e7SJames Wright }
110d4e423e7SJames Wright 
111d4e423e7SJames Wright PetscErrorCode FreestreamBCSetup(BCDefinition bc_def, ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx,
112d4e423e7SJames Wright                                  const StatePrimitive *reference) {
1130c373b74SJames Wright   Honee                honee = *(Honee *)ctx;
1140c373b74SJames Wright   MPI_Comm             comm  = honee->comm;
1150c373b74SJames Wright   Ceed                 ceed  = honee->ceed;
1169ed3d70dSJames Wright   FreestreamContext    freestream_ctx;
117e07531f7SJames Wright   CeedQFunctionContext freestream_qfctx;
118d4e423e7SJames Wright   RiemannFluxType      flux_type = RIEMANN_HLLC;
119c9f37605SMohammed Amin   Units                units     = honee->units;
120d4e423e7SJames Wright   FreestreamHoneeBCCtx freestream_honee_bc;
121d4e423e7SJames Wright   HoneeBCStruct        honee_bc;
1229ed3d70dSJames Wright 
1239ed3d70dSJames Wright   PetscFunctionBeginUser;
124d4e423e7SJames Wright   // Freestream inherits reference state. We re-dimensionalize so the defaults in -help will be visible in SI units.
125c9f37605SMohammed Amin   StatePrimitive Y_inf = {.pressure = reference->pressure / units->Pascal, .velocity = {0}, .temperature = reference->temperature / units->Kelvin};
126c9f37605SMohammed Amin   for (int i = 0; i < 3; i++) Y_inf.velocity[i] = reference->velocity[i] * units->second / units->meter;
1279ed3d70dSJames Wright 
1289ed3d70dSJames Wright   PetscOptionsBegin(comm, NULL, "Options for Freestream boundary condition", NULL);
1299ed3d70dSJames Wright   PetscCall(PetscOptionsEnum("-freestream_riemann", "Riemann solver to use in freestream boundary condition", NULL, RiemannSolverTypes,
130d4e423e7SJames Wright                              (PetscEnum)flux_type, (PetscEnum *)&flux_type, NULL));
1319ed3d70dSJames Wright   PetscCall(PetscOptionsScalar("-freestream_pressure", "Pressure at freestream condition", NULL, Y_inf.pressure, &Y_inf.pressure, NULL));
1329ed3d70dSJames Wright   PetscInt narray = 3;
1339ed3d70dSJames Wright   PetscCall(PetscOptionsScalarArray("-freestream_velocity", "Velocity at freestream condition", NULL, Y_inf.velocity, &narray, NULL));
1349ed3d70dSJames Wright   PetscCall(PetscOptionsScalar("-freestream_temperature", "Temperature at freestream condition", NULL, Y_inf.temperature, &Y_inf.temperature, NULL));
1359ed3d70dSJames Wright   PetscOptionsEnd();
1369ed3d70dSJames Wright 
137c9f37605SMohammed Amin   Y_inf.pressure *= units->Pascal;
138c9f37605SMohammed Amin   for (int i = 0; i < 3; i++) Y_inf.velocity[i] *= units->meter / units->second;
139c9f37605SMohammed Amin   Y_inf.temperature *= units->Kelvin;
1409ed3d70dSJames Wright 
1419ed3d70dSJames Wright   State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf);
1429ed3d70dSJames Wright 
1439ed3d70dSJames Wright   // -- Set freestream_ctx struct values
144*2d898fa6SJames Wright   PetscCall(PetscNew(&freestream_ctx));
1459ed3d70dSJames Wright   freestream_ctx->newtonian_ctx = *newtonian_ig_ctx;
1469ed3d70dSJames Wright   freestream_ctx->S_infty       = S_infty;
1479ed3d70dSJames Wright 
1480c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &freestream_qfctx));
149e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(freestream_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*freestream_ctx), freestream_ctx));
150e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(freestream_qfctx, CEED_MEM_HOST, FreeContextPetsc));
151d4e423e7SJames Wright 
152d4e423e7SJames Wright   PetscCall(PetscNew(&honee_bc));
153d4e423e7SJames Wright   PetscCall(PetscNew(&freestream_honee_bc));
154d4e423e7SJames Wright   freestream_honee_bc->flux_type = flux_type;
155d4e423e7SJames Wright   honee_bc->ctx                  = freestream_honee_bc;
156d4e423e7SJames Wright   honee_bc->DestroyCtx           = PetscCtxDestroyDefault;
157d4e423e7SJames Wright   honee_bc->honee                = honee;
1581abc2837SJames Wright   honee_bc->num_comps_jac_data   = honee->phys->implicit ? 5 : 0;
159d4e423e7SJames Wright   honee_bc->qfctx                = freestream_qfctx;
160d4e423e7SJames Wright   PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
161d4e423e7SJames Wright 
162d4e423e7SJames Wright   PetscCall(BCDefinitionSetIFunction(bc_def, FreestreamBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
163d4e423e7SJames Wright   PetscCall(BCDefinitionSetIJacobian(bc_def, FreestreamBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp));
1643e3353edSJames Wright 
1653e3353edSJames Wright   {
1663e3353edSJames Wright     PetscBool run_unit_tests = PETSC_FALSE;
1673e3353edSJames Wright 
1683e3353edSJames Wright     PetscCall(PetscOptionsGetBool(NULL, NULL, "-riemann_solver_unit_tests", &run_unit_tests, NULL));
1693e3353edSJames Wright     if (run_unit_tests) PetscCall(RiemannSolverUnitTests(newtonian_ig_ctx, 5e-7));
1703e3353edSJames Wright   }
1719ed3d70dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1729ed3d70dSJames Wright }
1739ed3d70dSJames Wright 
174224fc8c8SJames Wright // *****************************************************************************
175224fc8c8SJames Wright // Code for verifying the Riemann solver and Riemann Jacobian functions
176224fc8c8SJames Wright // *****************************************************************************
1773e3353edSJames Wright 
1783e3353edSJames Wright // @brief Calculate relative error, (A - B) / S
1793e3353edSJames Wright // If S < threshold, then set S=1
1803e3353edSJames Wright static inline CeedScalar RelativeError(CeedScalar S, CeedScalar A, CeedScalar B, CeedScalar threshold) {
1813e3353edSJames Wright   return (A - B) / (fabs(S) > threshold ? S : 1);
1823e3353edSJames Wright }
1833e3353edSJames Wright 
1843e3353edSJames Wright // @brief Check errors of a State vector and print if above tolerance
1853e3353edSJames Wright static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
1863e3353edSJames Wright                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
1873e3353edSJames Wright   CeedScalar relative_error[5];  // relative error
1883e3353edSJames Wright   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
1893e3353edSJames Wright 
1903e3353edSJames Wright   PetscFunctionBeginUser;
1913e3353edSJames Wright   relative_error[0] = RelativeError(Q_s[0], Q_a[0], Q_b[0], divisor_threshold);
1923e3353edSJames Wright   relative_error[4] = RelativeError(Q_s[4], Q_a[4], Q_b[4], divisor_threshold);
1933e3353edSJames Wright 
1943e3353edSJames Wright   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
1953e3353edSJames Wright   for (int i = 1; i < 4; i++) {
1963e3353edSJames Wright     relative_error[i] = RelativeError(u_magnitude, Q_a[i], Q_b[i], divisor_threshold);
1973e3353edSJames Wright   }
1983e3353edSJames Wright 
1993e3353edSJames Wright   if (fabs(relative_error[0]) >= rtol_0) {
2003e3353edSJames Wright     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
2013e3353edSJames Wright   }
2023e3353edSJames Wright   for (int i = 1; i < 4; i++) {
2033e3353edSJames Wright     if (fabs(relative_error[i]) >= rtol_u) {
2043e3353edSJames Wright       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
2053e3353edSJames Wright     }
2063e3353edSJames Wright   }
2073e3353edSJames Wright   if (fabs(relative_error[4]) >= rtol_4) {
2083e3353edSJames Wright     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
2093e3353edSJames Wright   }
2103e3353edSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2113e3353edSJames Wright }
2123e3353edSJames Wright 
2133e3353edSJames Wright // @brief Verify RiemannFlux_HLL_fwd function against finite-difference approximation
2143e3353edSJames Wright static PetscErrorCode TestRiemannHLL_fwd(NewtonianIdealGasContext gas, CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
2153e3353edSJames Wright   CeedScalar       eps = 4e-7;  // Finite difference step
2163e3353edSJames Wright   char             buf[128];
2173e3353edSJames Wright   const CeedScalar T           = 200;
2183e3353edSJames Wright   const CeedScalar rho         = 1.2;
2193e3353edSJames Wright   const CeedScalar p           = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
2203e3353edSJames Wright   const CeedScalar u_base      = 40;
2213e3353edSJames Wright   const CeedScalar u[3]        = {u_base, u_base * 1.1, u_base * 1.2};
2223e3353edSJames Wright   const CeedScalar Y0_left[5]  = {p, u[0], u[1], u[2], T};
2233e3353edSJames Wright   const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T};
2243e3353edSJames Wright   CeedScalar       normal[3]   = {1, 2, 3};
2253e3353edSJames Wright 
2263e3353edSJames Wright   PetscFunctionBeginUser;
2273e3353edSJames Wright   State left0  = StateFromY(gas, Y0_left);
2283e3353edSJames Wright   State right0 = StateFromY(gas, Y0_right);
22964667825SJames Wright   ScaleN(normal, 1 / Norm3(normal), 3);
2303e3353edSJames Wright 
2313e3353edSJames Wright   for (int i = 0; i < 10; i++) {
2323e3353edSJames Wright     CeedScalar dFlux[5] = {0.}, dFlux_fd[5] = {0.};
2333e3353edSJames Wright     {  // Calculate dFlux using *_fwd function
2343e3353edSJames Wright       CeedScalar dY_right[5] = {0};
2353e3353edSJames Wright       CeedScalar dY_left[5]  = {0};
2363e3353edSJames Wright 
2373e3353edSJames Wright       if (i < 5) {
2383e3353edSJames Wright         dY_left[i] = Y0_left[i];
2393e3353edSJames Wright       } else {
2403e3353edSJames Wright         dY_right[i % 5] = Y0_right[i % 5];
2413e3353edSJames Wright       }
2423e3353edSJames Wright       State dleft0  = StateFromY_fwd(gas, left0, dY_left);
2433e3353edSJames Wright       State dright0 = StateFromY_fwd(gas, right0, dY_right);
2443e3353edSJames Wright 
2453e3353edSJames Wright       StateConservative dFlux_state = RiemannFlux_HLL_fwd(gas, left0, dleft0, right0, dright0, normal);
2463e3353edSJames Wright       UnpackState_U(dFlux_state, dFlux);
2473e3353edSJames Wright     }
2483e3353edSJames Wright 
2493e3353edSJames Wright     {  // Calculate dFlux_fd via finite difference approximation
2503e3353edSJames Wright       CeedScalar Y1_left[5]  = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]};
2513e3353edSJames Wright       CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]};
2523e3353edSJames Wright       CeedScalar Flux0[5], Flux1[5];
2533e3353edSJames Wright 
2543e3353edSJames Wright       if (i < 5) {
2553e3353edSJames Wright         Y1_left[i] *= 1 + eps;
2563e3353edSJames Wright       } else {
2573e3353edSJames Wright         Y1_right[i % 5] *= 1 + eps;
2583e3353edSJames Wright       }
2593e3353edSJames Wright       State left1  = StateFromY(gas, Y1_left);
2603e3353edSJames Wright       State right1 = StateFromY(gas, Y1_right);
2613e3353edSJames Wright 
2623e3353edSJames Wright       StateConservative Flux0_state = RiemannFlux_HLL(gas, left0, right0, normal);
2633e3353edSJames Wright       StateConservative Flux1_state = RiemannFlux_HLL(gas, left1, right1, normal);
2643e3353edSJames Wright       UnpackState_U(Flux0_state, Flux0);
2653e3353edSJames Wright       UnpackState_U(Flux1_state, Flux1);
2663e3353edSJames Wright       for (int j = 0; j < 5; j++) dFlux_fd[j] = (Flux1[j] - Flux0[j]) / eps;
2673e3353edSJames Wright     }
2683e3353edSJames Wright 
2693e3353edSJames Wright     snprintf(buf, sizeof buf, "RiemannFlux_HLL i=%d: Flux", i);
2703e3353edSJames Wright     PetscCall(CheckQWithTolerance(dFlux_fd, dFlux, dFlux_fd, buf, rtol_0, rtol_u, rtol_4));
2713e3353edSJames Wright   }
2723e3353edSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2733e3353edSJames Wright }
2743e3353edSJames Wright 
2753e3353edSJames Wright // @brief Verify RiemannFlux_HLLC_fwd function against finite-difference approximation
2763e3353edSJames Wright static PetscErrorCode TestRiemannHLLC_fwd(NewtonianIdealGasContext gas, CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
2773e3353edSJames Wright   CeedScalar       eps = 4e-7;  // Finite difference step
2783e3353edSJames Wright   char             buf[128];
2793e3353edSJames Wright   const CeedScalar T           = 200;
2803e3353edSJames Wright   const CeedScalar rho         = 1.2;
2813e3353edSJames Wright   const CeedScalar p           = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
2823e3353edSJames Wright   const CeedScalar u_base      = 40;
2833e3353edSJames Wright   const CeedScalar u[3]        = {u_base, u_base * 1.1, u_base * 1.2};
2843e3353edSJames Wright   const CeedScalar Y0_left[5]  = {p, u[0], u[1], u[2], T};
2853e3353edSJames Wright   const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T};
2863e3353edSJames Wright   CeedScalar       normal[3]   = {1, 2, 3};
2873e3353edSJames Wright 
2883e3353edSJames Wright   PetscFunctionBeginUser;
2893e3353edSJames Wright   State left0  = StateFromY(gas, Y0_left);
2903e3353edSJames Wright   State right0 = StateFromY(gas, Y0_right);
29164667825SJames Wright   ScaleN(normal, 1 / Norm3(normal), 3);
2923e3353edSJames Wright 
2933e3353edSJames Wright   for (int i = 0; i < 10; i++) {
2943e3353edSJames Wright     CeedScalar dFlux[5] = {0.}, dFlux_fd[5] = {0.};
2953e3353edSJames Wright     {  // Calculate dFlux using *_fwd function
2963e3353edSJames Wright       CeedScalar dY_right[5] = {0};
2973e3353edSJames Wright       CeedScalar dY_left[5]  = {0};
2983e3353edSJames Wright 
2993e3353edSJames Wright       if (i < 5) {
3003e3353edSJames Wright         dY_left[i] = Y0_left[i];
3013e3353edSJames Wright       } else {
3023e3353edSJames Wright         dY_right[i % 5] = Y0_right[i % 5];
3033e3353edSJames Wright       }
3043e3353edSJames Wright       State dleft0  = StateFromY_fwd(gas, left0, dY_left);
3053e3353edSJames Wright       State dright0 = StateFromY_fwd(gas, right0, dY_right);
3063e3353edSJames Wright 
3073e3353edSJames Wright       StateConservative dFlux_state = RiemannFlux_HLLC_fwd(gas, left0, dleft0, right0, dright0, normal);
3083e3353edSJames Wright       UnpackState_U(dFlux_state, dFlux);
3093e3353edSJames Wright     }
3103e3353edSJames Wright 
3113e3353edSJames Wright     {  // Calculate dFlux_fd via finite difference approximation
3123e3353edSJames Wright       CeedScalar Y1_left[5]  = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]};
3133e3353edSJames Wright       CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]};
3143e3353edSJames Wright       CeedScalar Flux0[5], Flux1[5];
3153e3353edSJames Wright 
3163e3353edSJames Wright       if (i < 5) {
3173e3353edSJames Wright         Y1_left[i] *= 1 + eps;
3183e3353edSJames Wright       } else {
3193e3353edSJames Wright         Y1_right[i % 5] *= 1 + eps;
3203e3353edSJames Wright       }
3213e3353edSJames Wright       State left1  = StateFromY(gas, Y1_left);
3223e3353edSJames Wright       State right1 = StateFromY(gas, Y1_right);
3233e3353edSJames Wright 
3243e3353edSJames Wright       StateConservative Flux0_state = RiemannFlux_HLLC(gas, left0, right0, normal);
3253e3353edSJames Wright       StateConservative Flux1_state = RiemannFlux_HLLC(gas, left1, right1, normal);
3263e3353edSJames Wright       UnpackState_U(Flux0_state, Flux0);
3273e3353edSJames Wright       UnpackState_U(Flux1_state, Flux1);
3283e3353edSJames Wright       for (int j = 0; j < 5; j++) dFlux_fd[j] = (Flux1[j] - Flux0[j]) / eps;
3293e3353edSJames Wright     }
3303e3353edSJames Wright 
3313e3353edSJames Wright     snprintf(buf, sizeof buf, "RiemannFlux_HLLC i=%d: Flux", i);
3323e3353edSJames Wright     PetscCall(CheckQWithTolerance(dFlux_fd, dFlux, dFlux_fd, buf, rtol_0, rtol_u, rtol_4));
3333e3353edSJames Wright   }
3343e3353edSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3353e3353edSJames Wright }
3363e3353edSJames Wright 
3373e3353edSJames Wright // @brief Verify ComputeHLLSpeeds_Roe_fwd function against finite-difference approximation
3383e3353edSJames Wright static PetscErrorCode TestComputeHLLSpeeds_Roe_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) {
3393e3353edSJames Wright   CeedScalar       eps = 4e-7;  // Finite difference step
3403e3353edSJames Wright   char             buf[128];
3413e3353edSJames Wright   const CeedScalar T           = 200;
3423e3353edSJames Wright   const CeedScalar rho         = 1.2;
3433e3353edSJames Wright   const CeedScalar p           = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
3443e3353edSJames Wright   const CeedScalar u_base      = 40;
3453e3353edSJames Wright   const CeedScalar u[3]        = {u_base, u_base * 1.1, u_base * 1.2};
3463e3353edSJames Wright   const CeedScalar Y0_left[5]  = {p, u[0], u[1], u[2], T};
3473e3353edSJames Wright   const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T};
3483e3353edSJames Wright   CeedScalar       normal[3]   = {1, 2, 3};
3493e3353edSJames Wright 
3503e3353edSJames Wright   PetscFunctionBeginUser;
3513e3353edSJames Wright   State left0  = StateFromY(gas, Y0_left);
3523e3353edSJames Wright   State right0 = StateFromY(gas, Y0_right);
35364667825SJames Wright   ScaleN(normal, 1 / Norm3(normal), 3);
3543e3353edSJames Wright   CeedScalar u_left0  = Dot3(left0.Y.velocity, normal);
3553e3353edSJames Wright   CeedScalar u_right0 = Dot3(right0.Y.velocity, normal);
3563e3353edSJames Wright 
3573e3353edSJames Wright   for (int i = 0; i < 10; i++) {
3583e3353edSJames Wright     CeedScalar ds_left, ds_right, ds_left_fd, ds_right_fd;
3593e3353edSJames Wright     {  // Calculate ds_{left,right} using *_fwd function
3603e3353edSJames Wright       CeedScalar dY_right[5] = {0};
3613e3353edSJames Wright       CeedScalar dY_left[5]  = {0};
3623e3353edSJames Wright 
3633e3353edSJames Wright       if (i < 5) {
3643e3353edSJames Wright         dY_left[i] = Y0_left[i];
3653e3353edSJames Wright       } else {
3663e3353edSJames Wright         dY_right[i % 5] = Y0_right[i % 5];
3673e3353edSJames Wright       }
3683e3353edSJames Wright       State      dleft0   = StateFromY_fwd(gas, left0, dY_left);
3693e3353edSJames Wright       State      dright0  = StateFromY_fwd(gas, right0, dY_right);
3703e3353edSJames Wright       CeedScalar du_left  = Dot3(dleft0.Y.velocity, normal);
3713e3353edSJames Wright       CeedScalar du_right = Dot3(dright0.Y.velocity, normal);
3723e3353edSJames Wright 
3733e3353edSJames Wright       CeedScalar s_left, s_right;  // Throw away
3743e3353edSJames 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);
3753e3353edSJames Wright     }
3763e3353edSJames Wright 
3773e3353edSJames Wright     {  // Calculate ds_{left,right}_fd via finite difference approximation
3783e3353edSJames Wright       CeedScalar Y1_left[5]  = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]};
3793e3353edSJames Wright       CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]};
3803e3353edSJames Wright 
3813e3353edSJames Wright       if (i < 5) {
3823e3353edSJames Wright         Y1_left[i] *= 1 + eps;
3833e3353edSJames Wright       } else {
3843e3353edSJames Wright         Y1_right[i % 5] *= 1 + eps;
3853e3353edSJames Wright       }
3863e3353edSJames Wright       State      left1    = StateFromY(gas, Y1_left);
3873e3353edSJames Wright       State      right1   = StateFromY(gas, Y1_right);
3883e3353edSJames Wright       CeedScalar u_left1  = Dot3(left1.Y.velocity, normal);
3893e3353edSJames Wright       CeedScalar u_right1 = Dot3(right1.Y.velocity, normal);
3903e3353edSJames Wright 
3913e3353edSJames Wright       CeedScalar s_left0, s_right0, s_left1, s_right1;
3923e3353edSJames Wright       ComputeHLLSpeeds_Roe(gas, left0, u_left0, right0, u_right0, &s_left0, &s_right0);
3933e3353edSJames Wright       ComputeHLLSpeeds_Roe(gas, left1, u_left1, right1, u_right1, &s_left1, &s_right1);
3943e3353edSJames Wright       ds_left_fd  = (s_left1 - s_left0) / eps;
3953e3353edSJames Wright       ds_right_fd = (s_right1 - s_right0) / eps;
3963e3353edSJames Wright     }
3973e3353edSJames Wright 
3983e3353edSJames Wright     snprintf(buf, sizeof buf, "ComputeHLLSpeeds_Roe i=%d:", i);
3993e3353edSJames Wright     {
4003e3353edSJames Wright       CeedScalar divisor_threshold = 10 * CEED_EPSILON;
4013e3353edSJames Wright       CeedScalar ds_left_err, ds_right_err;
4023e3353edSJames Wright 
4033e3353edSJames Wright       ds_left_err  = RelativeError(ds_left_fd, ds_left, ds_left_fd, divisor_threshold);
4043e3353edSJames Wright       ds_right_err = RelativeError(ds_right_fd, ds_right, ds_right_fd, divisor_threshold);
4053e3353edSJames 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);
4063e3353edSJames 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);
4073e3353edSJames Wright     }
4083e3353edSJames Wright   }
4093e3353edSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4103e3353edSJames Wright }
4113e3353edSJames Wright 
4123e3353edSJames Wright // @brief Verify TotalSpecificEnthalpy_fwd function against finite-difference approximation
4133e3353edSJames Wright static PetscErrorCode TestTotalSpecificEnthalpy_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) {
4143e3353edSJames Wright   CeedScalar       eps = 4e-7;  // Finite difference step
4153e3353edSJames Wright   char             buf[128];
4163e3353edSJames Wright   const CeedScalar T      = 200;
4173e3353edSJames Wright   const CeedScalar rho    = 1.2;
4183e3353edSJames Wright   const CeedScalar p      = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
4193e3353edSJames Wright   const CeedScalar u_base = 40;
4203e3353edSJames Wright   const CeedScalar u[3]   = {u_base, u_base * 1.1, u_base * 1.2};
4213e3353edSJames Wright   const CeedScalar Y0[5]  = {p, u[0], u[1], u[2], T};
4223e3353edSJames Wright 
4233e3353edSJames Wright   PetscFunctionBeginUser;
4243e3353edSJames Wright   State state0 = StateFromY(gas, Y0);
4253e3353edSJames Wright 
4263e3353edSJames Wright   for (int i = 0; i < 5; i++) {
4273e3353edSJames Wright     CeedScalar dH, dH_fd;
4283e3353edSJames Wright     {  // Calculate dH using *_fwd function
4293e3353edSJames Wright       CeedScalar dY[5] = {0};
4303e3353edSJames Wright 
4313e3353edSJames Wright       dY[i]         = Y0[i];
4323e3353edSJames Wright       State dstate0 = StateFromY_fwd(gas, state0, dY);
4333e3353edSJames Wright       dH            = TotalSpecificEnthalpy_fwd(gas, state0, dstate0);
4343e3353edSJames Wright     }
4353e3353edSJames Wright 
4363e3353edSJames Wright     {  // Calculate dH_fd via finite difference approximation
4373e3353edSJames Wright       CeedScalar H0, H1;
4383e3353edSJames Wright       CeedScalar Y1[5] = {Y0[0], Y0[1], Y0[2], Y0[3], Y0[4]};
4393e3353edSJames Wright       Y1[i] *= 1 + eps;
4403e3353edSJames Wright       State state1 = StateFromY(gas, Y1);
4413e3353edSJames Wright 
4423e3353edSJames Wright       H0    = TotalSpecificEnthalpy(gas, state0);
4433e3353edSJames Wright       H1    = TotalSpecificEnthalpy(gas, state1);
4443e3353edSJames Wright       dH_fd = (H1 - H0) / eps;
4453e3353edSJames Wright     }
4463e3353edSJames Wright 
4473e3353edSJames Wright     snprintf(buf, sizeof buf, "TotalSpecificEnthalpy i=%d:", i);
4483e3353edSJames Wright     {
4493e3353edSJames Wright       CeedScalar divisor_threshold = 10 * CEED_EPSILON;
4503e3353edSJames Wright       CeedScalar dH_err;
4513e3353edSJames Wright 
4523e3353edSJames Wright       dH_err = RelativeError(dH_fd, dH, dH_fd, divisor_threshold);
4533e3353edSJames Wright       if (fabs(dH_err) >= rtol) printf("%s dH error %g (expected %.10e, got %.10e)\n", buf, dH_err, dH_fd, dH);
4543e3353edSJames Wright     }
4553e3353edSJames Wright   }
4563e3353edSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4573e3353edSJames Wright }
4583e3353edSJames Wright 
4593e3353edSJames Wright // @brief Verify RoeSetup_fwd function against finite-difference approximation
4603e3353edSJames Wright static PetscErrorCode TestRowSetup_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) {
4613e3353edSJames Wright   CeedScalar       eps = 4e-7;  // Finite difference step
4623e3353edSJames Wright   char             buf[128];
4633e3353edSJames Wright   const CeedScalar rho0[2] = {1.2, 1.4};
4643e3353edSJames Wright 
4653e3353edSJames Wright   PetscFunctionBeginUser;
4663e3353edSJames Wright   for (int i = 0; i < 2; i++) {
4673e3353edSJames Wright     RoeWeights dR, dR_fd;
4683e3353edSJames Wright     {  // Calculate using *_fwd function
4693e3353edSJames Wright       CeedScalar drho[5] = {0};
4703e3353edSJames Wright 
4713e3353edSJames Wright       drho[i] = rho0[i];
4723e3353edSJames Wright       dR      = RoeSetup_fwd(rho0[0], rho0[1], drho[0], drho[1]);
4733e3353edSJames Wright     }
4743e3353edSJames Wright 
4753e3353edSJames Wright     {  // Calculate via finite difference approximation
4763e3353edSJames Wright       RoeWeights R0, R1;
4773e3353edSJames Wright       CeedScalar rho1[5] = {rho0[0], rho0[1]};
4783e3353edSJames Wright       rho1[i] *= 1 + eps;
4793e3353edSJames Wright 
4803e3353edSJames Wright       R0          = RoeSetup(rho0[0], rho0[1]);
4813e3353edSJames Wright       R1          = RoeSetup(rho1[0], rho1[1]);
4823e3353edSJames Wright       dR_fd.left  = (R1.left - R0.left) / eps;
4833e3353edSJames Wright       dR_fd.right = (R1.right - R0.right) / eps;
4843e3353edSJames Wright     }
4853e3353edSJames Wright 
4863e3353edSJames Wright     snprintf(buf, sizeof buf, "RoeSetup i=%d:", i);
4873e3353edSJames Wright     {
4883e3353edSJames Wright       CeedScalar divisor_threshold = 10 * CEED_EPSILON;
4893e3353edSJames Wright       RoeWeights dR_err;
4903e3353edSJames Wright 
4913e3353edSJames Wright       dR_err.left  = RelativeError(dR_fd.left, dR.left, dR_fd.left, divisor_threshold);
4923e3353edSJames Wright       dR_err.right = RelativeError(dR_fd.right, dR.right, dR_fd.right, divisor_threshold);
4933e3353edSJames 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);
4943e3353edSJames 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);
4953e3353edSJames Wright     }
4963e3353edSJames Wright   }
4973e3353edSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4983e3353edSJames Wright }
4993e3353edSJames Wright 
5003e3353edSJames Wright // @brief Test Riemann solver related `*_fwd` functions via finite-difference approximation
5013e3353edSJames Wright static PetscErrorCode RiemannSolverUnitTests(NewtonianIdealGasContext gas, CeedScalar rtol) {
5023e3353edSJames Wright   PetscFunctionBeginUser;
5033e3353edSJames Wright   PetscCall(TestRiemannHLL_fwd(gas, rtol, rtol, rtol));
5043e3353edSJames Wright   PetscCall(TestRiemannHLLC_fwd(gas, rtol, rtol, rtol));
5053e3353edSJames Wright   PetscCall(TestComputeHLLSpeeds_Roe_fwd(gas, rtol));
5063e3353edSJames Wright   PetscCall(TestTotalSpecificEnthalpy_fwd(gas, rtol));
5073e3353edSJames Wright   PetscCall(TestRowSetup_fwd(gas, rtol));
5083e3353edSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5093e3353edSJames Wright }
510