PLUTO
set_indexes.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief Perform index permutation and set domain integration indexes.
5 
6  The function SetIndexes() performs two different tasks, depending on
7  the value of ::g_dir and the time-stepping algorithm:
8  - perform a cyclic permutation of the array indexes corresponding to
9  vector components (velocity, momentum and magnetic field);
10  - set the starting and final grid indexes (IBEG, IEND,...) of the
11  computational domain before commencing integration.
12 
13  The indexes coincide with the usual ones (IBEG, IEND,...) most of the time,
14  but can be expanded one further zone either in the transverse or normal
15  directions or both. This depends on the integration algorithm.
16 
17  With cell-centered fields, the following table holds:
18  \verbatim
19  RK CTU
20  +-------------------------------------------
21  Transverse++ | NO YES, at predictor step
22  |
23  Normal++ | NO NO
24  \endverbatim
25 
26  If constrained transport is enabled, then
27  \verbatim
28  RK CTU
29  +-------------------------------------------
30  Transverse++ | YES YES
31  |
32  Normal++ | NO YES, at predictor step
33  \endverbatim
34 
35  Predictor step (g_intStage = 1) in CTU requires extra transverse loop
36  to obtain transverse predictor in any case.
37  Also, with CTU+CT, one needs to expand the grid of one zone in the
38  \e normal direction as well.
39  This allows to computed fully corner coupled states in the boundary to
40  get electric field components during the constrained transport algorithm.
41 
42  \author A. Mignone (mignone@ph.unito.it)\n
43  \date Sep 17, 2012
44 */
45 /* ///////////////////////////////////////////////////////////////////// */
46 #include "pluto.h"
47 
48 /* ********************************************************************* */
49 void SetIndexes (Index *indx, Grid *grid)
50 /*!
51  * Set vector indices and integration index range.
52  *
53  * \param [out] indx pointer to an Index structure
54  * \param [in] grid pointer to an array of Grid structures
55  *
56  *********************************************************************** */
57 {
58  IBEG = grid[IDIR].lbeg; IEND = grid[IDIR].lend;
59  JBEG = grid[JDIR].lbeg; JEND = grid[JDIR].lend;
60  KBEG = grid[KDIR].lbeg; KEND = grid[KDIR].lend;
61 
62  if (g_dir == IDIR) { /* -- Order: X-Y-Z {in,t1,t2 = i,j,k} -- */
63 
64  EXPAND(VXn = MXn = VX1; ,
65  VXt = MXt = VX2; ,
66  VXb = MXb = VX3;)
67 #if PHYSICS == MHD || PHYSICS == RMHD
68  EXPAND(BXn = BX1; ,
69  BXt = BX2; ,
70  BXb = BX3;)
71 #endif
72 #if DUST == YES
73  EXPAND(VXn_D = MXn_D = VX1_D; ,
74  VXt_D = MXt_D = VX2_D; ,
75  VXb_D = MXb_D = VX3_D;)
76 #endif
77 
78  indx->ntot = grid[IDIR].np_tot;
79  indx->beg = IBEG; indx->end = IEND;
80  indx->t1_beg = JBEG; indx->t1_end = JEND;
81  indx->t2_beg = KBEG; indx->t2_end = KEND;
82 
83  }else if (g_dir == JDIR){ /* -- Order: Y-X-Z {in,t1,t2 = j,i,k} -- */
84 
85  EXPAND(VXn = MXn = VX2; ,
86  VXt = MXt = VX1; ,
87  VXb = MXb = VX3;)
88 #if PHYSICS == MHD || PHYSICS == RMHD
89  EXPAND(BXn = BX2; ,
90  BXt = BX1; ,
91  BXb = BX3;)
92 #endif
93 #if DUST == YES
94  EXPAND(VXn_D = MXn_D = VX2_D; ,
95  VXt_D = MXt_D = VX1_D; ,
96  VXb_D = MXb_D = VX3_D;)
97 #endif
98 
99  indx->ntot = grid[JDIR].np_tot;
100  indx->beg = JBEG; indx->end = JEND;
101  indx->t1_beg = IBEG; indx->t1_end = IEND;
102  indx->t2_beg = KBEG; indx->t2_end = KEND;
103 
104  }else if (g_dir == KDIR){ /* -- Order: Z-X-Y {in,t1,t2 = k,i,j} -- */
105 
106  VXn = MXn = VX3;
107  VXt = MXt = VX1;
108  VXb = MXb = VX2;
109  #if PHYSICS == MHD || PHYSICS == RMHD
110  BXn = BX3;
111  BXt = BX1;
112  BXb = BX2;
113  #endif
114 #if DUST == YES
115  EXPAND(VXn_D = MXn_D = VX3_D; ,
116  VXt_D = MXt_D = VX1_D; ,
117  VXb_D = MXb_D = VX2_D;)
118 #endif
119 
120  indx->ntot = grid[KDIR].np_tot;
121  indx->beg = KBEG; indx->end = KEND;
122  indx->t1_beg = IBEG; indx->t1_end = IEND;
123  indx->t2_beg = JBEG; indx->t2_end = JEND;
124 
125  }
126 
127 /* -------------------------------------------------------
128  Expand grid one further zone to account for proper
129  flux computation. This is necessary to obtain the EMF
130  in the boundary zones and to get transverse rhs for
131  corner coupled states.
132  ------------------------------------------------------- */
133 
134  #ifdef STAGGERED_MHD
135  D_EXPAND( ; ,
136  indx->t1_beg--; indx->t1_end++; ,
137  indx->t2_beg--; indx->t2_end++;)
138  #ifdef CTU
139  if (g_intStage == 1){
140  indx->beg--;
141  indx->end++;
142  D_EXPAND( ; ,
143  indx->t1_beg--; indx->t1_end++; ,
144  indx->t2_beg--; indx->t2_end++;)
145  }
146  #endif
147  #else
148  #ifdef CTU
149  if (g_intStage == 1){
150  #if (PARABOLIC_FLUX & EXPLICIT)
151  indx->beg--;
152  indx->end++;
153  #endif
154  D_EXPAND( ,
155  indx->t1_beg--; indx->t1_end++; ,
156  indx->t2_beg--; indx->t2_end++;)
157  }
158  #endif
159  #endif
160 }
161 
162 /* ********************************************************************* */
163 void ResetState (const Data *d, State_1D *state, Grid *grid)
164 /*!
165  * Initialize some of the elements of the State_1D structure to zero
166  * in order to speed up computations. These includes:
167  *
168  * - source term
169  * - left and right eigenvectors
170  * - the maximum eigenvalue ???
171  *
172  * \param [in] d pointer to Data structure
173  * \param [out] state pointer to a State_1D structure
174  * \param [in] grid pointer to an array of Grid structures
175  *
176  *********************************************************************** */
177 {
178  int i, j, k, nv;
179  double v[NVAR], lambda[NVAR];
180  double a;
181 
182  for (i = 0; i < NMAX_POINT; i++){
183  for (j = NVAR; j--; ) state->src[i][j] = 0.0;
184 
185  for (j = NFLX; j--; ){
186  for (k = NFLX; k--; ){
187  state->Lp[i][j][k] = state->Rp[i][j][k] = 0.0;
188  }}
189  }
190 
191 /* ---------------------------------------------------------
192  When using Finite Difference methods, we need to find,
193  for each characteristic k, its maximum over the
194  whole computational domain (LF global splitting).
195  --------------------------------------------------------- */
196 
197  #ifdef FINITE_DIFFERENCE
198  FD_GetMaxEigenvalues (d, state, grid);
199  #endif
200 }
static double a
Definition: init.c:135
#define VX2
Definition: mod_defs.h:29
void FD_GetMaxEigenvalues(const Data *d, State_1D *state, Grid *grid)
Definition: fd_flux.c:293
int MXn_D
Definition: globals.h:78
int ntot
Definition: structs.h:318
int VXn_D
Definition: globals.h:77
double *** Lp
Definition: structs.h:155
int t2_end
Definition: structs.h:320
#define YES
Definition: pluto.h:25
int t1_end
Definition: structs.h:319
int VXt_D
Definition: globals.h:77
int g_intStage
Gives the current integration stage of the time stepping method (predictor = 0, 1st corrector = 1...
Definition: globals.h:98
int lend
Local end index for the local array.
Definition: structs.h:118
int lbeg
Local start index for the local array.
Definition: structs.h:117
#define VX1
Definition: mod_defs.h:28
#define KDIR
Definition: pluto.h:195
int VXb_D
Definition: globals.h:77
int t2_beg
Definition: structs.h:320
int VXb
Definition: globals.h:73
int BXn
Definition: globals.h:75
int t1_beg
Definition: structs.h:319
int MXt
Definition: globals.h:74
int MXb_D
Definition: globals.h:78
double ** src
Definition: structs.h:160
void ResetState(const Data *d, State_1D *state, Grid *grid)
Definition: set_indexes.c:163
#define NFLX
Definition: mod_defs.h:32
#define IDIR
Definition: pluto.h:193
int g_dir
Specifies the current sweep or direction of integration.
Definition: globals.h:86
Definition: structs.h:78
#define PHYSICS
Definition: definitions_01.h:1
int j
Definition: analysis.c:2
void SetIndexes(Index *indx, Grid *grid)
Definition: set_indexes.c:49
int VXt
Definition: globals.h:73
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
int VXn
Definition: globals.h:73
#define DUST
Definition: pluto.h:333
int beg
Definition: structs.h:318
#define MHD
Definition: pluto.h:111
#define BX3
Definition: mod_defs.h:27
PLUTO main header file.
D_EXPAND(tot/[n]=(double) grid[IDIR].np_int_glob;, tot/[n]=(double) grid[JDIR].np_int_glob;, tot/[n]=(double) grid[KDIR].np_int_glob;)
Definition: analysis.c:27
Definition: structs.h:30
long int NMAX_POINT
Maximum number of points among the three directions, boundaries excluded.
Definition: globals.h:62
int i
Definition: analysis.c:2
Definition: structs.h:317
#define BX1
Definition: mod_defs.h:25
#define VX3
Definition: mod_defs.h:30
int MXt_D
Definition: globals.h:78
long int KBEG
Lower grid index of the computational domain in the the X3 direction for the local processor...
Definition: globals.h:43
double *** Rp
Left and right primitive eigenvectors.
Definition: structs.h:155
#define BX2
Definition: mod_defs.h:26
int np_tot
Total number of points in the local domain (boundaries included).
Definition: structs.h:100
long int KEND
Upper grid index of the computational domain in the the X3 direction for the local processor...
Definition: globals.h:45
#define JDIR
Definition: pluto.h:194
int BXt
Definition: globals.h:75
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
int BXb
Definition: globals.h:75
#define RMHD
Definition: pluto.h:112
long int JEND
Upper grid index of the computational domain in the the X2 direction for the local processor...
Definition: globals.h:41
long int IBEG
Lower grid index of the computational domain in the the X1 direction for the local processor...
Definition: globals.h:35
int MXb
Definition: globals.h:74