PLUTO
source.c
Go to the documentation of this file.
1 #include "pluto.h"
2 
3 #if DIVB_CONTROL == EIGHT_WAVES
4 /* ************************************************************************* */
5 void POWELL_DIVB_SOURCE(const State_1D *state, int is, int ie, Grid *grid)
6 /*
7  *
8  *
9  ************************************************************************** */
10 {
11  int i;
12  double divB;
13  double g_2, vB, *vc;
14  Grid *GG;
15 
16  GG = grid + g_dir;
17 
18 /* ------------------------------------------------------------
19  POWELL's divB monopole source term
20  ------------------------------------------------------------ */
21 
22  for (i = is; i <= ie; i++) {
23 
24  vc = state->vh[i];
25 
26 /* Make (divB) contribution from the normal direction (1) */
27 
28  divB = (state->vR[i][BXn] + state->vL[i][BXn])*GG->A[i] -
29  (state->vR[i - 1][BXn] + state->vL[i - 1][BXn])*GG->A[i - 1];
30 /*
31  divB = state->bn[i]*GG->A[i] - state->bn[i - 1]*GG->A[i - 1];
32 */
33  divB /= 2.0*GG->dV[i];
34 
35  vB = EXPAND(vc[VX1]*vc[BX1], + vc[VX2]*vc[BX2], + vc[VX3]*vc[BX3]);
36  g_2 = EXPAND(vc[VX1]*vc[VX1], + vc[VX2]*vc[VX2], + vc[VX3]*vc[VX3]);
37  g_2 = 1.0 - g_2;
38 
39  EXPAND (state->src[i][MX1] = -divB*(vc[BX1]*g_2 + vB*vc[VX1]); ,
40  state->src[i][MX2] = -divB*(vc[BX2]*g_2 + vB*vc[VX2]); ,
41  state->src[i][MX3] = -divB*(vc[BX3]*g_2 + vB*vc[VX3]);)
42 
43  state->src[i][ENG] = -divB*vB;
44 
45  EXPAND (state->src[i][BX1] = -divB*vc[VX1]; ,
46  state->src[i][BX2] = -divB*vc[VX2]; ,
47  state->src[i][BX3] = -divB*vc[VX3];)
48 
49  }
50 }
51 
52 /* ********************************************************************* */
53 void HLL_DIVB_SOURCE (const State_1D *state, double **Uhll, int beg,
54  int end, Grid *grid)
55 /*
56  *
57  * PURPOSE
58  *
59  * Include div.B source term to momentum, induction
60  * and energy equation. Used in conjunction with
61  * an HLL-type Riemann solver.
62  *
63  * LAST_MODIFIED
64  *
65  * Setp 18th 2008, by Andrea Mignone (mignone@to.astro.it)
66  *
67  *********************************************************************** */
68 {
69  int i, nv;
70  double vc[NVAR], *A, *src, *vm;
71  double r, s, vB;
72  static double *divB, *vp;
73  Grid *GG;
74 
75  if (divB == NULL){
76  divB = ARRAY_1D(NMAX_POINT, double);
77  vp = ARRAY_1D(NMAX_POINT, double);
78  }
79 
80  GG = grid + g_dir;
81  vm = vp - 1;
82  A = grid[g_dir].A;
83 
84 /* --------------------------------------------
85  Compute normal component of the field
86  -------------------------------------------- */
87 
88  for (i = beg - 1; i <= end; i++) {
89  vp[i] = state->flux[i][RHO] < 0.0 ? state->vR[i][BXn]: state->vL[i][BXn];
90  }
91 
92 /* --------------------------------------------
93  Compute div.B contribution from the normal
94  direction (1) in different geometries
95  -------------------------------------------- */
96 
97 
98  #if GEOMETRY == CARTESIAN
99 
100  for (i = beg; i <= end; i++) {
101  divB[i] = (vp[i] - vm[i])/GG->dx[i];
102  }
103 
104  #elif GEOMETRY == CYLINDRICAL
105 
106  if (g_dir == IDIR){ /* -- r -- */
107  for (i = beg; i <= end; i++) {
108  divB[i] = (vp[i]*A[i] - vm[i]*A[i - 1])/GG->dV[i];
109  }
110  }else if (g_dir == JDIR){ /* -- z -- */
111  for (i = beg; i <= end; i++) {
112  divB[i] = (vp[i] - vm[i])/GG->dx[i];
113  }
114  }
115 
116  #elif GEOMETRY == POLAR
117 
118  if (g_dir == IDIR){ /* -- r -- */
119  Ar = grid[IDIR].A;
120  for (i = beg; i <= end; i++) {
121  divB[i] = (vp[i]*A[i] - vm[i]*A[i - 1])/GG->dV[i];
122  }
123  }else if (g_dir == JDIR){ /* -- phi -- */
124  r = grid[IDIR].x[*g_i];
125  for (i = beg; i <= end; i++) {
126  divB[i] = (vp[i] - vm[i])/(r*GG->dx[i]);
127  }
128  }else if (g_dir == KDIR){ /* -- z -- */
129  for (i = beg; i <= end; i++) {
130  divB[i] = (vp[i] - vm[i])/GG->dx[i];
131  }
132  }
133 
134  #elif GEOMETRY == SPHERICAL
135 
136  if (g_dir == IDIR){ /* -- r -- */
137  for (i = beg; i <= end; i++) {
138  divB[i] = (vp[i]*A[i] - vm[i]*A[i - 1])/GG->dV[i];
139  }
140  }else if (g_dir == JDIR){ /* -- theta -- */
141  r = grid[IDIR].x[*g_i];
142  for (i = beg; i <= end; i++) {
143  divB[i] = (vp[i]*A[i] - vm[i]*A[i - 1])/(r*GG->dV[i]);
144  }
145  }else if (g_dir == KDIR){ /* -- phi -- */
146  r = grid[IDIR].x[*g_i];
147  s = sin(grid[JDIR].x[*g_j]);
148  for (i = beg; i <= end; i++) {
149  divB[i] = (vp[i] - vm[i])/(r*s*GG->dx[i]);
150  }
151  }
152 
153  #endif
154 
155  /* -----------------------------------------
156  compute total source terms
157  ----------------------------------------- */
158 
159  for (i = beg; i <= end; i++) {
160 
161  src = state->src[i];
162  for (nv = NFLX; nv--; ) vc[nv] = state->vh[i][nv];
163 
164  vB = EXPAND(vc[VX1]*vc[BX1], + vc[VX2]*vc[BX2], + vc[VX3]*vc[BX3]);
165  src[RHO] = 0.0;
166  EXPAND(src[MX1] = 0.0; ,
167  src[MX2] = 0.0; ,
168  src[MX3] = 0.0;)
169 
170  EXPAND(src[BX1] = -vc[VX1]*divB[i]; ,
171  src[BX2] = -vc[VX2]*divB[i]; ,
172  src[BX3] = -vc[VX3]*divB[i];)
173 
174  #if EOS != ISOTHERMAL
175  src[ENG] = 0.0;
176  #endif
177 
178  }
179 }
180 #endif
#define MX3
Definition: mod_defs.h:22
#define EOS
Definition: pluto.h:341
#define MX1
Definition: mod_defs.h:20
double ** vh
Primitive state at n+1/2 (only for one step method)
Definition: structs.h:162
#define VX2
Definition: mod_defs.h:29
double ** flux
upwind flux computed with the Riemann solver
Definition: structs.h:149
void HLL_DIVB_SOURCE(const State_1D *state, double **Uhll, int beg, int end, Grid *grid)
Definition: source.c:53
#define RHO
Definition: mod_defs.h:19
void POWELL_DIVB_SOURCE(const State_1D *state, int is, int ie, Grid *grid)
Definition: source.c:5
double * dV
Cell volume.
Definition: structs.h:86
double ** vR
Primitive variables to the right of the interface, .
Definition: structs.h:139
double * dx
Definition: structs.h:83
#define VX1
Definition: mod_defs.h:28
#define KDIR
Definition: pluto.h:195
int BXn
Definition: globals.h:75
int g_i
x1 grid index when sweeping along the x2 or x3 direction.
Definition: globals.h:82
double ** src
Definition: structs.h:160
#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
int g_j
x2 grid index when sweeping along the x1 or x3 direction.
Definition: globals.h:83
Definition: structs.h:78
#define MX2
Definition: mod_defs.h:21
double *** divB
Definition: analysis.c:4
double * x
Definition: structs.h:80
#define ISOTHERMAL
Definition: pluto.h:49
#define s
#define BX3
Definition: mod_defs.h:27
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
int i
Definition: analysis.c:2
#define BX1
Definition: mod_defs.h:25
#define VX3
Definition: mod_defs.h:30
double ** vL
Primitive variables to the left of the interface, .
Definition: structs.h:136
#define BX2
Definition: mod_defs.h:26
#define JDIR
Definition: pluto.h:194
#define NVAR
Definition: pluto.h:609
double * A
Right interface area, A[i] = .
Definition: structs.h:87