Shamrock 2025.10.0
Astrophysical Code
Loading...
Searching...
No Matches
iterative.hpp
Go to the documentation of this file.
1// -------------------------------------------------------//
2//
3// SHAMROCK code for hydrodynamics
4// Copyright (c) 2021-2026 Timothée David--Cléris <tim.shamrock@proton.me>
5// SPDX-License-Identifier: CeCILL Free Software License Agreement v2.1
6// Shamrock is licensed under the CeCILL 2.1 License, see LICENSE for more information
7//
8// -------------------------------------------------------//
9
10#pragma once
11
28#include "shambackends/math.hpp"
29#include "shambackends/sycl.hpp"
30
31namespace shammodels::gsph::riemann {
32
39 template<class Tscal>
41 Tscal p_star;
42 Tscal v_star;
43 };
44
68 template<class Tscal>
70 Tscal u_L,
71 Tscal rho_L,
72 Tscal p_L,
73 Tscal u_R,
74 Tscal rho_R,
75 Tscal p_R,
76 Tscal gamma,
77 Tscal tol = Tscal{1.0e-6},
78 u32 max_iter = 20) {
79
80 RiemannResult<Tscal> result;
81
82 // Safety check for non-physical values
83 const Tscal smallp = Tscal{1.0e-25};
84 const Tscal smallrho = Tscal{1.0e-25};
85
86 if (rho_L < smallrho || rho_R < smallrho || p_L < smallp || p_R < smallp) {
87 // Return acoustic approximation for near-vacuum
88 result.p_star = sycl::fmax(smallp, Tscal{0.5} * (p_L + p_R));
89 result.v_star = Tscal{0.5} * (u_L + u_R);
90 return result;
91 }
92
93 // Derived constants
94 const Tscal gm1 = gamma - Tscal{1};
95 const Tscal gp1 = gamma + Tscal{1};
96 const Tscal gamma1 = Tscal{0.5} * gp1 / gamma; // (gamma+1)/(2*gamma)
97
98 // Specific volumes
99 const Tscal V_L = Tscal{1} / rho_L;
100 const Tscal V_R = Tscal{1} / rho_R;
101
102 // Lagrangian sound speeds: c_lag = sqrt(gamma * p * rho)
103 const Tscal c_L = sycl::sqrt(gamma * p_L * rho_L);
104 const Tscal c_R = sycl::sqrt(gamma * p_R * rho_R);
105
106 // Initial guess for p_star using PVRS (Primitive Variable Riemann Solver)
107 // p_star = p_L + (p_R - p_L - c_R*(u_R - u_L)) * c_L / (c_L + c_R)
108 Tscal p_star = p_R - p_L - c_R * (u_R - u_L);
109 p_star = p_L + p_star * c_L / (c_L + c_R);
110 p_star = sycl::fmax(p_star, smallp);
111
112 // Newton-Raphson iteration
113 for (u32 iter = 0; iter < max_iter; ++iter) {
114 const Tscal p_star_old = p_star;
115
116 // Left wave impedance: W_L = c_L * sqrt(1 + gamma1*(p_star - p_L)/p_L)
117 Tscal W_L = Tscal{1} + gamma1 * (p_star - p_L) / p_L;
118 W_L = c_L * sycl::sqrt(sycl::fmax(W_L, smallp));
119
120 // Right wave impedance: W_R = c_R * sqrt(1 + gamma1*(p_star - p_R)/p_R)
121 Tscal W_R = Tscal{1} + gamma1 * (p_star - p_R) / p_R;
122 W_R = c_R * sycl::sqrt(sycl::fmax(W_R, smallp));
123
124 // Derivatives dW/dp for Newton-Raphson
125 // Z_L = -dW_L/dp * W_L (note the sign convention)
126 // Add smallp to denominator to prevent division by zero near vacuum/shock
127 Tscal Z_L = Tscal{4} * V_L * W_L * W_L;
128 Z_L = -Z_L * W_L / (Z_L - gp1 * (p_star - p_L) + smallp);
129
130 // Z_R = dW_R/dp * W_R
131 Tscal Z_R = Tscal{4} * V_R * W_R * W_R;
132 Z_R = Z_R * W_R / (Z_R - gp1 * (p_star - p_R) + smallp);
133
134 // Intermediate velocities from each side
135 // u*_L = u_L - (p* - p_L) / W_L
136 // u*_R = u_R + (p* - p_R) / W_R
137 const Tscal ustar_L = u_L - (p_star - p_L) / W_L;
138 const Tscal ustar_R = u_R + (p_star - p_R) / W_R;
139
140 // Newton-Raphson update: p_new = p - f(p)/f'(p)
141 // where f(p) = u*_R - u*_L (velocity mismatch)
142 // and f'(p) = du*_R/dp - du*_L/dp = 1/Z_R - 1/Z_L
143 const Tscal denom = Z_R - Z_L;
144 if (sycl::fabs(denom) > smallp) {
145 p_star = p_star + (ustar_R - ustar_L) * (Z_L * Z_R) / denom;
146 }
147 p_star = sycl::fmax(smallp, p_star);
148
149 // Check convergence
150 if (sycl::fabs(p_star - p_star_old) / p_star < tol) {
151 break;
152 }
153 }
154
155 // Recalculate wave impedances with final p_star
156 Tscal W_L = Tscal{1} + gamma1 * (p_star - p_L) / p_L;
157 W_L = c_L * sycl::sqrt(sycl::fmax(W_L, smallp));
158
159 Tscal W_R = Tscal{1} + gamma1 * (p_star - p_R) / p_R;
160 W_R = c_R * sycl::sqrt(sycl::fmax(W_R, smallp));
161
162 // Calculate final u_star (average of left and right estimates)
163 const Tscal ustar_L = u_L - (p_star - p_L) / W_L;
164 const Tscal ustar_R = u_R + (p_star - p_R) / W_R;
165 const Tscal u_star = Tscal{0.5} * (ustar_L + ustar_R);
166
167 result.p_star = p_star;
168 result.v_star = u_star;
169
170 return result;
171 }
172
189 template<class Tscal>
191 Tscal u_L, Tscal rho_L, Tscal p_L, Tscal u_R, Tscal rho_R, Tscal p_R, Tscal gamma) {
192
194 const Tscal smallval = Tscal{1.0e-25};
195
196 // Compute Eulerian sound speeds
197 const Tscal c_L = sycl::sqrt(gamma * p_L / sycl::fmax(rho_L, smallval));
198 const Tscal c_R = sycl::sqrt(gamma * p_R / sycl::fmax(rho_R, smallval));
199
200 // Roe averages for wave speed estimates
201 const Tscal sqrt_rho_L = sycl::sqrt(rho_L);
202 const Tscal sqrt_rho_R = sycl::sqrt(rho_R);
203 const Tscal roe_inv = Tscal{1} / (sqrt_rho_L + sqrt_rho_R + smallval);
204
205 const Tscal u_roe = (sqrt_rho_L * u_L + sqrt_rho_R * u_R) * roe_inv;
206 const Tscal c_roe = (sqrt_rho_L * c_L + sqrt_rho_R * c_R) * roe_inv;
207
208 // Wave speed estimates (following reference implementation)
209 const Tscal S_L = sycl::fmin(u_L - c_L, u_roe - c_roe);
210 const Tscal S_R = sycl::fmax(u_R + c_R, u_roe + c_roe);
211
212 // HLL flux formula (following reference g_fluid_force.cpp hll_solver)
213 // c1 = rho_L * (S_L - u_L)
214 // c2 = rho_R * (S_R - u_R)
215 // c3 = 1 / (c1 - c2)
216 // c4 = p_L - u_L * c1
217 // c5 = p_R - u_R * c2
218 // v* = (c5 - c4) * c3
219 // p* = (c1 * c5 - c2 * c4) * c3
220 const Tscal c1 = rho_L * (S_L - u_L);
221 const Tscal c2 = rho_R * (S_R - u_R);
222 const Tscal c3 = Tscal{1} / (c1 - c2 + smallval);
223 const Tscal c4 = p_L - u_L * c1;
224 const Tscal c5 = p_R - u_R * c2;
225
226 result.v_star = (c5 - c4) * c3;
227 result.p_star = sycl::fmax(smallval, (c1 * c5 - c2 * c4) * c3);
228
229 return result;
230 }
231
232} // namespace shammodels::gsph::riemann
std::uint32_t u32
32 bit unsigned integer
RiemannResult< Tscal > hllc_solver(Tscal u_L, Tscal rho_L, Tscal p_L, Tscal u_R, Tscal rho_R, Tscal p_R, Tscal gamma)
HLL approximate Riemann solver.
RiemannResult< Tscal > iterative_solver(Tscal u_L, Tscal rho_L, Tscal p_L, Tscal u_R, Tscal rho_R, Tscal p_R, Tscal gamma, Tscal tol=Tscal{1.0e-6}, u32 max_iter=20)
Iterative Riemann solver (van Leer 1997)
Definition iterative.hpp:69
Tscal v_star
Interface velocity (normal component)
Definition iterative.hpp:42