PLUTO  4.0
 All Data Structures Files Functions Variables Enumerations Macros Pages
Data Structures | Macros | Functions
mod_defs.h File Reference

Module header file for relativistic MHD (RMHD). More...

Go to the source code of this file.

Data Structures

struct  Map_param
 

Macros

#define SUBTRACT_DENSITY   YES
 

Functions

int ConsToPrim (double **, double **, int, int, unsigned char *)
 
void Enthalpy (double **, double *, int, int)
 
void Entropy (double **, double *, int, int)
 
int EntropySolve (Map_param *)
 
int EnergySolve (Map_param *)
 
int PressureFix (Map_param *)
 
void Flux (double **, double **, double *, double **, double *, int, int)
 
int MaxSignalSpeed (double **, double *, double *, double *, double *, int, int)
 
void PrimToCons (double **, double **, int, int)
 
int Magnetosonic (double *vp, double cs2, double h, double *lambda)
 
int QuarticSolve (double, double, double, double, double *)
 

Detailed Description

Set label, indexes and basic prototyping for the relativistic MHD module.

Authors
A. Mignone (migno.nosp@m.ne@p.nosp@m.h.uni.nosp@m.to.i.nosp@m.t)
C. Zanni (zanni.nosp@m.@oat.nosp@m.o.ina.nosp@m.f.it)
Date
Oct 5, 2012

Macro Definition Documentation

#define SUBTRACT_DENSITY   YES

By turning SUBTRACT_DENSITY to YES, we let PLUTO evolve the total energy minus the mass density contribution.

Function Documentation

int ConsToPrim ( double **  ucons,
double **  uprim,
int  ibeg,
int  iend,
unsigned char *  flag 
)

Convert from conservative to primitive variables.

Parameters
[in]uconsarray of conservative variables
[out]uprimarray of primitive variables
[in]begstarting index of computation
[in]endfinal index of computation
[out]flagarray of flags tagging zones where conversion went wrong.
Returns
Return (0) if conversion was succesful in every zone [ibeg,iend]. Otherwise, return a non-zero integer number giving the bit flag(s) turned on during the conversion process. In this case, flag contains the failure codes of those zones where where conversion did not go through.

Convert from conservative to primitive variables.

Parameters
[in]uconsarray of conservative variables
[out]uprimarray of primitive variables
[in]begstarting index of computation
[in]endfinal index of computation
[out]flagarray of flags tagging zones where conversion went wrong.
Returns
This function returns 0 on success. Otherwise, a non-zero integer number giving the bit flag turned on during the conversion.

Convert from conservative to primitive variables.

Parameters
[in]uconsarray of conservative variables
[out]uprimarray of primitive variables
[in]begstarting index of computation
[in]endfinal index of computation
[out]flagarray of flags tagging zones where conversion went wrong.
Returns
Return (0) if conversion was succesful in every zone [ibeg,iend]. Otherwise, return a non-zero integer number giving the bit flag(s) turned on during the conversion process. In this case, flag contains the failure codes of those zones where where conversion did not go through.

Convert from conservative to primitive variables.

Parameters
[in]uconsarray of conservative variables
[out]uprimarray of primitive variables
[in]begstarting index of computation
[in]endfinal index of computation
[out]flagarray of flags tagging zones where conversion went wrong.
Returns
This function returns 0 on success. Otherwise, a non-zero integer number giving the bit flag turned on during the conversion.
int EnergySolve ( Map_param par)

Use Newton algorithm to recover pressure p0 from conservative quantities D, m and E

Solve f(W) = 0, where f(W) is Equation (A4) or (A6).

Returns
Return (0) is successful, (1) otherwise.
void Enthalpy ( double **  v,
double *  h,
int  beg,
int  end 
)

Compute the enthalpy.

Parameters
[in]v1D array of primitive quantities
[in]h1D array of enthalpy values
[in]beginitial index of computation
[in]endfinal index of computation
Returns
This function has no return value.
void Entropy ( double **  v,
double *  s,
int  beg,
int  end 
)

Compute the entropy.

Parameters
[in]v1D array of primitive quantities
[in]s1D array of entropy values
[in]isinitial index of computation
[in]iefinal index of computation
Returns
This function has no return value.

Compute the entropy.

Parameters
[in]v1D array of primitive quantities
[in]s1D array of entropy values
[in]beginitial index of computation
[in]endfinal index of computation
Returns
This function has no return value.
int EntropySolve ( Map_param par)

Convert the conservative variables {D, m, sigma_c} (where sigma_c = D*sigma is the conserved entropy) to primitive variable using a Newton-Raphson/Bisection scheme.

LAST MODIFIED

February 5th (2012) by C. Zanni & A. Mignone (zanni.nosp@m.@oat.nosp@m.o.ina.nosp@m.f.it, migno.nosp@m.ne@o.nosp@m.ato.i.nosp@m.naf..nosp@m.it)

void Flux ( double **  ucons,
double **  uprim,
double *  h,
double **  fx,
real *  pr,
int  beg,
int  end 
)
Parameters
[in]u1D array of conserved quantities
[in]w1D array of primitive quantities
[in]a21D array of sound speeds
[out]fx1D array of fluxes (total pressure excluded)
[out]p1D array of pressure values
[in]beginitial index of computation
[in]endfinal index of computation
Returns
This function has no return value.
int Magnetosonic ( double *  vp,
double  cs2,
double  h,
double *  lambda 
)

