PLUTO
hll.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief HLL Riemann solver for RMHD.
5 
6  Solve the Riemann problem for the RMHD equations using the
7  single-state HLL solver by Toro.
8 
9  On input, this function takes left and right primitive state vectors
10  \c state->vL and \c state->vR at zone edge \c i+1/2;
11  On output, return flux and pressure vectors at the same interface
12  \c i+1/2 (note that the \c i refers to \c i+1/2).
13 
14  Also during this step, compute maximum wave propagation speed (cmax)
15  for explicit time step computation.
16 
17  \b Reference:
18  - "Riemann Solver and Numerical Methods for Fluid Dynamics"
19  by E.F. Toro (Chapter 10)
20  - "On Godunov-Type Method near Low Densities"
21  by B. Einfeldt, C.D. Munz, P.L. Roe, JCP 92, 273-295 (1991)
22 
23  \authors A. Mignone (mignone@ph.unito.it)
24  \date April 7, 2014
25 */
26 /* ///////////////////////////////////////////////////////////////////// */
27 #include"pluto.h"
28 
29 /* ********************************************************************* */
30 void HLL_Solver (const State_1D *state, int beg, int end,
31  real *cmax, Grid *grid)
32 /*!
33  * Solve Riemann problem for the adiabatic/isothermal MHD equations
34  * using the HLL Riemann solver.
35  *
36  * \param[in,out] state pointer to State_1D structure
37  * \param[in] beg initial grid index
38  * \param[out] end final grid index
39  * \param[out] cmax 1D array of maximum characteristic speeds
40  * \param[in] grid pointer to array of Grid structures.
41  *
42  *********************************************************************** */
43 {
44  int nv, i;
45  real scrh, bmin, bmax;
46  double *uL, *uR, *SL, *SR;
47  static double **fL, **fR;
48  static double *pL, *pR, *a2L, *a2R, *hL, *hR;
49  static double **Uhll;
50  static double **VL, **VR, **UL, **UR;
51 
52  if (fL == NULL){
53  fL = ARRAY_2D(NMAX_POINT, NFLX, double);
54  fR = ARRAY_2D(NMAX_POINT, NFLX, double);
55  pL = ARRAY_1D(NMAX_POINT, double);
56  pR = ARRAY_1D(NMAX_POINT, double);
57  SL = ARRAY_1D(NMAX_POINT, double);
58  SR = ARRAY_1D(NMAX_POINT, double);
59 
60  a2L = ARRAY_1D(NMAX_POINT, double);
61  a2R = ARRAY_1D(NMAX_POINT, double);
62 
63  hL = ARRAY_1D(NMAX_POINT, double);
64  hR = ARRAY_1D(NMAX_POINT, double);
65 
66  Uhll = ARRAY_2D(NMAX_POINT, NVAR, double); /* NVAR and not NFLX since it needs */
67  /* to be converted in HLL_DIVB_SOUCE */
68  #ifdef GLM_MHD
69  VL = ARRAY_2D(NMAX_POINT, NVAR, double);
70  VR = ARRAY_2D(NMAX_POINT, NVAR, double);
71  UL = ARRAY_2D(NMAX_POINT, NVAR, double);
72  UR = ARRAY_2D(NMAX_POINT, NVAR, double);
73  #endif
74  }
75 
76  #ifdef GLM_MHD
77  GLM_Solve (state, VL, VR, beg, end, grid);
78  PrimToCons (VL, UL, beg, end);
79  PrimToCons (VR, UR, beg, end);
80  #else
81  VL = state->vL; UL = state->uL;
82  VR = state->vR; UR = state->uR;
83  #endif
84 
85 /* ----------------------------------------------------
86  compute sound speed & fluxes at zone interfaces
87  ---------------------------------------------------- */
88 
89  SoundSpeed2 (VL, a2L, hL, beg, end, FACE_CENTER, grid);
90  SoundSpeed2 (VR, a2R, hR, beg, end, FACE_CENTER, grid);
91 
92  Flux (UL, VL, hL, fL, pL, beg, end);
93  Flux (UR, VR, hR, fR, pR, beg, end);
94 
95 /* ----------------------------------------
96  get max and min signal velocities
97  ---------------------------------------- */
98 
99  SL = state->SL; SR = state->SR;
100  HLL_Speed (VL, VR, a2L, a2R, hL, hR, SL, SR, beg, end);
101 
102  for (i = beg; i <= end; i++) {
103 
104  scrh = MAX(fabs(SL[i]), fabs(SR[i]));
105  cmax[i] = scrh;
106 
107  uL = UL[i]; uR = UR[i];
108 
109 /* ---- Compute fluxes ---- */
110 
111  bmin = MIN(0.0, SL[i]);
112  bmax = MAX(0.0, SR[i]);
113  scrh = 1.0/(bmax - bmin);
114  for (nv = NFLX; nv--; ){
115  state->flux[i][nv] = bmin*bmax*(uR[nv] - uL[nv])
116  + bmax*fL[i][nv] - bmin*fR[i][nv];
117  state->flux[i][nv] *= scrh;
118  }
119  state->press[i] = (bmax*pL[i] - bmin*pR[i])*scrh;
120  }
121 
122 /* --------------------------------------------------------
123  initialize source term
124  -------------------------------------------------------- */
125 
126  #if DIVB_CONTROL == EIGHT_WAVES
127 /*
128  POWELL_DIVB_SOURCE (state, beg, end, grid);
129 */
130 
131  /* ----------------------------------------------------
132  to avoid conversion problems in HLL_DIVB_SOURCE,
133  we use the HLL average provided by SR = -SL = 1
134  ---------------------------------------------------- */
135 /*
136  for (i = beg; i <= end; i++) {
137  uL = state->uL[i]; uR = state->uR[i];
138  for (nv = 0; nv < NFLX; nv++) {
139  Uhll[i][nv] = 0.5*(uR[nv] + uL[nv] + fL[i][nv] - fR[i][nv]);
140  }
141  Uhll[i][MXn] += (pL[i] - pR[i])*0.5;
142  for (nv = NFLX; nv < NVAR; nv++) Uhll[i][nv] = 0.0;
143  }
144 */
145  HLL_DIVB_SOURCE (state, Uhll, beg + 1, end, grid);
146  #endif
147 
148 }
#define MAX(a, b)
Definition: macros.h:101
void GLM_Solve(const State_1D *state, double **VL, double **VR, int beg, int end, Grid *grid)
Definition: glm.c:24
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 * SR
Rightmost velocity in the Riemann fan at i+1/2.
Definition: structs.h:167
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
void HLL_DIVB_SOURCE(const State_1D *, double **, int, int, Grid *)
Definition: source.c:53
#define MIN(a, b)
Definition: macros.h:104
#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 * SL
Leftmost velocity in the Riemann fan at i+1/2.
Definition: structs.h:166
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
void PrimToCons(double **uprim, double **ucons, int ibeg, int iend)
Definition: mappers.c:26
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
#define NVAR
Definition: pluto.h:609
double ** uL
same as vL, in conservative vars
Definition: structs.h:144