PLUTO
update_stage.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief Single stage integration for RK time stepping.
5 
6  Advance the equations in conservative form by taking a single stage
7  in the form
8  \f[
9  U \quad\Longrightarrow \quad U + \Delta t R(V)
10  \f]
11  where \c U is a 3D array of conservative variables, \c V is a 3D array
12  of primitive variables, \c R(V) is the right hand side containing
13  flux differences and source terms.
14  Note that \c U and \c V may \e not necessarily be the map of
15  each other, i.e., \c U is \e not \c U(V).
16  The right hand side can contain contributions from
17 
18  - the direction set by the global variable ::g_dir,
19  when DIMENSIONAL_SPLITTING == YES;
20  - all directions when DIMENSIONAL_SPLITTING == NO;
21 
22  When the integrator stage is the first one (predictor), this function
23  also computes the maximum of inverse time steps for hyperbolic and
24  parabolic terms (if the latters are included explicitly).
25 
26  \authors A. Mignone (mignone@ph.unito.it)\n
27  C. Zanni (zanni@oato.inaf.it)\n
28  P. Tzeferacos (petros.tzeferacos@ph.unito.it)
29  T. Matsakos
30  \date March 02, 2014
31 */
32 /* ///////////////////////////////////////////////////////////////////// */
33 #include "pluto.h"
34 static void SaveAMRFluxes (const State_1D *, double **, int, int, Grid *);
35 static intList TimeStepIndexList();
36 
37 /* ********************************************************************* */
38 void UpdateStage(const Data *d, Data_Arr UU, double **aflux,
39  Riemann_Solver *Riemann, double dt, Time_Step *Dts,
40  Grid *grid)
41 /*!
42  *
43  * \param [in,out] d pointer to PLUTO Data structure
44  * \param [in,out] UU data array containing conservative variables
45  * at the previous time step to be updated
46  * \param [out] aflux interface fluxes needed for refluxing operations
47  * (only with AMR)
48  * \param [in] Riemann pointer to a Riemann solver function
49  * \param [in] dt the time step for the current update step
50  * \param [in,out] Dts pointer to time step structure
51  * \param [in] grid pointer to array of Grid structures
52  *********************************************************************** */
53 {
54  int i, j, k;
55  int nv, dir, beg_dir, end_dir;
56  int *ip;
57  double *inv_dl, dl2;
58  static double ***T, ***C_dt[NVAR], **dcoeff;
59  static State_1D state;
60  Index indx;
61  intList cdt_list;
62 
63  #if DIMENSIONAL_SPLITTING == YES
64  beg_dir = end_dir = g_dir;
65  #else
66  beg_dir = 0;
67  end_dir = DIMENSIONS-1;
68  #endif
69 
70  cdt_list = TimeStepIndexList();
71 
72 /* --------------------------------------------------------------
73  1. Allocate memory
74  -------------------------------------------------------------- */
75 
76  if (state.v == NULL){
77  MakeState (&state);
78  #if (PARABOLIC_FLUX & EXPLICIT)
79  dcoeff = ARRAY_2D(NMAX_POINT, NVAR, double);
80  #endif
81  }
82 
83 /* --------------------------------------------------------------
84  2. Reset arrays.
85  C_dt is an array used to store the inverse time step for
86  advection and diffusion.
87  We use C_dt[RHO] for advection,
88  C_dt[MX1] for viscosity,
89  C_dt[BX1...BX3] for resistivity and
90  C_dt[ENG] for thermal conduction.
91  --------------------------------------------------------------- */
92 
93  #if DIMENSIONAL_SPLITTING == NO
94  if (C_dt[RHO] == NULL){
95  FOR_EACH(nv, 0, (&cdt_list)) {
96  C_dt[nv] = ARRAY_3D(NX3_MAX, NX2_MAX, NX1_MAX, double);
97  }
98  }
99 
100  if (g_intStage == 1) KTOT_LOOP(k) JTOT_LOOP(j){
101  FOR_EACH(nv, 0, (&cdt_list)) {
102  memset ((void *)C_dt[nv][k][j],'\0', NX1_TOT*sizeof(double));
103  }
104  }
105  #endif
106 
107 /* ------------------------------------------------
108  2a. Compute current arrays
109  ------------------------------------------------ */
110 
111  #if (RESISTIVITY == EXPLICIT) && (defined STAGGERED_MHD)
112  GetCurrent (d, -1, grid);
113  #endif
114 
115 /* ------------------------------------------------
116  2b. Compute Temperature array
117  ------------------------------------------------ */
118 
119  #if THERMAL_CONDUCTION == EXPLICIT
120  if (T == NULL) T = ARRAY_3D(NX3_MAX, NX2_MAX, NX1_MAX, double);
121  TOT_LOOP(k,j,i) T[k][j][i] = d->Vc[PRS][k][j][i]/d->Vc[RHO][k][j][i];
122  #endif
123 
124 /* ----------------------------------------------------------------
125  3. Main loop on directions
126  ---------------------------------------------------------------- */
127 
128  for (dir = beg_dir; dir <= end_dir; dir++){
129 
130  g_dir = dir;
131  SetIndexes (&indx, grid); /* -- set normal and transverse indices -- */
132  ResetState (d, &state, grid);
133 
134  #if (RESISTIVITY == EXPLICIT) && !(defined STAGGERED_MHD)
135  GetCurrent(d, dir, grid);
136  #endif
137 
138  TRANSVERSE_LOOP(indx,ip,i,j,k){
139  g_i = i; g_j = j; g_k = k;
140  for ((*ip) = 0; (*ip) < indx.ntot; (*ip)++) {
141  VAR_LOOP(nv) state.v[(*ip)][nv] = d->Vc[nv][k][j][i];
142  state.flag[*ip] = d->flag[k][j][i];
143  #ifdef STAGGERED_MHD
144  state.bn[(*ip)] = d->Vs[g_dir][k][j][i];
145  #endif
146  }
147  CheckNaN (state.v, 0, indx.ntot-1,0);
148  States (&state, indx.beg - 1, indx.end + 1, grid);
149  Riemann (&state, indx.beg - 1, indx.end, Dts->cmax, grid);
150  #ifdef STAGGERED_MHD
151  CT_StoreEMF (&state, indx.beg - 1, indx.end, grid);
152  #endif
153  #if (PARABOLIC_FLUX & EXPLICIT)
154  ParabolicFlux(d->Vc, d->J, T, &state, dcoeff, indx.beg-1, indx.end, grid);
155  #endif
156  #if UPDATE_VECTOR_POTENTIAL == YES
157  VectorPotentialUpdate (d, NULL, &state, grid);
158  #endif
159  #ifdef SHEARINGBOX
160  SB_SaveFluxes (&state, grid);
161  #endif
162  RightHandSide (&state, Dts, indx.beg, indx.end, dt, grid);
163 
164  /* -- update: U = U + dt*R -- */
165 
166  #ifdef CHOMBO
167  for ((*ip) = indx.beg; (*ip) <= indx.end; (*ip)++) {
168  VAR_LOOP(nv) UU[nv][k][j][i] += state.rhs[*ip][nv];
169  }
170  SaveAMRFluxes (&state, aflux, indx.beg-1, indx.end, grid);
171  #else
172  for ((*ip) = indx.beg; (*ip) <= indx.end; (*ip)++) {
173  VAR_LOOP(nv) UU[k][j][i][nv] += state.rhs[*ip][nv];
174  }
175  #endif
176 
177  if (g_intStage > 1) continue;
178 
179  /* -- compute inverse dt coefficients when g_intStage = 1 -- */
180 
181  inv_dl = GetInverse_dl(grid);
182  for ((*ip) = indx.beg; (*ip) <= indx.end; (*ip)++) {
183  #if DIMENSIONAL_SPLITTING == NO
184 
185  #if !GET_MAX_DT
186  C_dt[0][k][j][i] += 0.5*( Dts->cmax[(*ip)-1]
187  + Dts->cmax[*ip])*inv_dl[*ip];
188  #endif
189  #if (PARABOLIC_FLUX & EXPLICIT)
190  dl2 = 0.5*inv_dl[*ip]*inv_dl[*ip];
191  FOR_EACH(nv, 1, (&cdt_list)) {
192  C_dt[nv][k][j][i] += (dcoeff[*ip][nv]+dcoeff[(*ip)-1][nv])*dl2;
193  }
194  #endif
195 
196  #elif DIMENSIONAL_SPLITTING == YES
197 
198  #if !GET_MAX_DT
199  Dts->inv_dta = MAX(Dts->inv_dta, Dts->cmax[*ip]*inv_dl[*ip]);
200  #endif
201  #if (PARABOLIC_FLUX & EXPLICIT)
202  dl2 = inv_dl[*ip]*inv_dl[*ip];
203  FOR_EACH(nv, 1, (&cdt_list)) {
204  Dts->inv_dtp = MAX(Dts->inv_dtp, dcoeff[*ip][nv]*dl2);
205  }
206  #endif
207  #endif
208  }
209  }
210  }
211 
212 /* -------------------------------------------------------------------
213  4. Additional terms here
214  ------------------------------------------------------------------- */
215 
216  #if (ENTROPY_SWITCH) && (RESISTIVITY == EXPLICIT)
217  EntropyOhmicHeating(d, UU, dt, grid);
218  #endif
219 
220  #ifdef SHEARINGBOX
221  SB_CorrectFluxes (UU, 0.0, dt, grid);
222  #endif
223 
224  #ifdef STAGGERED_MHD
225  CT_Update(d, d->Vs, dt, grid);
226  #endif
227 
228  #if DIMENSIONAL_SPLITTING == YES
229  return;
230  #endif
231 
232 /* -------------------------------------------------------------------
233  5. Reduce dt for dimensionally unsplit schemes.
234  ------------------------------------------------------------------- */
235 
236  if (g_intStage > 1) return;
237 
238  for (k = KBEG; k <= KEND; k++){ g_k = k;
239  for (j = JBEG; j <= JEND; j++){ g_j = j;
240  for (i = IBEG; i <= IEND; i++){
241  #if !GET_MAX_DT
242  Dts->inv_dta = MAX(Dts->inv_dta, C_dt[0][k][j][i]);
243  #endif
244  #if (PARABOLIC_FLUX & EXPLICIT)
245  FOR_EACH(nv, 1, (&cdt_list)) {
246  Dts->inv_dtp = MAX(Dts->inv_dtp, C_dt[nv][k][j][i]);
247  }
248  #endif
249  }
250  }}
251  #if !GET_MAX_DT
252  Dts->inv_dta /= (double)DIMENSIONS;
253  #endif
254  #if (PARABOLIC_FLUX & EXPLICIT)
255  Dts->inv_dtp /= (double)DIMENSIONS;
256  #endif
257 
258 }
259 
260 #ifdef CHOMBO
261 /* ********************************************************************* */
262 void SaveAMRFluxes (const State_1D *state, double **aflux,
263  int beg, int end, Grid *grid)
264 /*!
265  * Rebuild fluxes in a way suitable for AMR operation
266  * by adding pressure and multiplying by area.
267  *
268  *********************************************************************** */
269 {
270  int i, j, k, nv, *in;
271  int nxf, nyf, nzf;
272  int nxb, nyb, nzb;
273  long int indf;
274  double wflux, r, area;
275 
276  for (i = beg; i <= end; i++) state->flux[i][MXn] += state->press[i];
277 
278 #if (GEOMETRY == CARTESIAN) && (CH_SPACEDIM > 1)
279  if ((g_dir == IDIR) && (g_stretch_fact != 1.)) {
280  for (i = beg; i <= end; i++) {
281  NVAR_LOOP(nv) state->flux[i][nv] *= g_stretch_fact;
282  }
283  }
284  #if (CH_SPACEDIM == 3)
285  if ((g_dir == JDIR) && (g_x3stretch != 1.)) {
286  for (i = beg; i <= end; i++) {
287  NVAR_LOOP(nv) state->flux[i][nv] *= g_x3stretch;
288  }
289  }
290  if ((g_dir == KDIR) && (g_x2stretch != 1.)) {
291  for (i = beg; i <= end; i++) {
292  NVAR_LOOP(nv) state->flux[i][nv] *= g_x2stretch;
293  }
294  }
295  #endif
296 #endif
297 
298 #if GEOMETRY == CYLINDRICAL
299  if (g_dir == IDIR){
300  for (i = beg; i <= end; i++) {
301  NVAR_LOOP(nv) {
302  state->flux[i][nv] *= grid[IDIR].A[i];
303  #if CH_SPACEDIM > 1
304  state->flux[i][nv] *= g_x2stretch;
305  #endif
306  }
307  }
308  }else{
309  area = fabs(grid[IDIR].x[g_i]);
310  for (i = beg; i <= end; i++) {
311  NVAR_LOOP(nv) state->flux[i][nv] *= area;
312  }
313  }
314 #endif
315 
316 #if GEOMETRY == SPHERICAL
317  if (g_dir == IDIR){
318  #if CH_SPACEDIM > 1
319  area = grid[JDIR].dV[g_j]/g_level_dx;
320  #if CH_SPACEDIM == 3
321  area *= g_x3stretch;
322  #endif
323  #endif
324  for (i = beg; i <= end; i++) {
325  NVAR_LOOP(nv) {
326  state->flux[i][nv] *= grid[IDIR].A[i];
327  #if CH_SPACEDIM > 1
328  state->flux[i][nv] *= area;
329  #endif
330  }
331  #if (COMPONENTS == 3) && CHOMBO_CONS_AM
332  state->flux[i][iMPHI] *= grid[IDIR].xr[i]*sin(grid[JDIR].x[g_j]);
333  #endif
334  }
335  }
336  if (g_dir == JDIR){
337  area = fabs(grid[IDIR].x[g_i]);
338  #if CHOMBO_LOGR == YES
339  area *= grid[IDIR].dx[g_i]/g_level_dx;
340  #endif
341  #if CH_SPACEDIM == 3
342  area *= g_x3stretch;
343  #endif
344  for (i = beg; i <= end; i++) {
345  NVAR_LOOP(nv) {
346  state->flux[i][nv] *= grid[JDIR].A[i]*area;
347  }
348  #if (COMPONENTS == 3) && CHOMBO_CONS_AM
349  state->flux[i][iMPHI] *= grid[IDIR].x[g_i]*sin(grid[JDIR].xr[i]);
350  #endif
351  }
352  }
353  if (g_dir == KDIR){
354  for (i = beg; i <= end; i++) {
355  area = g_x2stretch* fabs(grid[IDIR].x[g_i]);
356  #if CHOMBO_LOGR == YES
357  area *= grid[IDIR].dx[g_i]/g_level_dx;
358  #endif
359  NVAR_LOOP(nv) {
360  state->flux[i][nv] *= area;
361  }
362  #if (COMPONENTS == 3) && CHOMBO_CONS_AM
363  state->flux[i][iMPHI] *= grid[IDIR].x[g_i]*sin(grid[JDIR].x[g_j]);
364  #endif
365  }
366  }
367 #endif
368 
369 #if GEOMETRY == POLAR
370  if (g_dir == IDIR){
371  for (i = beg; i <= end; i++) {
372  NVAR_LOOP(nv) {
373  state->flux[i][nv] *= grid[IDIR].A[i];
374  #if CH_SPACEDIM > 1
375  state->flux[i][nv] *= g_x2stretch;
376  #endif
377  #if CH_SPACEDIM == 3
378  state->flux[i][nv] *= g_x3stretch;
379  #endif
380  }
381  #if (COMPONENTS > 1) && CHOMBO_CONS_AM
382  state->flux[i][iMPHI] *= grid[IDIR].xr[i];
383  #endif
384  }
385  }
386  if (g_dir == JDIR){
387  area = g_x3stretch;
388  #if CHOMBO_LOGR == YES
389  area *= grid[IDIR].dx[g_i]/g_level_dx;
390  #endif
391  #if CH_SPACEDIM == 3
392  area *= g_x3stretch;
393  #endif
394  if (area != 1.) {
395  for (i = beg; i <= end; i++) {
396  NVAR_LOOP(nv) {
397  state->flux[i][nv] *= area;
398  }
399  }
400  }
401  #if (COMPONENTS > 1) && CHOMBO_CONS_AM
402  for (i = beg; i <= end; i++) state->flux[i][iMPHI] *= grid[IDIR].x[g_i];
403  #endif
404  }
405  if (g_dir == KDIR){
406  for (i = beg; i <= end; i++) {
407  area = g_x2stretch*fabs(grid[IDIR].x[g_i]);
408  #if CHOMBO_LOGR == YES
409  area *= grid[IDIR].dx[g_i]/g_level_dx;
410  #endif
411  NVAR_LOOP(nv) {
412  state->flux[i][nv] *= area;
413  }
414  #if (COMPONENTS > 1) && CHOMBO_CONS_AM
415  state->flux[i][iMPHI] *= grid[IDIR].x[g_i];
416  #endif
417  }
418  }
419 #endif
420 
421 /* -- store fluxes for re-fluxing operation -- */
422 
423  nxf = grid[IDIR].np_int + (g_dir == IDIR);
424  nyf = grid[JDIR].np_int + (g_dir == JDIR);
425  nzf = grid[KDIR].np_int + (g_dir == KDIR);
426 
427  nxb = grid[IDIR].lbeg - (g_dir == IDIR);
428  nyb = grid[JDIR].lbeg - (g_dir == JDIR);
429  nzb = grid[KDIR].lbeg - (g_dir == KDIR);
430 
431 #if TIME_STEPPING == RK2
432  wflux = 0.5;
433 #else
434  wflux = 1.0;
435 #endif
436 
437  i = g_i; j = g_j; k = g_k;
438  if (g_dir == IDIR) in = &i;
439  if (g_dir == JDIR) in = &j;
440  if (g_dir == KDIR) in = &k;
441 
442  for ((*in) = beg; (*in) <= end; (*in)++) {
443  #if HAVE_ENERGY && ENTROPY_SWITCH
444  state->flux[*in][ENTR] = 0.0;
445  #endif
446  NVAR_LOOP(nv) {
447  indf = nv*nzf*nyf*nxf + (k - nzb)*nyf*nxf + (j - nyb)*nxf + (i - nxb);
448  aflux[g_dir][indf] = wflux*state->flux[(*in)][nv];
449  }
450  }
451 
452 }
453 #endif
454 
455 /* ********************************************************************* */
457 /*!
458  * Return the intList of inverse time step indices.
459  *
460  *********************************************************************** */
461 {
462  int i = 0;
463  intList cdt;
464 
465  cdt.indx[i++] = RHO;
466  #if VISCOSITY == EXPLICIT
467  cdt.indx[i++] = MX1;
468  #endif
469  #if RESISTIVITY == EXPLICIT
470  EXPAND(cdt.indx[i++] = BX1; ,
471  cdt.indx[i++] = BX2; ,
472  cdt.indx[i++] = BX3;)
473  #endif
474  #if THERMAL_CONDUCTION == EXPLICIT
475  cdt.indx[i++] = ENG;
476  #endif
477  cdt.nvar = i;
478 
479  return cdt;
480 }
static void SaveAMRFluxes(const State_1D *, double **, int, int, Grid *)
#define MAX(a, b)
Definition: macros.h:101
#define iMPHI
Definition: mod_defs.h:71
tuple T
Definition: Sph_disk.py:33
void CT_StoreEMF(const State_1D *state, int beg, int end, Grid *grid)
Definition: ct_emf.c:35
double ** v
Cell-centered primitive varables at the base time level, v[i] = .
Definition: structs.h:134
double **** Vs
The main four-index data array used for face-centered staggered magnetic fields.
Definition: structs.h:43
void ResetState(const Data *, State_1D *, Grid *)
Definition: set_indexes.c:163
#define MX1
Definition: mod_defs.h:20
void EntropyOhmicHeating(const Data *, Data_Arr, double, Grid *)
double **** Data_Arr
Definition: pluto.h:492
double **** J
Electric current defined as curl(B).
Definition: structs.h:54
void GetCurrent(const Data *d, int dir, Grid *grid)
Definition: res_functions.c:97
double ** flux
upwind flux computed with the Riemann solver
Definition: structs.h:149
int ntot
Definition: structs.h:318
void Riemann_Solver(const State_1D *, int, int, double *, Grid *)
Definition: pluto.h:489
double ** rhs
Conservative right hand side.
Definition: structs.h:163
double * xr
Definition: structs.h:81
void States(const State_1D *, int, int, Grid *)
Definition: plm_states.c:430
#define RHO
Definition: mod_defs.h:19
int g_intStage
Gives the current integration stage of the time stepping method (predictor = 0, 1st corrector = 1...
Definition: globals.h:98
double **** Vc
The main four-index data array used for cell-centered primitive variables.
Definition: structs.h:31
#define KTOT_LOOP(k)
Definition: macros.h:40
int lbeg
Local start index for the local array.
Definition: structs.h:117
double * GetInverse_dl(const Grid *)
Definition: set_geometry.c:205
void SB_CorrectFluxes(Data_Arr U, double t, double dt, Grid *grid)
Definition: sb_flux.c:105
double * dV
Cell volume.
Definition: structs.h:86
void CT_Update(const Data *d, Data_Arr Bs, double dt, Grid *grid)
Definition: ct.c:29
#define ARRAY_3D(nx, ny, nz, type)
Definition: prototypes.h:172
double * dx
Definition: structs.h:83
double * cmax
Maximum signal velocity for hyperbolic eqns.
Definition: structs.h:214
#define FOR_EACH(nv, beg, list)
Definition: macros.h:89
#define KDIR
Definition: pluto.h:195
unsigned char *** flag
Pointer to a 3D array setting useful integration flags that are retrieved during integration.
Definition: structs.h:55
int nvar
Number of variables.
Definition: structs.h:330
void MakeState(State_1D *)
Definition: tools.c:51
double inv_dtp
Inverse of diffusion (parabolic) time step .
Definition: structs.h:217
#define NX3_MAX
Definition: pluto.h:743
double inv_dta
Inverse of advection (hyperbolic) time step, .
Definition: structs.h:215
int g_i
x1 grid index when sweeping along the x2 or x3 direction.
Definition: globals.h:82
int CheckNaN(double **, int, int, int)
Definition: tools.c:19
#define VAR_LOOP(n)
Definition: macros.h:226
#define NX2_MAX
Definition: pluto.h:742
#define TOT_LOOP(k, j, i)
Definition: macros.h:44
#define IDIR
Definition: pluto.h:193
int g_dir
Specifies the current sweep or direction of integration.
Definition: globals.h:86
#define NVAR_LOOP(n)
Definition: pluto.h:618
int g_j
x2 grid index when sweeping along the x1 or x3 direction.
Definition: globals.h:83
Definition: structs.h:78
int j
Definition: analysis.c:2
if(divB==NULL)
Definition: analysis.c:10
unsigned char * flag
Definition: structs.h:168
long int IEND
Upper grid index of the computational domain in the the X1 direction for the local processor...
Definition: globals.h:37
int k
Definition: analysis.c:2
int MXn
Definition: globals.h:74
void SetIndexes(Index *indx, Grid *grid)
Definition: set_indexes.c:49
double * x
Definition: structs.h:80
void RightHandSide(const State_1D *state, Time_Step *Dts, int beg, int end, double dt, Grid *grid)
Definition: rhs.c:88
int g_k
x3 grid index when sweeping along the x1 or x2 direction.
Definition: globals.h:84
#define GEOMETRY
Definition: definitions_01.h:4
#define CARTESIAN
Definition: pluto.h:33
int beg
Definition: structs.h:318
#define BX3
Definition: mod_defs.h:27
PLUTO main header file.
Definition: structs.h:30
long int NMAX_POINT
Maximum number of points among the three directions, boundaries excluded.
Definition: globals.h:62
#define JTOT_LOOP(j)
Definition: macros.h:39
int i
Definition: analysis.c:2
int indx[2046]
Array of integers containg variables indices.
Definition: structs.h:329
double * bn
Face magentic field, bn = bx(i+1/2)
Definition: structs.h:165
Definition: structs.h:317
void SB_SaveFluxes(State_1D *state, Grid *grid)
Definition: sb_flux.c:57
#define BX1
Definition: mod_defs.h:25
void ParabolicFlux(Data_Arr V, Data_Arr J, double ***T, const State_1D *state, double **dcoeff, int beg, int end, Grid *grid)
double * press
Upwind pressure term computed with the Riemann solver.
Definition: structs.h:164
long int KBEG
Lower grid index of the computational domain in the the X3 direction for the local processor...
Definition: globals.h:43
#define ARRAY_2D(nx, ny, type)
Definition: prototypes.h:171
static intList TimeStepIndexList()
Definition: update_stage.c:456
#define BX2
Definition: mod_defs.h:26
long int KEND
Upper grid index of the computational domain in the the X3 direction for the local processor...
Definition: globals.h:45
void VectorPotentialUpdate(const Data *d, const void *vp, const State_1D *state, const Grid *grid)
#define JDIR
Definition: pluto.h:194
int end
Definition: structs.h:318
long int JBEG
Lower grid index of the computational domain in the the X2 direction for the local processor...
Definition: globals.h:39
#define NVAR
Definition: pluto.h:609
void UpdateStage(const Data *d, Data_Arr UU, double **aflux, Riemann_Solver *Riemann, double dt, Time_Step *Dts, Grid *grid)
Definition: update_stage.c:38
#define NX1_MAX
Definition: pluto.h:741
long int JEND
Upper grid index of the computational domain in the the X2 direction for the local processor...
Definition: globals.h:41
#define TRANSVERSE_LOOP(indx, ip, i, j, k)
Definition: macros.h:54
long int NX1_TOT
Total number of zones in the X1 direction (boundaries included) for the local processor.
Definition: globals.h:55
long int IBEG
Lower grid index of the computational domain in the the X1 direction for the local processor...
Definition: globals.h:35
double * A
Right interface area, A[i] = .
Definition: structs.h:87
#define DIMENSIONS
Definition: definitions_01.h:2
int np_int
Total number of points in the local domain (boundaries excluded).
Definition: structs.h:102