PLUTO
tvdlf.c File Reference
#include "pluto.h"
Include dependency graph for tvdlf.c:

Go to the source code of this file.

Functions

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

Function Documentation

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

Definition at line 4 of file tvdlf.c.

39 {
40  int nv, i;
41  double *uL, *uR, cL, cR;
42  static double **fL, **fR;
43  static double *pR, *pL, *a2L, *a2R, *hL, *hR;
44  static double *cminL, *cmaxL;
45  static double *cminR, *cmaxR;
46 
47  if (fR == NULL){
48  fR = ARRAY_2D(NMAX_POINT, NFLX, double);
49  fL = ARRAY_2D(NMAX_POINT, NFLX, double);
50  pR = ARRAY_1D(NMAX_POINT, double);
51  pL = ARRAY_1D(NMAX_POINT, double);
52  cminL = ARRAY_1D(NMAX_POINT, double);
53  cmaxL = ARRAY_1D(NMAX_POINT, double);
54  cminR = ARRAY_1D(NMAX_POINT, double);
55  cmaxR = ARRAY_1D(NMAX_POINT, double);
56 
57  a2R = ARRAY_1D(NMAX_POINT, double);
58  a2L = ARRAY_1D(NMAX_POINT, double);
59 
60  hR = ARRAY_1D(NMAX_POINT, double);
61  hL = ARRAY_1D(NMAX_POINT, double);
62  }
63 
64 /* ----------------------------------------------------
65  compute sound speed & fluxes at zone interfaces
66  ---------------------------------------------------- */
67 
68  SoundSpeed2 (state->vL, a2L, hL, beg, end, FACE_CENTER, grid);
69  SoundSpeed2 (state->vR, a2R, hR, beg, end, FACE_CENTER, grid);
70 
71  Flux (state->uL, state->vL, a2L, fL, pL, beg, end);
72  Flux (state->uR, state->vR, a2R, fR, pR, beg, end);
73 
74  MaxSignalSpeed (state->vL, a2L, cminL, cmaxL, beg, end);
75  MaxSignalSpeed (state->vR, a2R, cminR, cmaxR, beg, end);
76 
77  for (i = beg; i <= end; i++) {
78  cL = MAX(fabs(cminL[i]), fabs(cmaxL[i]));
79  cR = MAX(fabs(cminR[i]), fabs(cmaxR[i]));
80  cmax[i] = MAX(cL, cR);
81 
82  uL = state->uL[i];
83  uR = state->uR[i];
84  for (nv = NFLX; nv--; ) {
85  state->flux[i][nv] = 0.5*(fL[i][nv] + fR[i][nv] - cmax[i]*(uR[nv] - uL[nv]));
86  }
87  state->press[i] = 0.5*(pL[i] + pR[i]);
88  }
89 }
#define MAX(a, b)
Definition: macros.h:101
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
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
#define FACE_CENTER
Definition: pluto.h:206
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
int i
Definition: analysis.c:2
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
void MaxSignalSpeed(double **v, double *cs2, double *cmin, double *cmax, int beg, int end)
Definition: eigenv.c:34
#define ARRAY_2D(nx, ny, type)
Definition: prototypes.h:171
double ** uL
same as vL, in conservative vars
Definition: structs.h:144

Here is the call graph for this function: