Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
riemann.h
Go to the documentation of this file.
1/** @file riemann.h
2 */
3/* A collection of Riemann solvers for the Saint-Venant system
4 *
5 * References:
6 * [1] Kurganov, A., & Levy, D. (2002). Central-upwind schemes for the
7 * Saint-Venant system. Mathematical Modelling and Numerical
8 * Analysis, 36(3), 397-425.
9 */
10
11#define SQRT3 1.73205080756888
12#define epsilon 1e-30
13
14void kinetic (double hm, double hp, double um, double up, double Delta,
15 double * fh, double * fq, double * dtmax)
16{
17 double ci = sqrt(G*hm/2.);
18 double Mp = max(um + ci*SQRT3, 0.);
19 double Mm = max(um - ci*SQRT3, 0.);
20 double cig = ci/(6.*G*SQRT3);
21 *fh = cig*3.*(Mp*Mp - Mm*Mm);
22 *fq = cig*2.*(Mp*Mp*Mp - Mm*Mm*Mm);
23 if (Mp > 0.) {
24 double dt = CFL*Delta/Mp;
25 if (dt < *dtmax)
26 *dtmax = dt;
27 }
28
29 ci = sqrt(G*hp/2.);
30 Mp = min(up + ci*SQRT3, 0.);
31 Mm = min(up - ci*SQRT3, 0.);
32 cig = ci/(6.*G*SQRT3);
33 *fh += cig*3.*(Mp*Mp - Mm*Mm);
34 *fq += cig*2.*(Mp*Mp*Mp - Mm*Mm*Mm);
35 if (Mm < - epsilon) {
36 double dt = CFL*Delta/-Mm;
37 if (dt < *dtmax)
38 *dtmax = dt;
39 }
40}
41
42void kurganov (double hm, double hp, double um, double up, double Delta,
43 double * fh, double * fq, double * dtmax)
44{
45 double cp = sqrt(G*hp), cm = sqrt(G*hm);
46 double ap = max(up + cp, um + cm); ap = max(ap, 0.);
47 double am = min(up - cp, um - cm); am = min(am, 0.);
48 double qm = hm*um, qp = hp*up;
49 double a = max(ap, -am);
50 if (a > epsilon) {
51 *fh = (ap*qm - am*qp + ap*am*(hp - hm))/(ap - am); // (4.5) of [1]
52 *fq = (ap*(qm*um + G*sq(hm)/2.) - am*(qp*up + G*sq(hp)/2.) +
53 ap*am*(qp - qm))/(ap - am);
54 double dt = CFL*Delta/a;
55 if (dt < *dtmax)
56 *dtmax = dt;
57 }
58 else
59 *fh = 0., *fq = 0.;
60}
61
62void hllc (double hm, double hp, double um, double up, double Delta,
63 double * fh, double * fq, double * dtmax)
64{
65 double cm = sqrt (G*hm), cp = sqrt (G*hp);
66 double ustar = (um + up)/2. + cm - cp;
67 double cstar = (cm + cp)/2. + (um - up)/4.;
68 double SL = hm == 0. ? up - 2.*cp : min (um - cm, ustar - cstar);
69 double SR = hp == 0. ? um + 2.*cm : max (up + cp, ustar + cstar);
70
71 if (0. <= SL) {
72 *fh = um*hm;
73 *fq = hm*(um*um + G*hm/2.);
74 }
75 else if (0. >= SR) {
76 *fh = up*hp;
77 *fq = hp*(up*up + G*hp/2.);
78 }
79 else {
80 double fhm = um*hm;
81 double fum = hm*(um*um + G*hm/2.);
82 double fhp = up*hp;
83 double fup = hp*(up*up + G*hp/2.);
84 *fh = (SR*fhm - SL*fhp + SL*SR*(hp - hm))/(SR - SL);
85 *fq = (SR*fum - SL*fup + SL*SR*(hp*up - hm*um))/(SR - SL);
86 }
87
88 double a = max(fabs(SL), fabs(SR));
89 if (a > epsilon) {
90 double dt = CFL*Delta/a;
91 if (dt < *dtmax)
92 *dtmax = dt;
93 }
94}
double dtmax
The timestep is computed using the CFL condition on the face velocity field.
Definition all-mach.h:105
const vector a
Definition all-mach.h:59
double G
Definition atmosphere.h:12
int min
Definition balance.h:192
int x
Definition common.h:76
const scalar cm[]
Definition common.h:397
static number sq(number x)
Definition common.h:11
double dt
size_t max
Definition mtrace.h:8
#define epsilon
Definition riemann.h:12
void kinetic(double hm, double hp, double um, double up, double Delta, double *fh, double *fq, double *dtmax)
Definition riemann.h:14
void hllc(double hm, double hp, double um, double up, double Delta, double *fh, double *fq, double *dtmax)
Definition riemann.h:62
#define SQRT3
Definition riemann.h:11
void kurganov(double hm, double hp, double um, double up, double Delta, double *fh, double *fq, double *dtmax)
Definition riemann.h:42
double CFL
Definition utils.h:10