Find the four magneto-sonic speeds (fast and slow) for the relativistic MHD equations (RMHD). We follow the notations introduced in Eq. (16) in

Del Zanna, Bucciantini and Londrillo, A&A, 400, 397–413, 2003

and also Eq. (57-58) of Mignone & Bodo, MNRAS, 2006 for the degenerate cases.

The quartic equation is solved analytically and the solver (function 'QuarticSolve') assumes all roots are real (which should be always the case here, since eigenvalues are recovered from primitive variables). The coefficients of the quartic were found through the following MAPLE script:

 *   ------------------------------------------------
restart;

u[0] := gamma;
u[x] := gamma*v[x];
b[0] := gamma*vB;
b[x] := B[x]/gamma + b[0]*v[x];
wt   := w + b^2;

#############################################################
"fdZ will be equation (16) of Del Zanna 2003 times w_{tot}";

e2  := c[s]^2 + b^2*(1 - c[s]^2)/wt;
fdZ := (1-e2)*(u[0]*lambda - u[x])^4 + (1-lambda^2)*
       (c[s]^2*(b[0]*lambda - b[x])^2/wt - e2*(u[0]*lambda - u[x])^2);

########################################################
"print the coefficients of the quartic polynomial fdZ";

coeff(fMB,lambda,4);
coeff(fMB,lambda,3);
coeff(fMB,lambda,2);
coeff(fMB,lambda,1);
coeff(fMB,lambda,0);

fdZ := fdZ*wt;  

######################################################
"fMB will be equation (56) of Mignone & Bodo (2006)";

a  := gamma*(lambda-v[x]);
BB := b[x] - lambda*b[0];
fMB := w*(1-c[s]^2)*a^4 - (1-lambda^2)*((b^2+w*c[s]^2)*a^2-c[s]^2*BB^2);

######################################
"check that fdZ = fMB";

simplify(fdZ - fMB);

########################################################
"print the coefficients of the quartic polynomial fMB";

coeff(fMB,lambda,4);
coeff(fMB,lambda,3);
coeff(fMB,lambda,2);
coeff(fMB,lambda,1);
coeff(fMB,lambda,0);


###############################################
"Degenerate case 2: Bx = 0, Eq (58) in MB06";

B[x] := 0;

a2 := w*(c[s]^2 + gamma^2*(1-c[s]^2)) + Q; 
a1 := -2*w*gamma^2*v[x]*(1-c[s]^2);
a0 := w*(-c[s]^2 + gamma^2*v[x]^2*(1-c[s]^2))-Q;

"delc is a good way to evaluate the determinant so that it is > 0";
del  := a1^2 - 4*a2*a0;
delc := 4*(w*c[s]^2 + Q)*(w*(1-c[s]^2)*gamma^2*(1-v[x]^2) + (w*c[s]^2+Q));
simplify(del-delc);

############################################################################
"check the correctness of the coefficients of Eq. (58) in MB06";

Q  := b^2 - c[s]^2*vB^2;
divide(fMB, (lambda-v[x])^2,'fMB2');
simplify(a2 - coeff(fMB2,lambda,2)/gamma^2);
simplify(a1 - coeff(fMB2,lambda,1)/gamma^2);
simplify(a0 - coeff(fMB2,lambda,0)/gamma^2);
int MaxSignalSpeed ( double **  v,
double *  a2,
double *  h,
double *  cmin,
double *  cmax,
int  beg,
int  end 
)

Return the rightmost (cmax) and leftmost (cmin) wave propagation speed in the Riemann fan

int PressureFix ( Map_param par)

Fix p to a small value, solve for the square of velocity by using secant algorithm applied to Eq (9) of Mignone, Plewa & Bodo (2005). This step involved re-computing W at each step of the iteration. Once the root has been found, we recompute total energy E. Return 0 if succesful, 1 otherwise.

void PrimToCons ( double **  uprim,
double **  ucons,
int  ibeg,
int  iend 
)

Convert primitive variables to conservative variables.

Parameters
[in]uprimarray of primitive variables
[out]uconsarray of conservative variables
[in]begstarting index of computation
[in]endfinal index of computation

Convert primitive variables in conservative variables.

Parameters
[in]uprimarray of primitive variables
[out]uconsarray of conservative variables
[in]begstarting index of computation
[in]endfinal index of computation

Convert primitive variables to conservative variables.

Parameters
[in]uprimarray of primitive variables
[out]uconsarray of conservative variables
[in]begstarting index of computation
[in]endfinal index of computation

Convert primitive variables in conservative variables.

Parameters
[in]uprimarray of primitive variables
[out]uconsarray of conservative variables
[in]begstarting index of computation
[in]endfinal index of computation
int QuarticSolve ( double  b,
double  c,
double  d,
double  e,
double *  z 
)

Solve a quartic equation in the form

\[ z^4 + bz^3 + cz^2 + dz + e = 0 \]

For its purpose, it is assumed that ALL roots are double. This makes things faster.

Parameters
[in]bcoefficient of the quartic
[in]ccoefficient of the quartic
[in]dcoefficient of the quartic
[in]ecoefficient of the quartic
[out]zvector containing the (double) roots of the quartic

Reference:

http://www.1728.com/quartic2.htm

Returns
Return 0 on success.