PLUTO
init.c File Reference

Hydrostatic atmosphere in a smoothed gravitational potential. More...

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

Go to the source code of this file.

Functions

void Init (double *us, double x1, double x2, double x3)
 
void Analysis (const Data *d, Grid *grid)
 
void UserDefBoundary (const Data *d, RBox *box, int side, Grid *grid)
 

Variables

static real acf
 
static real bcf
 
static real ccf
 

Detailed Description

Hydrostatic atmosphere in a smoothed gravitational potential.

Set initial conditions for an hydrostatic atmosphere in cylindrical coordinates. Gravity is given as

\[ g = \left\{\begin{array}{ll} -1/R^2 & \quad{\rm for} \quad R > 1 \\ \noalign{\medskip} aR + bR^2 + cR^3 & \quad{\rm for} \quad R < 1 \end{array}\right.\,, \]

where $R = \sqrt{r^2 + z^2}$ is the spherical radius. The coefficients $a$, $b$ and $c$ are chosen to guarantee continuity of $g$, its first and second derivative (optionally). Density and pressure are tied by the isothermal condition $P = \rho/a$ so that the hydrostatic condition is

\[ \frac{1}{\rho} \frac{d\rho}{dr} = ag\,, \]

with the normalization $\rho = 1$ at $R = 1$.

The runtime parameters that are read from pluto.ini are

  • g_inputParam[ALPHA]: sets the value of $a$;

Configurations:

  • #1-4: 2D cylindrical
  • #5: 3D cartesian

References:

Author
A. Mignone (migno.nosp@m.ne@p.nosp@m.h.uni.nosp@m.to.i.nosp@m.t)
Date
July 09, 2014

Definition in file init.c.

Function Documentation

void Analysis ( const Data d,
Grid grid 
)

Perform runtime data analysis.

Parameters
[in]dthe PLUTO Data structure
[in]gridpointer to array of Grid structures

Definition at line 93 of file init.c.

98 {
99 
100 }
void Init ( double *  us,
double  x1,
double  x2,
double  x3 
)

The Init() function can be used to assign initial conditions as as a function of spatial position.

Parameters
[out]va pointer to a vector of primitive variables
[in]x1coordinate point in the 1st dimension
[in]x2coordinate point in the 2nd dimension
[in]x3coordinate point in the 3rdt dimension

The meaning of x1, x2 and x3 depends on the geometry:

\[ \begin{array}{cccl} x_1 & x_2 & x_3 & \mathrm{Geometry} \\ \noalign{\medskip} \hline x & y & z & \mathrm{Cartesian} \\ \noalign{\medskip} R & z & - & \mathrm{cylindrical} \\ \noalign{\medskip} R & \phi & z & \mathrm{polar} \\ \noalign{\medskip} r & \theta & \phi & \mathrm{spherical} \end{array} \]

Variable names are accessed by means of an index v[nv], where nv = RHO is density, nv = PRS is pressure, nv = (VX1, VX2, VX3) are the three components of velocity, and so forth.

Definition at line 43 of file init.c.

49 {
50  double rs, scrh;
51 
52 /* ------------------------------------------
53  with this choice g, g' and g''
54  will be continuous at R = 1
55  ----------------------------------------- */
56 
57  acf = -10.0;
58  bcf = 15.0;
59  ccf = -6.0;
60 
61 /* -----------------------------------
62  with this choice g and g'
63  will be continuous
64  ---------------------------------- */
65 
66  acf = -3.0;
67  bcf = 2.0;
68  ccf = 0.0;
69 
70  #if GEOMETRY == CARTESIAN
71  rs = sqrt(x1*x1 + x2*x2 + x3*x3);
72  #elif GEOMETRY == CYLINDRICAL
73  rs = sqrt(x1*x1 + x2*x2);
74  #elif GEOMETRY == SPHERICAL
75  rs = sqrt(x1*x1);
76  #endif
77 
78  if (rs > 1.0){
79  us[RHO] = exp(g_inputParam[ALPHA]*(1.0/rs - 1.0));
80  }else{
81  scrh = 0.5*acf*(rs*rs - 1.0) + 1.0/3.0*bcf*(rs*rs*rs - 1.0)
82  + 0.25*ccf*(rs*rs*rs*rs - 1.0);
83  us[RHO] = exp(scrh*g_inputParam[ALPHA]);
84  }
85 
86  us[VX1] = 0.0;
87  us[VX2] = 0.0;
88  us[VX3] = 0.0;
89  us[PRS] = us[RHO]/g_inputParam[ALPHA];
90  us[TRC] = 0.0;
91 }
#define VX2
Definition: mod_defs.h:29
#define RHO
Definition: mod_defs.h:19
tuple scrh
Definition: configure.py:200
#define TRC
Definition: pluto.h:581
#define VX1
Definition: mod_defs.h:28
static real ccf
Definition: init.c:40
double g_inputParam[32]
Array containing the user-defined parameters.
Definition: globals.h:131
static real bcf
Definition: init.c:40
#define ALPHA
static real acf
Definition: init.c:40
#define VX3
Definition: mod_defs.h:30
void UserDefBoundary ( const Data d,
RBox box,
int  side,
Grid grid 
)

Assign user-defined boundary conditions.

Parameters
[in,out]dpointer to the PLUTO data structure containing cell-centered primitive quantities (d->Vc) and staggered magnetic fields (d->Vs, when used) to be filled.
[in]boxpointer to a RBox structure containing the lower and upper indices of the ghost zone-centers/nodes or edges at which data values should be assigned.
[in]sidespecifies the boundary side where ghost zones need to be filled. It can assume the following pre-definite values: X1_BEG, X1_END, X2_BEG, X2_END, X3_BEG, X3_END. The special value side == 0 is used to control a region inside the computational domain.
[in]gridpointer to an array of Grid structures.

Assign user-defined boundary conditions in the lower boundary ghost zones. The profile is top-hat:

\[ V_{ij} = \left\{\begin{array}{ll} V_{\rm jet} & \quad\mathrm{for}\quad r_i < 1 \\ \noalign{\medskip} \mathrm{Reflect}(V) & \quad\mathrm{otherwise} \end{array}\right. \]

where $ V_{\rm jet} = (\rho,v,p)_{\rm jet} = (1,M,1/\Gamma)$ and M is the flow Mach number (the unit velocity is the jet sound speed, so $ v = M$).

Assign user-defined boundary conditions:

  • left side (x beg): constant shocked values
  • bottom side (y beg): constant shocked values for x < 1/6 and reflective boundary otherwise.
  • top side (y end): time-dependent boundary: for $ x < x_s(t) = 10 t/\sin\alpha + 1/6 + 1.0/\tan\alpha $ we use fixed (post-shock) values. Unperturbed values otherwise.

Definition at line 102 of file init.c.

106 {
107  int i, j, k;
108  double *x1, *x2, *x3;
109  double rs;
110 
111  x1 = grid[IDIR].xgc;
112  x2 = grid[JDIR].xgc;
113  x3 = grid[KDIR].xgc;
114 
115  if (side == X1_END) {
116 
117  X1_END_LOOP(k,j,i){
118  #if GEOMETRY == CARTESIAN
119  rs = sqrt(x1[i]*x1[i] + x2[j]*x2[j] + x3[k]*x3[k]);
120  #elif GEOMETRY == CYLINDRICAL
121  rs = sqrt(x1[i]*x1[i] + x2[j]*x2[j]);
122  #elif GEOMETRY == SPHERICAL
123  rs = sqrt(x1[i]*x1[i]);
124  #endif
125  d->Vc[RHO][k][j][i] = exp(g_inputParam[ALPHA]*(1.0/rs - 1.0));
126  EXPAND(d->Vc[VX1][k][j][i] = 0.0; ,
127  d->Vc[VX2][k][j][i] = 0.0; ,
128  d->Vc[VX3][k][j][i] = 0.0;)
129  d->Vc[PRS][k][j][i] = exp(g_inputParam[ALPHA]*(1.0/rs - 1.0))/g_inputParam[ALPHA];
130  }
131 
132  } else if (side == X2_END) {
133 
134  X2_END_LOOP(k,j,i){
135  #if GEOMETRY == CARTESIAN
136  rs = sqrt(x1[i]*x1[i] + x2[j]*x2[j] + x3[k]*x3[k]);
137  #elif GEOMETRY == CYLINDRICAL
138  rs = sqrt(x1[i]*x1[i] + x2[j]*x2[j]);
139  #elif GEOMETRY == SPHERICAL
140  rs = sqrt(x1[i]*x1[i]);
141  #endif
142  d->Vc[RHO][k][j][i] = exp(g_inputParam[ALPHA]*(1.0/rs - 1.0));
143  EXPAND(d->Vc[VX1][k][j][i] = 0.0; ,
144  d->Vc[VX2][k][j][i] = 0.0; ,
145  d->Vc[VX3][k][j][i] = 0.0;)
146  d->Vc[PRS][k][j][i] = exp(g_inputParam[ALPHA]*(1.0/rs - 1.0))/g_inputParam[ALPHA];
147  }
148 
149  } else if (side == X3_END) { /* Only Cartesian */
150 
151  X3_END_LOOP(k,j,i){
152  rs = sqrt(x1[i]*x1[i] + x2[j]*x2[j] + x3[k]*x3[k]);
153  d->Vc[RHO][k][j][i] = exp(g_inputParam[ALPHA]*(1.0/rs - 1.0));
154  EXPAND(d->Vc[VX1][k][j][i] = 0.0; ,
155  d->Vc[VX2][k][j][i] = 0.0; ,
156  d->Vc[VX3][k][j][i] = 0.0;)
157  d->Vc[PRS][k][j][i] = exp(g_inputParam[ALPHA]*(1.0/rs - 1.0))/g_inputParam[ALPHA];
158  }
159  }
160 }
#define VX2
Definition: mod_defs.h:29
#define RHO
Definition: mod_defs.h:19
double **** Vc
The main four-index data array used for cell-centered primitive variables.
Definition: structs.h:31
#define X3_END_LOOP(k, j, i)
Definition: macros.h:52
#define X1_END
Boundary region at X1 end.
Definition: pluto.h:147
#define VX1
Definition: mod_defs.h:28
#define KDIR
Definition: pluto.h:195
#define X2_END_LOOP(k, j, i)
Definition: macros.h:51
#define IDIR
Definition: pluto.h:193
#define X2_END
Boundary region at X2 end.
Definition: pluto.h:149
double g_inputParam[32]
Array containing the user-defined parameters.
Definition: globals.h:131
int j
Definition: analysis.c:2
if(divB==NULL)
Definition: analysis.c:10
int k
Definition: analysis.c:2
#define X1_END_LOOP(k, j, i)
Definition: macros.h:50
double * xgc
Cell volumetric centroid (!= x when geometry != CARTESIAN).
Definition: structs.h:84
#define ALPHA
#define X3_END
Boundary region at X3 end.
Definition: pluto.h:151
int i
Definition: analysis.c:2
#define VX3
Definition: mod_defs.h:30
#define JDIR
Definition: pluto.h:194

Variable Documentation

real acf
static

Definition at line 40 of file init.c.

real bcf
static

Definition at line 40 of file init.c.

real ccf
static

Definition at line 40 of file init.c.