Lines Matching +full:- +full:f

1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
6 // The left and right states are specified from the perspective of an outward-facing normal vector …
10 // |------------| /
12 // | Left |----> Right
14 // |------------|
17 // Much of the work references Eleuterio F. Toro's "Riemann Solvers and Numerical Methods for Fluid…
37 …CeedScalar r_right = (sqrt_left / (2 * sqrt_right * square_sum_root)) * drho_right - (sqrt_right /… in RoeSetup_fwd()
38 …CeedScalar r_left = (sqrt_right / (2 * sqrt_left * square_sum_root)) * drho_left - (sqrt_left / (… in RoeSetup_fwd()
58 …F_hll[i] = (s_right * F_left[i] - s_left * F_right[i] + s_left * s_right * (U_right[i] - U_left[i]… in Flux_HLL()
60 StateConservative F = { in Flux_HLL() local
65 return F; in Flux_HLL()
83 const CeedScalar S_diff = S_r - S_l; in Flux_HLL_fwd()
85 dF_hll[i] += (S_l * (-F_l[i] + F_r[i] + S_l * U_l[i] - S_l * U_r[i]) / Square(S_diff)) * dS_r; in Flux_HLL_fwd()
86 dF_hll[i] += (S_r * (F_l[i] - F_r[i] - S_r * U_l[i] + S_r * U_r[i]) / Square(S_diff)) * dS_l; in Flux_HLL_fwd()
87 dF_hll[i] += (S_r * dF_l[i] - S_l * dF_r[i] + S_r * S_l * (dU_r[i] - dU_l[i])) / S_diff; in Flux_HLL_fwd()
103 // Roe average eigenvalues for left and right non-linear waves. in ComputeHLLSpeeds_Roe()
110 CeedScalar a_roe = sqrt((gamma - 1) * (H_roe - 0.5 * Square(u_roe))); in ComputeHLLSpeeds_Roe()
113 *s_left = u_roe - a_roe; in ComputeHLLSpeeds_Roe()
125 // Roe average eigenvalues for left and right non-linear waves. in ComputeHLLSpeeds_Roe_fwd()
137 CeedScalar a_roe = sqrt((gamma - 1) * (H_roe - 0.5 * Square(u_roe))); in ComputeHLLSpeeds_Roe_fwd()
138 …CeedScalar da_roe = 0.5 * sqrt((gamma - 1) / (H_roe - 0.5 * Square(u_roe))) * dH_roe; // (da/dH) … in ComputeHLLSpeeds_Roe_fwd()
139 …da_roe -= 0.5 * sqrt(gamma - 1) * u_roe / sqrt(H_roe - 0.5 * Square(u_roe)) * du_roe; // (da/du) … in ComputeHLLSpeeds_Roe_fwd()
141 *s_left = u_roe - a_roe; in ComputeHLLSpeeds_Roe_fwd()
142 *ds_left = du_roe - da_roe; in ComputeHLLSpeeds_Roe_fwd()
150 // The left and right states are specified from the perspective of an outward-facing normal vector …
180 // @brief Forward-mode Derivative of Harten Lax VanLeer (HLL) approximate Riemann solver.
217 CeedScalar fact = side.U.density * (s_side - u_side) / (s_side - s_star); in RiemannFlux_HLLC_Star()
218 CeedScalar denom = side.U.density * (s_side - u_side); in RiemannFlux_HLLC_Star()
223 side.Y.velocity[0] + (s_star - u_side) * normal[0], in RiemannFlux_HLLC_Star()
224 side.Y.velocity[1] + (s_star - u_side) * normal[1], in RiemannFlux_HLLC_Star()
225 side.Y.velocity[2] + (s_star - u_side) * normal[2], in RiemannFlux_HLLC_Star()
228 + (s_star - u_side) * (s_star + side.Y.pressure / denom) in RiemannFlux_HLLC_Star()
230 return StateConservativeAXPBYPCZ(1, F_side, s_side * fact, star, -s_side, side.U); in RiemannFlux_HLLC_Star()
237 CeedScalar fact = side.U.density * (s_side - u_side) / (s_side - s_star); in RiemannFlux_HLLC_Star_fwd()
238 …CeedScalar dfact = (side.U.density * (ds_side - du_side) + dside.U.density * (s_side - u_side)) / … in RiemannFlux_HLLC_Star_fwd()
239 - fact / (s_side - s_star) * (ds_side - ds_star); in RiemannFlux_HLLC_Star_fwd()
240 CeedScalar denom = side.U.density * (s_side - u_side); in RiemannFlux_HLLC_Star_fwd()
241 CeedScalar ddenom = side.U.density * (ds_side - du_side) + dside.U.density * (s_side - u_side); in RiemannFlux_HLLC_Star_fwd()
246 side.Y.velocity[0] + (s_star - u_side) * normal[0], in RiemannFlux_HLLC_Star_fwd()
247 side.Y.velocity[1] + (s_star - u_side) * normal[1], in RiemannFlux_HLLC_Star_fwd()
248 side.Y.velocity[2] + (s_star - u_side) * normal[2], in RiemannFlux_HLLC_Star_fwd()
251 + (s_star - u_side) * (s_star + side.Y.pressure / denom) in RiemannFlux_HLLC_Star_fwd()
256 dside.Y.velocity[0] + (ds_star - du_side) * normal[0], in RiemannFlux_HLLC_Star_fwd()
257 dside.Y.velocity[1] + (ds_star - du_side) * normal[1], in RiemannFlux_HLLC_Star_fwd()
258 dside.Y.velocity[2] + (ds_star - du_side) * normal[2], in RiemannFlux_HLLC_Star_fwd()
260 … dside.U.E_total / side.U.density - side.U.E_total / Square(side.U.density) * dside.U.density // in RiemannFlux_HLLC_Star_fwd()
261 + (ds_star - du_side) * (s_star + side.Y.pressure / denom) // in RiemannFlux_HLLC_Star_fwd()
262 …+ (s_star - u_side) * (ds_star + dside.Y.pressure / denom - side.Y.pressure / Square(denom) * dden… in RiemannFlux_HLLC_Star_fwd()
265 …const CeedScalar a[] = {1, ds_side * fact + s_side * dfact, s_side * fact, -ds_side, -s_sid… in RiemannFlux_HLLC_Star_fwd()
282 …CeedScalar numer = right.Y.pressure - left.Y.pressure + rhou_left * (s_left - u_left) - rhou_righ… in RiemannFlux_HLLC()
283 CeedScalar denom = left.U.density * (s_left - u_left) - right.U.density * (s_right - u_right); in RiemannFlux_HLLC()
316 CeedScalar numer = right.Y.pressure - left.Y.pressure // in RiemannFlux_HLLC_fwd()
317 + rhou_left * (s_left - u_left) // in RiemannFlux_HLLC_fwd()
318 - rhou_right * (s_right - u_right); in RiemannFlux_HLLC_fwd()
319 CeedScalar dnumer = dright.Y.pressure - dleft.Y.pressure // in RiemannFlux_HLLC_fwd()
320 + rhou_left * (ds_left - du_left) + drhou_left * (s_left - u_left) // in RiemannFlux_HLLC_fwd()
321 - rhou_right * (ds_right - du_right) - drhou_right * (s_right - u_right); in RiemannFlux_HLLC_fwd()
322 CeedScalar denom = left.U.density * (s_left - u_left) - right.U.density * (s_right - u_right); in RiemannFlux_HLLC_fwd()
323 CeedScalar ddenom = left.U.density * (ds_left - du_left) + dleft.U.density * (s_left - u_left) // in RiemannFlux_HLLC_fwd()
324- right.U.density * (ds_right - du_right) - dright.U.density * (s_right - u_right); in RiemannFlux_HLLC_fwd()
326 CeedScalar ds_star = dnumer / denom - numer / Square(denom) * ddenom; in RiemannFlux_HLLC_fwd()