24#ifndef DUMUX_FLUX_SHALLOW_WATER_EXACT_RIEMANN_HH
25#define DUMUX_FLUX_SHALLOW_WATER_EXACT_RIEMANN_HH
31namespace ShallowWater{
33template<
typename Scalar>
80 constexpr int maxsteps = 200;
81 constexpr Scalar tol = 1.0E-12;
89 const auto cl = sqrt(grav * max(dl, 0.0));
90 const auto cr = sqrt(grav * max(dr, 0.0));
91 const auto dcrit = (ur - ul) - 2.0*(cl + cr);
94 if ((dl <= 0) || (dr <= 0) || (dcrit >= 0))
99 const auto shr = ur + cr;
112 const auto str = ur - 2.0*cr;
115 sol.
velocityX = (ur - 2.0 *cr + 2.0*s)/3.0;
116 const auto c = (-ur + 2.0*cr +s )/3.0;
136 const auto shl = ul-cl;
147 const auto stl = ul + 2.0 *cl;
151 sol.
velocityX = (ul + 2.0 *cl + 2.0*s)/3.0;
152 const auto c = (ul + 2.0*cl - s)/3.0;
169 const auto shl = ul - cl;
170 const auto shr = ur + cr;
171 const auto ssl = ul + 2.0 * cl;
172 const auto ssr = ur - 2.0 * cr;
181 if((s > shl) && (s <= ssl))
183 sol.
velocityX = (ul + 2.0 * cl + 2.0*s)/3.0;
184 const auto c = (ul + 2.0*cl - s)/3.0;
189 if((s > ssl) && (s <= ssr))
196 if ((s > ssr ) && (s <= shr))
198 sol.
velocityX = (ur - 2.0*cr + 2.0*s)/3.0;
199 const auto c = (-ur + 2.0*cr + s)/3.0;
218 const auto tmp = 0.5 * (cl + cr) - 0.25 * (ul - ur);
219 auto ds = (1.0/grav) * tmp * tmp;
222 const auto dmin = min(dl, dr);
226 const auto gel = sqrt(0.5*grav * (ds+dl)/(ds*dl));
227 const auto ger = sqrt(0.5*grav * (ds+dr)/(ds*dr));
228 ds = (gel * dl + ger * dr - (ur-ul))/(gel + ger);
235 Scalar fl = 0.0, fr = 0.0, fld = 0.0, frd = 0.0;
237 for (
int i=0; i <= maxsteps; ++i)
244 const auto c = sqrt(grav * ds);
251 const auto ges = sqrt(0.5 * grav * (ds+dl)/(ds*dl));
252 fl = (ds - dl) * ges;
253 fld = ges - 0.25 * grav * (ds - dl)/(ges * ds * ds);
259 const auto c = sqrt(grav * ds);
266 const auto ges = sqrt(0.5 * grav * (ds+dr)/(ds*dr));
267 fr = (ds - dr) * ges;
268 frd = ges - 0.25 * grav * (ds - dr)/(ges * ds * ds);
271 ds -= (fl + fr + ur - ul)/(fld + frd);
272 const auto cha = abs(ds - d0)/(0.5*(ds + d0));
290 const auto c = sqrt(grav * ds);
296 const auto ges = sqrt(0.5*grav*(ds + dl)/(ds*dl));
297 fl = (ds - dl) * ges;
302 const auto c = sqrt(grav * ds);
308 const auto ges = sqrt(0.5*grav*(ds + dr)/(ds*dr));
309 fr = (ds - dr) * ges;
313 const auto us = 0.5*(ul + ur) + 0.5*(fr - fl);
314 const auto cs = sqrt(grav*ds);
324 const auto ql = sqrt((ds + dl)*ds/(2.0*dl*dl));
325 const auto sl = ul - cl*ql;
345 const auto shl = ul-cl;
356 const auto stl = us -cs;
361 sol.
velocityX = (ul + 2.0*cl + 2.0 * s)/3.0;
362 const auto c = (ul + 2.0* cl - s)/3.0;
384 const auto qr = sqrt((ds + dr)*ds /(2.0*dr*dr));
385 const auto sr = ur + cr*qr;
406 const auto shr = ur + cr;
417 const auto str = us + cs;
422 sol.
velocityX = (ur - 2.0*cr + 2.0*s)/3.0;
423 const auto c = (-ur + 2.0*cr + s)/3.0;
RiemannSolution< Scalar > exactRiemann(const Scalar dl, const Scalar dr, const Scalar ul, const Scalar ur, const Scalar vl, const Scalar vr, const Scalar grav, const Scalar s=0.0)
Exact Riemann solver for the shallow water equations.
Definition: exactriemann.hh:66
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: exactriemann.hh:34
Scalar waterDepth
Definition: exactriemann.hh:36
std::array< Scalar, 3 > flux
Definition: exactriemann.hh:35
Scalar velocityY
Definition: exactriemann.hh:38
Scalar velocityX
Definition: exactriemann.hh:37