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