xref: /honee/problems/bc_freestream.c (revision 43f5052dfd2ac7a767d1919bf17877b95c0f9699)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Utility functions for setting up Freestream boundary condition
6 
7 #include "../qfunctions/bc_freestream.h"
8 
9 #include <ceed.h>
10 #include <petscdm.h>
11 
12 #include <navierstokes.h>
13 #include "../qfunctions/newtonian_types.h"
14 
15 static const char *const RiemannSolverTypes[] = {"HLL", "HLLC", "RiemannSolverTypes", "RIEMANN_", NULL};
16 
17 static PetscErrorCode RiemannSolverUnitTests(NewtonianIdealGasContext gas, CeedScalar rtol);
18 
19 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference) {
20   Honee                honee = *(Honee *)ctx;
21   MPI_Comm             comm  = honee->comm;
22   Ceed                 ceed  = honee->ceed;
23   FreestreamContext    freestream_ctx;
24   CeedQFunctionContext freestream_qfctx;
25   RiemannFluxType      riemann = RIEMANN_HLLC;
26   PetscScalar          meter   = honee->units->meter;
27   PetscScalar          second  = honee->units->second;
28   PetscScalar          Kelvin  = honee->units->Kelvin;
29   PetscScalar          Pascal  = honee->units->Pascal;
30 
31   PetscFunctionBeginUser;
32   // Freestream inherits reference state. We re-dimensionalize so the defaults
33   // in -help will be visible in SI units.
34   StatePrimitive Y_inf = {.pressure = reference->pressure / Pascal, .velocity = {0}, .temperature = reference->temperature / Kelvin};
35   for (int i = 0; i < 3; i++) Y_inf.velocity[i] = reference->velocity[i] * second / meter;
36 
37   PetscOptionsBegin(comm, NULL, "Options for Freestream boundary condition", NULL);
38   PetscCall(PetscOptionsEnum("-freestream_riemann", "Riemann solver to use in freestream boundary condition", NULL, RiemannSolverTypes,
39                              (PetscEnum)riemann, (PetscEnum *)&riemann, NULL));
40   PetscCall(PetscOptionsScalar("-freestream_pressure", "Pressure at freestream condition", NULL, Y_inf.pressure, &Y_inf.pressure, NULL));
41   PetscInt narray = 3;
42   PetscCall(PetscOptionsScalarArray("-freestream_velocity", "Velocity at freestream condition", NULL, Y_inf.velocity, &narray, NULL));
43   PetscCall(PetscOptionsScalar("-freestream_temperature", "Temperature at freestream condition", NULL, Y_inf.temperature, &Y_inf.temperature, NULL));
44   PetscOptionsEnd();
45 
46   switch (honee->phys->state_var) {
47     case STATEVAR_CONSERVATIVE:
48       switch (riemann) {
49         case RIEMANN_HLL:
50           problem->apply_freestream.qf_func_ptr          = Freestream_Conserv_HLL;
51           problem->apply_freestream.qf_loc               = Freestream_Conserv_HLL_loc;
52           problem->apply_freestream_jacobian.qf_func_ptr = Freestream_Jacobian_Conserv_HLL;
53           problem->apply_freestream_jacobian.qf_loc      = Freestream_Jacobian_Conserv_HLL_loc;
54           break;
55         case RIEMANN_HLLC:
56           problem->apply_freestream.qf_func_ptr          = Freestream_Conserv_HLLC;
57           problem->apply_freestream.qf_loc               = Freestream_Conserv_HLLC_loc;
58           problem->apply_freestream_jacobian.qf_func_ptr = Freestream_Jacobian_Conserv_HLLC;
59           problem->apply_freestream_jacobian.qf_loc      = Freestream_Jacobian_Conserv_HLLC_loc;
60           break;
61       }
62       break;
63     case STATEVAR_PRIMITIVE:
64       switch (riemann) {
65         case RIEMANN_HLL:
66           problem->apply_freestream.qf_func_ptr          = Freestream_Prim_HLL;
67           problem->apply_freestream.qf_loc               = Freestream_Prim_HLL_loc;
68           problem->apply_freestream_jacobian.qf_func_ptr = Freestream_Jacobian_Prim_HLL;
69           problem->apply_freestream_jacobian.qf_loc      = Freestream_Jacobian_Prim_HLL_loc;
70           break;
71         case RIEMANN_HLLC:
72           problem->apply_freestream.qf_func_ptr          = Freestream_Prim_HLLC;
73           problem->apply_freestream.qf_loc               = Freestream_Prim_HLLC_loc;
74           problem->apply_freestream_jacobian.qf_func_ptr = Freestream_Jacobian_Prim_HLLC;
75           problem->apply_freestream_jacobian.qf_loc      = Freestream_Jacobian_Prim_HLLC_loc;
76           break;
77       }
78       break;
79     case STATEVAR_ENTROPY:
80       switch (riemann) {
81         case RIEMANN_HLL:
82           problem->apply_freestream.qf_func_ptr          = Freestream_Entropy_HLL;
83           problem->apply_freestream.qf_loc               = Freestream_Entropy_HLL_loc;
84           problem->apply_freestream_jacobian.qf_func_ptr = Freestream_Jacobian_Entropy_HLL;
85           problem->apply_freestream_jacobian.qf_loc      = Freestream_Jacobian_Entropy_HLL_loc;
86           break;
87         case RIEMANN_HLLC:
88           problem->apply_freestream.qf_func_ptr          = Freestream_Entropy_HLLC;
89           problem->apply_freestream.qf_loc               = Freestream_Entropy_HLLC_loc;
90           problem->apply_freestream_jacobian.qf_func_ptr = Freestream_Jacobian_Entropy_HLLC;
91           problem->apply_freestream_jacobian.qf_loc      = Freestream_Jacobian_Entropy_HLLC_loc;
92           break;
93       }
94       break;
95   }
96 
97   Y_inf.pressure *= Pascal;
98   for (int i = 0; i < 3; i++) Y_inf.velocity[i] *= meter / second;
99   Y_inf.temperature *= Kelvin;
100 
101   State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf);
102 
103   // -- Set freestream_ctx struct values
104   PetscCall(PetscCalloc1(1, &freestream_ctx));
105   freestream_ctx->newtonian_ctx = *newtonian_ig_ctx;
106   freestream_ctx->S_infty       = S_infty;
107 
108   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &freestream_qfctx));
109   PetscCallCeed(ceed, CeedQFunctionContextSetData(freestream_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*freestream_ctx), freestream_ctx));
110   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(freestream_qfctx, CEED_MEM_HOST, FreeContextPetsc));
111   problem->apply_freestream.qfctx = freestream_qfctx;
112   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(freestream_qfctx, &problem->apply_freestream_jacobian.qfctx));
113 
114   {
115     PetscBool run_unit_tests = PETSC_FALSE;
116 
117     PetscCall(PetscOptionsGetBool(NULL, NULL, "-riemann_solver_unit_tests", &run_unit_tests, NULL));
118     if (run_unit_tests) PetscCall(RiemannSolverUnitTests(newtonian_ig_ctx, 5e-7));
119   }
120   PetscFunctionReturn(PETSC_SUCCESS);
121 }
122 
123 // *****************************************************************************
124 // Code for verifying the Riemann solver and Riemann Jacobian functions
125 // *****************************************************************************
126 
127 // @brief Calculate relative error, (A - B) / S
128 // If S < threshold, then set S=1
129 static inline CeedScalar RelativeError(CeedScalar S, CeedScalar A, CeedScalar B, CeedScalar threshold) {
130   return (A - B) / (fabs(S) > threshold ? S : 1);
131 }
132 
133 // @brief Check errors of a State vector and print if above tolerance
134 static PetscErrorCode CheckQWithTolerance(const CeedScalar Q_s[5], const CeedScalar Q_a[5], const CeedScalar Q_b[5], const char *name,
135                                           PetscReal rtol_0, PetscReal rtol_u, PetscReal rtol_4) {
136   CeedScalar relative_error[5];  // relative error
137   CeedScalar divisor_threshold = 10 * CEED_EPSILON;
138 
139   PetscFunctionBeginUser;
140   relative_error[0] = RelativeError(Q_s[0], Q_a[0], Q_b[0], divisor_threshold);
141   relative_error[4] = RelativeError(Q_s[4], Q_a[4], Q_b[4], divisor_threshold);
142 
143   CeedScalar u_magnitude = sqrt(Square(Q_s[1]) + Square(Q_s[2]) + Square(Q_s[3]));
144   for (int i = 1; i < 4; i++) {
145     relative_error[i] = RelativeError(u_magnitude, Q_a[i], Q_b[i], divisor_threshold);
146   }
147 
148   if (fabs(relative_error[0]) >= rtol_0) {
149     printf("%s[0] error %g (expected %.10e, got %.10e)\n", name, relative_error[0], Q_s[0], Q_a[0]);
150   }
151   for (int i = 1; i < 4; i++) {
152     if (fabs(relative_error[i]) >= rtol_u) {
153       printf("%s[%d] error %g (expected %.10e, got %.10e)\n", name, i, relative_error[i], Q_s[i], Q_a[i]);
154     }
155   }
156   if (fabs(relative_error[4]) >= rtol_4) {
157     printf("%s[4] error %g (expected %.10e, got %.10e)\n", name, relative_error[4], Q_s[4], Q_a[4]);
158   }
159   PetscFunctionReturn(PETSC_SUCCESS);
160 }
161 
162 // @brief Verify RiemannFlux_HLL_fwd function against finite-difference approximation
163 static PetscErrorCode TestRiemannHLL_fwd(NewtonianIdealGasContext gas, CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
164   CeedScalar       eps = 4e-7;  // Finite difference step
165   char             buf[128];
166   const CeedScalar T           = 200;
167   const CeedScalar rho         = 1.2;
168   const CeedScalar p           = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
169   const CeedScalar u_base      = 40;
170   const CeedScalar u[3]        = {u_base, u_base * 1.1, u_base * 1.2};
171   const CeedScalar Y0_left[5]  = {p, u[0], u[1], u[2], T};
172   const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T};
173   CeedScalar       normal[3]   = {1, 2, 3};
174 
175   PetscFunctionBeginUser;
176   State left0  = StateFromY(gas, Y0_left);
177   State right0 = StateFromY(gas, Y0_right);
178   ScaleN(normal, 1 / Norm3(normal), 3);
179 
180   for (int i = 0; i < 10; i++) {
181     CeedScalar dFlux[5] = {0.}, dFlux_fd[5] = {0.};
182     {  // Calculate dFlux using *_fwd function
183       CeedScalar dY_right[5] = {0};
184       CeedScalar dY_left[5]  = {0};
185 
186       if (i < 5) {
187         dY_left[i] = Y0_left[i];
188       } else {
189         dY_right[i % 5] = Y0_right[i % 5];
190       }
191       State dleft0  = StateFromY_fwd(gas, left0, dY_left);
192       State dright0 = StateFromY_fwd(gas, right0, dY_right);
193 
194       StateConservative dFlux_state = RiemannFlux_HLL_fwd(gas, left0, dleft0, right0, dright0, normal);
195       UnpackState_U(dFlux_state, dFlux);
196     }
197 
198     {  // Calculate dFlux_fd via finite difference approximation
199       CeedScalar Y1_left[5]  = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]};
200       CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]};
201       CeedScalar Flux0[5], Flux1[5];
202 
203       if (i < 5) {
204         Y1_left[i] *= 1 + eps;
205       } else {
206         Y1_right[i % 5] *= 1 + eps;
207       }
208       State left1  = StateFromY(gas, Y1_left);
209       State right1 = StateFromY(gas, Y1_right);
210 
211       StateConservative Flux0_state = RiemannFlux_HLL(gas, left0, right0, normal);
212       StateConservative Flux1_state = RiemannFlux_HLL(gas, left1, right1, normal);
213       UnpackState_U(Flux0_state, Flux0);
214       UnpackState_U(Flux1_state, Flux1);
215       for (int j = 0; j < 5; j++) dFlux_fd[j] = (Flux1[j] - Flux0[j]) / eps;
216     }
217 
218     snprintf(buf, sizeof buf, "RiemannFlux_HLL i=%d: Flux", i);
219     PetscCall(CheckQWithTolerance(dFlux_fd, dFlux, dFlux_fd, buf, rtol_0, rtol_u, rtol_4));
220   }
221   PetscFunctionReturn(PETSC_SUCCESS);
222 }
223 
224 // @brief Verify RiemannFlux_HLLC_fwd function against finite-difference approximation
225 static PetscErrorCode TestRiemannHLLC_fwd(NewtonianIdealGasContext gas, CeedScalar rtol_0, CeedScalar rtol_u, CeedScalar rtol_4) {
226   CeedScalar       eps = 4e-7;  // Finite difference step
227   char             buf[128];
228   const CeedScalar T           = 200;
229   const CeedScalar rho         = 1.2;
230   const CeedScalar p           = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
231   const CeedScalar u_base      = 40;
232   const CeedScalar u[3]        = {u_base, u_base * 1.1, u_base * 1.2};
233   const CeedScalar Y0_left[5]  = {p, u[0], u[1], u[2], T};
234   const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T};
235   CeedScalar       normal[3]   = {1, 2, 3};
236 
237   PetscFunctionBeginUser;
238   State left0  = StateFromY(gas, Y0_left);
239   State right0 = StateFromY(gas, Y0_right);
240   ScaleN(normal, 1 / Norm3(normal), 3);
241 
242   for (int i = 0; i < 10; i++) {
243     CeedScalar dFlux[5] = {0.}, dFlux_fd[5] = {0.};
244     {  // Calculate dFlux using *_fwd function
245       CeedScalar dY_right[5] = {0};
246       CeedScalar dY_left[5]  = {0};
247 
248       if (i < 5) {
249         dY_left[i] = Y0_left[i];
250       } else {
251         dY_right[i % 5] = Y0_right[i % 5];
252       }
253       State dleft0  = StateFromY_fwd(gas, left0, dY_left);
254       State dright0 = StateFromY_fwd(gas, right0, dY_right);
255 
256       StateConservative dFlux_state = RiemannFlux_HLLC_fwd(gas, left0, dleft0, right0, dright0, normal);
257       UnpackState_U(dFlux_state, dFlux);
258     }
259 
260     {  // Calculate dFlux_fd via finite difference approximation
261       CeedScalar Y1_left[5]  = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]};
262       CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]};
263       CeedScalar Flux0[5], Flux1[5];
264 
265       if (i < 5) {
266         Y1_left[i] *= 1 + eps;
267       } else {
268         Y1_right[i % 5] *= 1 + eps;
269       }
270       State left1  = StateFromY(gas, Y1_left);
271       State right1 = StateFromY(gas, Y1_right);
272 
273       StateConservative Flux0_state = RiemannFlux_HLLC(gas, left0, right0, normal);
274       StateConservative Flux1_state = RiemannFlux_HLLC(gas, left1, right1, normal);
275       UnpackState_U(Flux0_state, Flux0);
276       UnpackState_U(Flux1_state, Flux1);
277       for (int j = 0; j < 5; j++) dFlux_fd[j] = (Flux1[j] - Flux0[j]) / eps;
278     }
279 
280     snprintf(buf, sizeof buf, "RiemannFlux_HLLC i=%d: Flux", i);
281     PetscCall(CheckQWithTolerance(dFlux_fd, dFlux, dFlux_fd, buf, rtol_0, rtol_u, rtol_4));
282   }
283   PetscFunctionReturn(PETSC_SUCCESS);
284 }
285 
286 // @brief Verify ComputeHLLSpeeds_Roe_fwd function against finite-difference approximation
287 static PetscErrorCode TestComputeHLLSpeeds_Roe_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) {
288   CeedScalar       eps = 4e-7;  // Finite difference step
289   char             buf[128];
290   const CeedScalar T           = 200;
291   const CeedScalar rho         = 1.2;
292   const CeedScalar p           = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
293   const CeedScalar u_base      = 40;
294   const CeedScalar u[3]        = {u_base, u_base * 1.1, u_base * 1.2};
295   const CeedScalar Y0_left[5]  = {p, u[0], u[1], u[2], T};
296   const CeedScalar Y0_right[5] = {1.2 * p, 1.2 * u[0], 1.2 * u[1], 1.2 * u[2], 1.2 * T};
297   CeedScalar       normal[3]   = {1, 2, 3};
298 
299   PetscFunctionBeginUser;
300   State left0  = StateFromY(gas, Y0_left);
301   State right0 = StateFromY(gas, Y0_right);
302   ScaleN(normal, 1 / Norm3(normal), 3);
303   CeedScalar u_left0  = Dot3(left0.Y.velocity, normal);
304   CeedScalar u_right0 = Dot3(right0.Y.velocity, normal);
305 
306   for (int i = 0; i < 10; i++) {
307     CeedScalar ds_left, ds_right, ds_left_fd, ds_right_fd;
308     {  // Calculate ds_{left,right} using *_fwd function
309       CeedScalar dY_right[5] = {0};
310       CeedScalar dY_left[5]  = {0};
311 
312       if (i < 5) {
313         dY_left[i] = Y0_left[i];
314       } else {
315         dY_right[i % 5] = Y0_right[i % 5];
316       }
317       State      dleft0   = StateFromY_fwd(gas, left0, dY_left);
318       State      dright0  = StateFromY_fwd(gas, right0, dY_right);
319       CeedScalar du_left  = Dot3(dleft0.Y.velocity, normal);
320       CeedScalar du_right = Dot3(dright0.Y.velocity, normal);
321 
322       CeedScalar s_left, s_right;  // Throw away
323       ComputeHLLSpeeds_Roe_fwd(gas, left0, dleft0, u_left0, du_left, right0, dright0, u_right0, du_right, &s_left, &ds_left, &s_right, &ds_right);
324     }
325 
326     {  // Calculate ds_{left,right}_fd via finite difference approximation
327       CeedScalar Y1_left[5]  = {Y0_left[0], Y0_left[1], Y0_left[2], Y0_left[3], Y0_left[4]};
328       CeedScalar Y1_right[5] = {Y0_right[0], Y0_right[1], Y0_right[2], Y0_right[3], Y0_right[4]};
329 
330       if (i < 5) {
331         Y1_left[i] *= 1 + eps;
332       } else {
333         Y1_right[i % 5] *= 1 + eps;
334       }
335       State      left1    = StateFromY(gas, Y1_left);
336       State      right1   = StateFromY(gas, Y1_right);
337       CeedScalar u_left1  = Dot3(left1.Y.velocity, normal);
338       CeedScalar u_right1 = Dot3(right1.Y.velocity, normal);
339 
340       CeedScalar s_left0, s_right0, s_left1, s_right1;
341       ComputeHLLSpeeds_Roe(gas, left0, u_left0, right0, u_right0, &s_left0, &s_right0);
342       ComputeHLLSpeeds_Roe(gas, left1, u_left1, right1, u_right1, &s_left1, &s_right1);
343       ds_left_fd  = (s_left1 - s_left0) / eps;
344       ds_right_fd = (s_right1 - s_right0) / eps;
345     }
346 
347     snprintf(buf, sizeof buf, "ComputeHLLSpeeds_Roe i=%d:", i);
348     {
349       CeedScalar divisor_threshold = 10 * CEED_EPSILON;
350       CeedScalar ds_left_err, ds_right_err;
351 
352       ds_left_err  = RelativeError(ds_left_fd, ds_left, ds_left_fd, divisor_threshold);
353       ds_right_err = RelativeError(ds_right_fd, ds_right, ds_right_fd, divisor_threshold);
354       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);
355       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);
356     }
357   }
358   PetscFunctionReturn(PETSC_SUCCESS);
359 }
360 
361 // @brief Verify TotalSpecificEnthalpy_fwd function against finite-difference approximation
362 static PetscErrorCode TestTotalSpecificEnthalpy_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) {
363   CeedScalar       eps = 4e-7;  // Finite difference step
364   char             buf[128];
365   const CeedScalar T      = 200;
366   const CeedScalar rho    = 1.2;
367   const CeedScalar p      = (HeatCapacityRatio(gas) - 1) * rho * gas->cv * T;
368   const CeedScalar u_base = 40;
369   const CeedScalar u[3]   = {u_base, u_base * 1.1, u_base * 1.2};
370   const CeedScalar Y0[5]  = {p, u[0], u[1], u[2], T};
371 
372   PetscFunctionBeginUser;
373   State state0 = StateFromY(gas, Y0);
374 
375   for (int i = 0; i < 5; i++) {
376     CeedScalar dH, dH_fd;
377     {  // Calculate dH using *_fwd function
378       CeedScalar dY[5] = {0};
379 
380       dY[i]         = Y0[i];
381       State dstate0 = StateFromY_fwd(gas, state0, dY);
382       dH            = TotalSpecificEnthalpy_fwd(gas, state0, dstate0);
383     }
384 
385     {  // Calculate dH_fd via finite difference approximation
386       CeedScalar H0, H1;
387       CeedScalar Y1[5] = {Y0[0], Y0[1], Y0[2], Y0[3], Y0[4]};
388       Y1[i] *= 1 + eps;
389       State state1 = StateFromY(gas, Y1);
390 
391       H0    = TotalSpecificEnthalpy(gas, state0);
392       H1    = TotalSpecificEnthalpy(gas, state1);
393       dH_fd = (H1 - H0) / eps;
394     }
395 
396     snprintf(buf, sizeof buf, "TotalSpecificEnthalpy i=%d:", i);
397     {
398       CeedScalar divisor_threshold = 10 * CEED_EPSILON;
399       CeedScalar dH_err;
400 
401       dH_err = RelativeError(dH_fd, dH, dH_fd, divisor_threshold);
402       if (fabs(dH_err) >= rtol) printf("%s dH error %g (expected %.10e, got %.10e)\n", buf, dH_err, dH_fd, dH);
403     }
404   }
405   PetscFunctionReturn(PETSC_SUCCESS);
406 }
407 
408 // @brief Verify RoeSetup_fwd function against finite-difference approximation
409 static PetscErrorCode TestRowSetup_fwd(NewtonianIdealGasContext gas, CeedScalar rtol) {
410   CeedScalar       eps = 4e-7;  // Finite difference step
411   char             buf[128];
412   const CeedScalar rho0[2] = {1.2, 1.4};
413 
414   PetscFunctionBeginUser;
415   for (int i = 0; i < 2; i++) {
416     RoeWeights dR, dR_fd;
417     {  // Calculate using *_fwd function
418       CeedScalar drho[5] = {0};
419 
420       drho[i] = rho0[i];
421       dR      = RoeSetup_fwd(rho0[0], rho0[1], drho[0], drho[1]);
422     }
423 
424     {  // Calculate via finite difference approximation
425       RoeWeights R0, R1;
426       CeedScalar rho1[5] = {rho0[0], rho0[1]};
427       rho1[i] *= 1 + eps;
428 
429       R0          = RoeSetup(rho0[0], rho0[1]);
430       R1          = RoeSetup(rho1[0], rho1[1]);
431       dR_fd.left  = (R1.left - R0.left) / eps;
432       dR_fd.right = (R1.right - R0.right) / eps;
433     }
434 
435     snprintf(buf, sizeof buf, "RoeSetup i=%d:", i);
436     {
437       CeedScalar divisor_threshold = 10 * CEED_EPSILON;
438       RoeWeights dR_err;
439 
440       dR_err.left  = RelativeError(dR_fd.left, dR.left, dR_fd.left, divisor_threshold);
441       dR_err.right = RelativeError(dR_fd.right, dR.right, dR_fd.right, divisor_threshold);
442       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);
443       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);
444     }
445   }
446   PetscFunctionReturn(PETSC_SUCCESS);
447 }
448 
449 // @brief Test Riemann solver related `*_fwd` functions via finite-difference approximation
450 static PetscErrorCode RiemannSolverUnitTests(NewtonianIdealGasContext gas, CeedScalar rtol) {
451   PetscFunctionBeginUser;
452   PetscCall(TestRiemannHLL_fwd(gas, rtol, rtol, rtol));
453   PetscCall(TestRiemannHLLC_fwd(gas, rtol, rtol, rtol));
454   PetscCall(TestComputeHLLSpeeds_Roe_fwd(gas, rtol));
455   PetscCall(TestTotalSpecificEnthalpy_fwd(gas, rtol));
456   PetscCall(TestRowSetup_fwd(gas, rtol));
457   PetscFunctionReturn(PETSC_SUCCESS);
458 }
459