PLUTO
hll.c
Go to the documentation of this file.
1 #include"pluto.h"
2 
3 /* ********************************************************************* */
4 void HLL_Solver (const State_1D *state, int beg, int end,
5  real *cmax, Grid *grid)
6 /*
7  *
8  * NAME
9  *
10  * HLL_SOLVER
11  *
12  *
13  * PURPOSE
14  *
15  * Solve riemann problem for the Euler equations
16  * using th HLLE solver;
17  *
18  * Reference: --
19  *
20  *
21  * LAST_MODIFIED
22  *
23  * April 4th 2006, by Andrea Mignone (mignone@to.astro.it)
24  *
25  *
26  **************************************************************************** */
27 {
28  int nv, i;
29  real scrh;
30  double *uL, *uR;
31  static real **fL, **fR;
32  static real *SL, *SR, *pL, *pR;
33  static double *a2L, *a2R, *hL, *hR;
34 
35  if (fL == NULL){
36  fL = ARRAY_2D(NMAX_POINT, NFLX, double);
37  fR = ARRAY_2D(NMAX_POINT, NFLX, double);
38 
39  pR = ARRAY_1D(NMAX_POINT, double);
40  pL = ARRAY_1D(NMAX_POINT, double);
41  SR = ARRAY_1D(NMAX_POINT, double);
42  SL = ARRAY_1D(NMAX_POINT, double);
43 
44  hR = ARRAY_1D(NMAX_POINT, double);
45  hL = ARRAY_1D(NMAX_POINT, double);
46 
47  a2R = ARRAY_1D(NMAX_POINT, double);
48  a2L = ARRAY_1D(NMAX_POINT, double);
49  }
50 /* ----------------------------------------------------
51  compute sound speed & fluxes at zone interfaces
52  ---------------------------------------------------- */
53 
54  SoundSpeed2 (state->vL, a2L, hL, beg, end, FACE_CENTER, grid);
55  SoundSpeed2 (state->vR, a2R, hR, beg, end, FACE_CENTER, grid);
56 
57  Flux (state->uL, state->vL, a2L, fL, pL, beg, end);
58  Flux (state->uR, state->vR, a2R, fR, pR, beg, end);
59 
60  HLL_Speed (state->vL, state->vR, a2L, a2R, SL, SR, beg, end);
61 
62  for (i = beg; i <= end; i++) {
63 
64  scrh = MAX(fabs(SL[i]), fabs(SR[i]));
65  cmax[i] = scrh;
66 
67 /* --------------------------------------------------
68  compute HLL flux
69  -------------------------------------------------- */
70 
71  if (SL[i] >= 0.0){
72 
73  for (nv = NFLX; nv--; ) state->flux[i][nv] = fL[i][nv];
74  state->press[i] = pL[i];
75 
76  }else if (SR[i] <= 0.0){
77 
78  for (nv = NFLX; nv--; ) state->flux[i][nv] = fR[i][nv];
79  state->press[i] = pR[i];
80 
81  }else{
82 
83  uL = state->uL[i];
84  uR = state->uR[i];
85 
86  scrh = 1.0/(SR[i] - SL[i]);
87  for (nv = NFLX; nv--; ){
88  state->flux[i][nv] = SL[i]*SR[i]*(uR[nv] - uL[nv])
89  + SR[i]*fL[i][nv] - SL[i]*fR[i][nv];
90  state->flux[i][nv] *= scrh;
91  }
92  state->press[i] = (SR[i]*pL[i] - SL[i]*pR[i])*scrh;
93  }
94  }
95 }
#define MAX(a, b)
Definition: macros.h:101
double ** flux
upwind flux computed with the Riemann solver
Definition: structs.h:149
double real
Definition: pluto.h:488
void Flux(double **u, double **w, double *a2, double **fx, double *p, int beg, int end)
Definition: fluxes.c:23
tuple scrh
Definition: configure.py:200
double ** vR
Primitive variables to the right of the interface, .
Definition: structs.h:139
void SoundSpeed2(double **v, double *cs2, double *h, int beg, int end, int pos, Grid *grid)
Definition: eos.c:16
#define NFLX
Definition: mod_defs.h:32
void HLL_Solver(const State_1D *state, int beg, int end, double *cmax, Grid *grid)
Definition: hll.c:30
#define FACE_CENTER
Definition: pluto.h:206
Definition: structs.h:78
double ** uR
same as vR, in conservative vars
Definition: structs.h:145
PLUTO main header file.
#define ARRAY_1D(nx, type)
Definition: prototypes.h:170
long int NMAX_POINT
Maximum number of points among the three directions, boundaries excluded.
Definition: globals.h:62
int i
Definition: analysis.c:2
void HLL_Speed(double **vL, double **vR, double *a2L, double *a2R, double *SL, double *SR, int beg, int end)
Definition: hll_speed.c:24
double ** vL
Primitive variables to the left of the interface, .
Definition: structs.h:136
double * press
Upwind pressure term computed with the Riemann solver.
Definition: structs.h:164
#define ARRAY_2D(nx, ny, type)
Definition: prototypes.h:171
double ** uL
same as vL, in conservative vars
Definition: structs.h:144