PLUTO
hll.c File Reference

HLL Riemann solver for MHD. More...

#include "pluto.h"
Include dependency graph for hll.c:

Go to the source code of this file.

Functions

void HLL_Solver (const State_1D *state, int beg, int end, real *cmax, Grid *grid)
 

Detailed Description

HLL Riemann solver for MHD.

Solve the Riemann problem for the MHD equations using the single-state HLL solver by Toro.

On input, this function takes left and right primitive state vectors state->vL and state->vR at zone edge i+1/2; On output, return flux and pressure vectors at the same interface i+1/2 (note that the i refers to i+1/2).

Also during this step, compute maximum wave propagation speed (cmax) for explicit time step computation.

Reference:

  • "Riemann Solver and Numerical Methods for Fluid Dynamics" by E.F. Toro (Chapter 10)
  • "On Godunov-Type Method near Low Densities" by B. Einfeldt, C.D. Munz, P.L. Roe, JCP 92, 273-295 (1991)
Authors
A. Mignone (migno.nosp@m.ne@p.nosp@m.h.uni.nosp@m.to.i.nosp@m.t)
Date
March 23, 2012

Definition in file hll.c.

Function Documentation

void HLL_Solver ( const State_1D state,
int  beg,
int  end,
real cmax,
Grid grid 
)

Solve Riemann problem for the adiabatic/isothermal MHD equations using the HLL Riemann solver.

Parameters
[in,out]statepointer to State_1D structure
[in]beginitial grid index
[out]endfinal grid index
[out]cmax1D array of maximum characteristic speeds
[in]gridpointer to array of Grid structures.

Definition at line 30 of file hll.c.

43 {
44  int nv, i, xdface;
45  double scrh;
46  double *vL, *vR, *uL, *uR, *SR, *SL;
47  static double **fL, **fR, **Uhll;
48  static double **VL, **VR, **UL, **UR;
49  static double *pL, *pR, *a2L, *a2R;
50  double **bgf;
51 
52  if (fL == NULL){
53  fL = ARRAY_2D(NMAX_POINT, NFLX, double);
54  fR = ARRAY_2D(NMAX_POINT, NFLX, double);
55 
56  pL = ARRAY_1D(NMAX_POINT, double);
57  pR = ARRAY_1D(NMAX_POINT, double);
58  a2L = ARRAY_1D(NMAX_POINT, double);
59  a2R = ARRAY_1D(NMAX_POINT, double);
60 
61  Uhll = ARRAY_2D(NMAX_POINT, NFLX, double);
62  #ifdef GLM_MHD
63  VL = ARRAY_2D(NMAX_POINT, NVAR, double);
64  VR = ARRAY_2D(NMAX_POINT, NVAR, double);
65  UL = ARRAY_2D(NMAX_POINT, NVAR, double);
66  UR = ARRAY_2D(NMAX_POINT, NVAR, double);
67  #endif
68  }
69 
70  #if BACKGROUND_FIELD == YES
71  bgf = GetBackgroundField (beg, end, FACE_CENTER, grid);
72  #endif
73 
74 /* ------------------------------------------------
75  Solve 2x2 Riemann problem with GLM Cleaning
76  ------------------------------------------------ */
77 
78  #ifdef GLM_MHD
79  GLM_Solve (state, VL, VR, beg, end, grid);
80  PrimToCons (VL, UL, beg, end);
81  PrimToCons (VR, UR, beg, end);
82  #else
83  VL = state->vL; UL = state->uL;
84  VR = state->vR; UR = state->uR;
85  #endif
86 
87 /* ----------------------------------------------------
88  compute sound speed & fluxes at zone interfaces
89  ---------------------------------------------------- */
90 
91  SoundSpeed2 (VL, a2L, NULL, beg, end, FACE_CENTER, grid);
92  SoundSpeed2 (VR, a2R, NULL, beg, end, FACE_CENTER, grid);
93 
94  Flux (UL, VL, a2L, bgf, fL, pL, beg, end);
95  Flux (UR, VR, a2R, bgf, fR, pR, beg, end);
96 
97 /* ----------------------------------------
98  get max and min signal velocities
99  ---------------------------------------- */
100 
101  SL = state->SL; SR = state->SR;
102  HLL_Speed (VL, VR, a2L, a2R, bgf, SL, SR, beg, end);
103 
104 /* ----------------------------------------
105  compute HLL flux
106  ---------------------------------------- */
107 
108  for (i = beg; i <= end; i++) {
109 
110  scrh = MAX(fabs(SL[i]), fabs(SR[i]));
111  cmax[i] = scrh;
112 
113  if (SL[i] > 0.0){
114 
115  for (nv = 0; nv < NFLX; nv++) {
116  state->flux[i][nv] = fL[i][nv];
117  }
118  state->press[i] = pL[i];
119 
120  }else if (SR[i] < 0.0){
121 
122  for (nv = 0; nv < NFLX; nv++) {
123  state->flux[i][nv] = fR[i][nv];
124  }
125  state->press[i] = pR[i];
126 
127  }else{
128 
129  uL = UL[i]; uR = UR[i];
130 
131  scrh = 1.0 / (SR[i] - SL[i]);
132 
133  for (nv = 0; nv < NFLX; nv++) {
134  state->flux[i][nv] = SL[i]*SR[i]*(uR[nv] - uL[nv]) +
135  SR[i]*fL[i][nv] - SL[i]*fR[i][nv];
136  state->flux[i][nv] *= scrh;
137  }
138  state->press[i] = (SR[i]*pL[i] - SL[i]*pR[i])*scrh;
139  }
140 
141  }
142 
143 
144 /* -----------------------------------------------------
145  initialize source term
146  ----------------------------------------------------- */
147 
148  #if DIVB_CONTROL == EIGHT_WAVES
149 /*
150  ROE_DIVB_SOURCE (state, beg, end, grid);
151 */
152 /*
153  for (i = beg; i <= end; i++) {
154  uR = state->uR[i]; uL = state->uL[i];
155  scrh = 1.0 / (SR[i] - SL[i]);
156  for (nv = 0; nv < NFLX; nv++) {
157  Uhll[i][nv] = SR[i]*uR[nv] - SL[i]*uL[nv] +
158  fL[i][nv] - fR[i][nv];
159  Uhll[i][nv] *= scrh;
160  }
161  Uhll[i][MXn] += (pL[i] - pR[i])*scrh;
162  }
163 */
164  HLL_DivBSource (state, Uhll, beg + 1, end, grid);
165 
166  #endif
167 }
#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
void HLL_DivBSource(const State_1D *state, double **Uhll, int beg, int end, Grid *grid)
Definition: source.c:179
double ** flux
upwind flux computed with the Riemann solver
Definition: structs.h:149
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
double ** GetBackgroundField(int beg, int end, int where, Grid *grid)
Definition: bckgrnd_field.c:5
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
#define FACE_CENTER
Definition: pluto.h:206
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
#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

Here is the call graph for this function: