LORENE
Lorene::Eos_bf_poly Class Reference

Analytic equation of state for two fluids (relativistic case). More...

#include <eos_bifluid.h>

Inheritance diagram for Lorene::Eos_bf_poly:
Lorene::Eos_bifluid Lorene::Eos_bf_poly_newt

Public Member Functions

 Eos_bf_poly (double kappa1, double kappa2, double kappa3, double beta)
 Standard constructor. More...
 
 Eos_bf_poly (double gamma1, double gamma2, double gamma3, double gamma4, double gamma5, double gamma6, double kappa1, double kappa2, double kappa3, double beta, double mass1=1, double mass2=1, double relax=0.5, double precis=1.e-9, double ecart=1.e-8)
 Standard constructor with all parameters. More...
 
 Eos_bf_poly (const Eos_bf_poly &)
 Copy constructor. More...
 
virtual ~Eos_bf_poly ()
 Destructor. More...
 
void operator= (const Eos_bf_poly &)
 Assignment to another Eos_bf_poly. More...
 
virtual bool operator== (const Eos_bifluid &) const
 Comparison operator (egality) More...
 
virtual bool operator!= (const Eos_bifluid &) const
 Comparison operator (difference) More...
 
virtual int identify () const
 Returns a number to identify the sub-classe of Eos_bifluid
the object belongs to. More...
 
double get_gam1 () const
 Returns the adiabatic index $\gamma_1$. More...
 
double get_gam2 () const
 Returns the adiabatic index $\gamma_2$. More...
 
double get_gam3 () const
 Returns the adiabatic index $\gamma_3$. More...
 
double get_gam4 () const
 Returns the adiabatic index $\gamma_4$. More...
 
double get_gam5 () const
 Returns the adiabatic index $\gamma_5$. More...
 
double get_gam6 () const
 Returns the adiabatic index $\gamma_6$. More...
 
double get_kap1 () const
 Returns the pressure coefficient $\kappa_1$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$. More...
 
double get_kap2 () const
 Returns the pressure coefficient $\kappa_2$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$. More...
 
double get_kap3 () const
 Returns the pressure coefficient $\kappa_3$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$. More...
 
double get_beta () const
 Returns the coefficient $\beta$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$. More...
 
int get_typeos () const
 
virtual void sauve (FILE *) const
 Save in a file. More...
 
virtual bool nbar_ent_p (const double ent1, const double ent2, const double delta2, double &nbar1, double &nbar2) const
 Computes both baryon densities from the log-enthalpies. More...
 
virtual double nbar_ent_p1 (const double ent1) const
 Computes baryon density out of the log-enthalpy asuming that only fluid 1 is present. More...
 
virtual double nbar_ent_p2 (const double ent2) const
 Computes baryon density out of the log-enthalpy assuming that only fluid 2 is present. More...
 
virtual double ener_nbar_p (const double nbar1, const double nbar2, const double delta2) const
 Computes the total energy density from the baryonic densities and the relative velocity. More...
 
virtual double press_nbar_p (const double nbar1, const double nbar2, const double delta2) const
 Computes the pressure from the baryonic densities and the relative velocity. More...
 
virtual Eostrans2Eos () const
 Makes a translation from Eos_bifluid to Eos . More...
 
virtual double get_K11 (const double n1, const double n2, const double delta2) const
 Computes the derivative of the energy with respect to (baryonic density 1) $^2$. More...
 
virtual double get_K12 (const double n1, const double n2, const double delta2) const
 Computes the derivative of the energy with respect to $x^2=n_1n_2\Gamma_\Delta$. More...
 
virtual double get_K22 (const double n1, const double n2, const double delta2) const
 Computes the derivative of the energy/(baryonic density 2) $^2$. More...
 
string get_name () const
 Returns the EOS name. More...
 
double get_m1 () const
 Return the individual particule mass $m_1$. More...
 
double get_m2 () const
 Return the individual particule mass $m_2$. More...
 
virtual void calcule_tout (const Cmp &ent1, const Cmp &ent2, const Cmp &delta2, Cmp &nbar1, Cmp &nbar2, Cmp &ener, Cmp &press, int nzet, int l_min=0) const
 General computational method for Cmp 's, it computes both baryon densities, energy and pressure profiles. More...
 
void nbar_ent (const Cmp &ent1, const Cmp &ent2, const Cmp &delta2, Cmp &nbar1, Cmp &nbar2, int nzet, int l_min=0) const
 Computes both baryon density fields from the log-enthalpy fields and the relative velocity. More...
 
Cmp ener_ent (const Cmp &ent1, const Cmp &ent2, const Cmp &delta2, int nzet, int l_min=0) const
 Computes the total energy density from the log-enthalpy fields and the relative velocity. More...
 
Cmp press_ent (const Cmp &ent1, const Cmp &ent2, const Cmp &delta2, int nzet, int l_min=0) const
 Computes the pressure from the log-enthalpy fields and the relative velocity. More...
 
Cmp get_Knn (const Cmp &nbar1, const Cmp &nbar2, const Cmp &x2, int nzet, int l_min=0) const
 Computes the derivatives of the energy/(baryonic density 1) $^2$. More...
 
Cmp get_Kpp (const Cmp &nbar1, const Cmp &nbar2, const Cmp &x2, int nzet, int l_min=0) const
 Computes the derivatives of the energy/(baryonic density 2) $^2$. More...
 
Cmp get_Knp (const Cmp &nbar1, const Cmp &nbar2, const Cmp &x2, int nzet, int l_min=0) const
 Computes the derivatives of the energy with respect to $x^2=n_1n_2\Gamma_\Delta^2$. More...
 
void calcule (const Cmp &nbar1, const Cmp &nbar2, const Cmp &x2, int nzet, int l_min, double(Eos_bifluid::*fait)(double, double, double) const, Cmp &resu) const
 General computational method for Cmp 's ( $K^{ij}$'s). More...
 

Static Public Member Functions

static Eos_bifluideos_from_file (FILE *)
 Construction of an EOS from a binary file. More...
 
static Eos_bifluideos_from_file (const char *fname)
 Construction of an EOS from a formatted file. More...
 
static Eos_bifluideos_from_file (ifstream &)
 Construction of an EOS from a formatted file. More...
 

Protected Member Functions

 Eos_bf_poly (FILE *)
 Constructor from a binary file (created by the function sauve(FILE*) ). More...
 
 Eos_bf_poly (const char *fname)
 Constructor from a formatted file. More...
 
void set_auxiliary ()
 Computes the auxiliary quantities gam1m1 , gam2m1 and gam3m1. More...
 
void determine_type ()
 Determines the type of the analytical EOS (see typeos ) More...
 
virtual ostream & operator>> (ostream &) const
 Operator >> More...
 

Protected Attributes

double gam1
 Adiabatic indexes $\gamma_1$, see Eq.~eeosbfpolye}. More...
 
double gam2
 Adiabatic indexes $\gamma_2$, see Eq.~eeosbfpolye}. More...
 
double gam3
 Adiabatic indexes $\gamma_3$, see Eq.~eeosbfpolye}. More...
 
double gam4
 Adiabatic indexes $\gamma_4$, see Eq.~eeosbfpolye}. More...
 
double gam5
 Adiabatic indexes $\gamma_5$, see Eq.~eeosbfpolye}. More...
 
double gam6
 Adiabatic indexes $\gamma_6$, see Eq.~eeosbfpolye}. More...
 
double kap1
 Pressure coefficient $\kappa_1$ , see Eq. More...
 
double kap2
 Pressure coefficient $\kappa_2$ , see Eq. More...
 
double kap3
 Pressure coefficient $\kappa_3$ , see Eq. More...
 
double beta
 Coefficient $\beta$ , see Eq. More...
 
double gam1m1
 $\gamma_1-1$ More...
 
double gam2m1
 $\gamma_2-1$ More...
 
double gam34m1
 $\gamma_3+\gamma_4-1$ More...
 
double gam56m1
 $\gamma_5+\gamma_6-1$ More...
 
int typeos
 The bi-fluid analytical EOS type: More...
 
double relax
 Parameters needed for some inversions of the EOS. More...
 
double precis
 contains the precision required in zerosec_b More...
 
double ecart
 contains the precision required in the relaxation nbar_ent_p More...
 
string name
 EOS name. More...
 
double m_1
 Individual particle mass $m_1$
[unit: $m_B = 1.66\ 10^{-27} \ {\rm kg}$]. More...
 
double m_2
 Individual particle mass $m_2$
[unit: $m_B = 1.66\ 10^{-27} \ {\rm kg}$]. More...
 

Friends

Eos_bifluidEos_bifluid::eos_from_file (FILE *)
 The construction functions from a file. More...
 
Eos_bifluidEos_bifluid::eos_from_file (const char *fname)
 

Detailed Description

Analytic equation of state for two fluids (relativistic case).

This equation of state (EOS) corresponds to two types of relativistic particles of rest mass is $m_1$ and $m_2$, whose total energy density $\cal{E}$ is related to their numerical densities $n_1$, $n_2$ and relative velocity

\[ \Gamma_\Delta = -g_{\alpha\beta} u^\alpha_{\rm n} u^\beta_{\rm p} \label{e:defgamamdelta} \]

( $u^\alpha_{\rm n}$ and $u^\alpha_{\rm p}$ being the 4-velocities of both fluids), by

\[ \label{eeosbfpolye} {\cal{E}} = \frac{1}{2}\kappa_1 n_1^{\gamma_1} + m_1 c^2 \, n_1 \ + \frac{1}{2}\kappa_2 n_2^{\gamma_2} + m_2 c^2 \, n_2 \ + \kappa_3 n_1^{\gamma_3} n_2^{\gamma_4} \ + \Delta^2 \beta n_1^{\gamma_5} n_2^{\gamma_6}\ . \]

The relativistic (i.e. including rest mass energy) chemical potentials are then

\[ \mu_1 := {{\rm d}{\cal{E}} \over {\rm d}n_1} = \frac{1}{2}\gamma_1\kappa_1 n_1^{\gamma_1-1} + m_1 c^2 + \gamma_3 \kappa_3 n_1^{\gamma_3-1} n_2^{\gamma_4} + \Delta^2 \gamma_5 \beta n_1^{\gamma_5-1} n_2^{\gamma_6}\ , \]

\[ \mu_2 := {{\rm d}{\cal{E}} \over {\rm d}n_2} = \frac{1}{2}\gamma_2\kappa_2 n_2^{\gamma_2-1} + m_2 c^2 + \gamma_4 \kappa_3 n_1^{\gamma_3} n_2^{\gamma_4-1} + \Delta^2 \gamma_6 \beta n_1^{\gamma_5} n_2^{\gamma_6-1} \ . \]

The pressure is given by the (zero-temperature) First Law of Thermodynamics: $p = \mu_1 n_1 + \mu_2 n_2 - {\cal E}$, so that

\[ p = \frac{1}{2} (\gamma_1 -1)\kappa_1 n_1^{\gamma_1} + \frac{1}{2}(\gamma_2-1)\kappa_2 n_2^{\gamma_2} + (\gamma_3 +\gamma_4 -1)\kappa_3 n_1^{\gamma_3}n_2^{\gamma_4} + \Delta^2 \beta \left( (\gamma_5 + \gamma_6 - 1) n_1^{\gamma_5} n_2^{\gamma_6} \right) \ . \]

The log-enthalpies are defined as the logarithm of the ratio of the enthalpy per particle (see Eq.~eeosbfdefent}) by the particle rest mass energy :

\[\label{eeosbfentanal} H_i := c^2 \ln \left( \frac{f_i}{m_ic^2} \right) \ . \]

From this system, the particle densities are obtained in term of the log-enthalpies. (The system (eeosbfentanal}) is a linear one if $\gamma_1 = \gamma_2 = 2$ and $\gamma_3 = \gamma_4 = \gamma_5 = \gamma_6 = 1$).

The energy density $\cal E$and pressure p can then be obtained as functions of baryonic densities.()

Definition at line 717 of file eos_bifluid.h.

Constructor & Destructor Documentation

◆ Eos_bf_poly() [1/5]

Lorene::Eos_bf_poly::Eos_bf_poly ( double  kappa1,
double  kappa2,
double  kappa3,
double  beta 
)

Standard constructor.

The adiabatic indexes $\gamma_1$ and $\gamma_2$ are set to 2. All other adiabatic indexes $\gamma_i$, $i\in 3\dots 6$ are set to 1. The individual particle masses $m_1$ and $m_2$ are set to the mean baryon mass $m_B = 1.66\ 10^{-27} \ {\rm kg}$. The inversion parameters are set to their default values (see hereafter the consrtuctor with all parameters).

Parameters
kappa1pressure coefficient $\kappa_1$
kappa2pressure coefficient $\kappa_2$
kappa3pressure coefficient $\kappa_3$
betacoefficient in the entrainment term $\beta$
(cf. Eq.~(eeosbfpolye})) [unit: $\rho_{\rm nuc} c^2$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$

Definition at line 160 of file eos_bf_poly.C.

References determine_type(), and set_auxiliary().

◆ Eos_bf_poly() [2/5]

Lorene::Eos_bf_poly::Eos_bf_poly ( double  gamma1,
double  gamma2,
double  gamma3,
double  gamma4,
double  gamma5,
double  gamma6,
double  kappa1,
double  kappa2,
double  kappa3,
double  beta,
double  mass1 = 1,
double  mass2 = 1,
double  relax = 0.5,
double  precis = 1.e-9,
double  ecart = 1.e-8 
)

Standard constructor with all parameters.

Parameters
gamma1adiabatic index $\gamma_1$
gamma2adiabatic index $\gamma_2$
gamma3adiabatic index $\gamma_3$
gamma4adiabatic index $\gamma_4$
gamma5adiabatic index $\gamma_5$
gamma6adiabatic index $\gamma_6$ (cf. Eq.~(eeosbfpolye}))
kappa1pressure coefficient $\kappa_1$
kappa2pressure coefficient $\kappa_2$
kappa3pressure coefficient $\kappa_3$
betacoefficient in the entrainment term $\beta$
(cf. Eq.~(eeosbfpolye})) [unit: $\rho_{\rm nuc} c^2$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$
mass1individual particule mass $m_1$ (neutrons)
mass2individual particule mass $m_2$ (protons)
[unit: $m_B = 1.66\ 10^{-27} \ {\rm kg}$]
relaxrelaxation parameter (see par_inv)
precisprecision parameter for zerosec_b (see par_inv)
relaxprecision parameter for relaxation procedure (see par_inv)

Definition at line 174 of file eos_bf_poly.C.

References determine_type(), and set_auxiliary().

◆ Eos_bf_poly() [3/5]

Lorene::Eos_bf_poly::Eos_bf_poly ( const Eos_bf_poly eosi)

Copy constructor.

Definition at line 192 of file eos_bf_poly.C.

References set_auxiliary().

◆ Eos_bf_poly() [4/5]

Lorene::Eos_bf_poly::Eos_bf_poly ( FILE *  fich)
protected

Constructor from a binary file (created by the function sauve(FILE*) ).

This constructor is protected because any EOS construction from a binary file must be done via the function Eos_bifluid::eos_from_file(FILE*) .

Definition at line 208 of file eos_bf_poly.C.

References beta, determine_type(), ecart, Lorene::fread_be(), gam1, gam2, gam3, gam4, gam5, gam6, kap1, kap2, kap3, precis, relax, and set_auxiliary().

◆ Eos_bf_poly() [5/5]

Lorene::Eos_bf_poly::Eos_bf_poly ( const char *  fname)
protected

Constructor from a formatted file.

This constructor is protected because any EOS construction from a formatted file must be done via the function Eos_bifluid::eos_from_file(const char*) .

Definition at line 234 of file eos_bf_poly.C.

References beta, determine_type(), ecart, gam1, gam2, gam3, gam4, gam5, gam6, kap1, kap2, kap3, precis, Lorene::read_variable(), relax, set_auxiliary(), and typeos.

◆ ~Eos_bf_poly()

Lorene::Eos_bf_poly::~Eos_bf_poly ( )
virtual

Destructor.

Definition at line 288 of file eos_bf_poly.C.

Member Function Documentation

◆ calcule()

void Lorene::Eos_bifluid::calcule ( const Cmp nbar1,
const Cmp nbar2,
const Cmp x2,
int  nzet,
int  l_min,
double(Eos_bifluid::*)(double, double, double) const  fait,
Cmp resu 
) const
inherited

General computational method for Cmp 's ( $K^{ij}$'s).

Parameters
nbar1[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 1 at which the derivatives are to be computed.
nbar2[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 2 at which the derivatives are to be computed
x2[input, unit $n_{\rm nuc}^2c^2$] relative velocity $\times$both densities at which the derivative is to be computed
nzet[input] number of domains where resu is to be computed.
l_min[input] index of the innermost domain is which resu is to be computed [default value: 0]; resu is computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.
fait[input] pointer on the member function of class Eos_bifluid which performs the pointwise calculation.
resu[output] result of the computation.

Definition at line 661 of file eos_bifluid.C.

References Lorene::Cmp::get_etat().

◆ calcule_tout()

void Lorene::Eos_bifluid::calcule_tout ( const Cmp ent1,
const Cmp ent2,
const Cmp delta2,
Cmp nbar1,
Cmp nbar2,
Cmp ener,
Cmp press,
int  nzet,
int  l_min = 0 
) const
virtualinherited

General computational method for Cmp 's, it computes both baryon densities, energy and pressure profiles.

Parameters
ent1[input] the first log-enthalpy field $H_1$.
ent2[input] the second log-enthalpy field $H_2$.
delta2[input] the relative velocity field $\Delta^2 $
nbar1[output] baryonic density of the first fluid
nbar2[output] baryonic density of the second fluid [unit: $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$]
ener[output] total energy density $\cal E$ of both fluids together
press[output] pressure p of both fluids together
nzet[input] number of domains where resu is to be computed.
l_min[input] index of the innermost domain is which resu is to be computed [default value: 0]; resu is computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.

Definition at line 286 of file eos_bifluid.C.

References Lorene::Cmp::get_etat().

◆ determine_type()

void Lorene::Eos_bf_poly::determine_type ( )
protected

Determines the type of the analytical EOS (see typeos )

Definition at line 337 of file eos_bf_poly.C.

References gam1, gam2, gam3, gam4, gam5, gam6, and typeos.

◆ ener_ent()

Cmp Lorene::Eos_bifluid::ener_ent ( const Cmp ent1,
const Cmp ent2,
const Cmp delta2,
int  nzet,
int  l_min = 0 
) const
inherited

Computes the total energy density from the log-enthalpy fields and the relative velocity.

Parameters
ent1[input, unit: $c^2$] log-enthalpy $H_1$
ent2[input, unit: $c^2$] log-enthalpy $H_2$
delta2[input, unit: $c^2$] relative velocity $\Delta^2$
nzetnumber of domains where the energy density is to be computed.
l_minindex of the innermost domain is which the energy density is to be computed [default value: 0]; the energy density is computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.
Returns
energy density field [unit: $\rho_{\rm nuc} c^2$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$

Definition at line 504 of file eos_bifluid.C.

References Lorene::Cmp::get_etat(), and Lorene::Cmp::get_mp().

◆ ener_nbar_p()

double Lorene::Eos_bf_poly::ener_nbar_p ( const double  nbar1,
const double  nbar2,
const double  delta2 
) const
virtual

Computes the total energy density from the baryonic densities and the relative velocity.

Parameters
nbar1[input] baryonic density of the first fluid
nbar2[input] baryonic density of the second fluid [unit: $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$]
delta2[input, unit: $c^2$] relative velocity $\Delta^2$
Returns
energy density $\cal E$ [unit: $\rho_{\rm nuc} c^2$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 892 of file eos_bf_poly.C.

References beta, gam1, gam2, gam3, gam4, gam5, gam6, kap1, kap2, kap3, Lorene::Eos_bifluid::m_1, Lorene::Eos_bifluid::m_2, and Lorene::pow().

◆ eos_from_file() [1/3]

Eos_bifluid * Lorene::Eos_bifluid::eos_from_file ( FILE *  fich)
staticinherited

Construction of an EOS from a binary file.

The file must have been created by the function sauve(FILE*) .

Definition at line 111 of file eos_bf_file.C.

References Lorene::fread_be().

◆ eos_from_file() [2/3]

Eos_bifluid * Lorene::Eos_bifluid::eos_from_file ( const char *  fname)
staticinherited

Construction of an EOS from a formatted file.

The following field has to be present:\ ident: [int] identifying the type of 2-fluid EOS 1 = relativistic polytropic EOS (class Eos_bf_poly ). \ 2 = Newtonian polytropic EOS (class Eos_bf_poly_newt ).

Definition at line 148 of file eos_bf_file.C.

References Lorene::read_variable().

◆ eos_from_file() [3/3]

Eos_bifluid * Lorene::Eos_bifluid::eos_from_file ( ifstream &  fich)
staticinherited

Construction of an EOS from a formatted file.

The fist line of the file must start by the EOS number, according to the following conventions:

  • 1 = 2-fluid relativistic polytropic EOS (class Eos_bf_poly ).
  • 2 = 2-fluid Newtonian polytropic EOS (class Eos_bf_poly_newt ).
  • 3 = 2-fluid tabulated EOS (class Eos_bf_tabul). The second line in the file should contain a name given by the user to the EOS. The following lines should contain the EOS parameters (one parameter per line), in the same order than in the class declaration.

Definition at line 190 of file eos_bf_file.C.

◆ get_beta()

double Lorene::Eos_bf_poly::get_beta ( ) const
inline

Returns the coefficient $\beta$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$.

Definition at line 953 of file eos_bifluid.h.

References beta.

◆ get_gam1()

double Lorene::Eos_bf_poly::get_gam1 ( ) const
inline

Returns the adiabatic index $\gamma_1$.

Definition at line 914 of file eos_bifluid.h.

References gam1.

◆ get_gam2()

double Lorene::Eos_bf_poly::get_gam2 ( ) const
inline

Returns the adiabatic index $\gamma_2$.

Definition at line 917 of file eos_bifluid.h.

References gam2.

◆ get_gam3()

double Lorene::Eos_bf_poly::get_gam3 ( ) const
inline

Returns the adiabatic index $\gamma_3$.

Definition at line 920 of file eos_bifluid.h.

References gam3.

◆ get_gam4()

double Lorene::Eos_bf_poly::get_gam4 ( ) const
inline

Returns the adiabatic index $\gamma_4$.

Definition at line 923 of file eos_bifluid.h.

References gam4.

◆ get_gam5()

double Lorene::Eos_bf_poly::get_gam5 ( ) const
inline

Returns the adiabatic index $\gamma_5$.

Definition at line 926 of file eos_bifluid.h.

References gam5.

◆ get_gam6()

double Lorene::Eos_bf_poly::get_gam6 ( ) const
inline

Returns the adiabatic index $\gamma_6$.

Definition at line 929 of file eos_bifluid.h.

References gam6.

◆ get_K11()

double Lorene::Eos_bf_poly::get_K11 ( const double  n1,
const double  n2,
const double  delta2 
) const
virtual

Computes the derivative of the energy with respect to (baryonic density 1) $^2$.

Parameters
n1[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density of fluid 1 at which the derivative is to be computed
n2[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density of fluid 2 at which the derivative is to be computed
x[input, unit $n_{\rm nuc}^2c^2$] relative Lorentz factor $\times$both densities at which the derivative is to be computed
Returns
derivative $K^{11}=2\frac{\partial{\cal{E}}}{\partial{n_1^2}}$

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 933 of file eos_bf_poly.C.

References beta, gam1, gam3, gam4, gam5, gam6, kap1, kap3, Lorene::Eos_bifluid::m_1, and Lorene::pow().

◆ get_K12()

double Lorene::Eos_bf_poly::get_K12 ( const double  n1,
const double  n2,
const double  delta2 
) const
virtual

Computes the derivative of the energy with respect to $x^2=n_1n_2\Gamma_\Delta$.

Parameters
n1[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density of fluid 1 at which the derivative is to be computed
n2[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density of fluid 2 at which the derivative is to be computed
x[input, unit $n_{\rm nuc}^2c^2$] relative Lorentz factor $\times$both densities at which the derivative is to be computed
Returns
derivative $K^{12}=\frac{\partial {\cal E}}{\partial (n_1n_2\Gamma_\Delta)}$

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 963 of file eos_bf_poly.C.

References beta, gam5, gam6, and Lorene::pow().

◆ get_K22()

double Lorene::Eos_bf_poly::get_K22 ( const double  n1,
const double  n2,
const double  delta2 
) const
virtual

Computes the derivative of the energy/(baryonic density 2) $^2$.

Parameters
n1[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density of fluid 1 at which the derivative is to be computed
n2[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density of fluid 2 at which the derivative is to be computed
x[input, unit $n_{\rm nuc}^2c^2$] relative Lorentz factor $\times$both densities at which the derivative is to be computed
Returns
derivative $K^{22} = 2\frac{\partial {\cal E}}{\partial n_2^2}$

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 948 of file eos_bf_poly.C.

References beta, gam2, gam3, gam4, gam5, gam6, kap2, kap3, Lorene::Eos_bifluid::m_2, and Lorene::pow().

◆ get_kap1()

double Lorene::Eos_bf_poly::get_kap1 ( ) const
inline

Returns the pressure coefficient $\kappa_1$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$.

Definition at line 935 of file eos_bifluid.h.

References kap1.

◆ get_kap2()

double Lorene::Eos_bf_poly::get_kap2 ( ) const
inline

Returns the pressure coefficient $\kappa_2$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$.

Definition at line 941 of file eos_bifluid.h.

References kap2.

◆ get_kap3()

double Lorene::Eos_bf_poly::get_kap3 ( ) const
inline

Returns the pressure coefficient $\kappa_3$
[unit: $\rho_{\rm nuc} c^2 $], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$.

Definition at line 947 of file eos_bifluid.h.

References kap3.

◆ get_Knn()

Cmp Lorene::Eos_bifluid::get_Knn ( const Cmp nbar1,
const Cmp nbar2,
const Cmp x2,
int  nzet,
int  l_min = 0 
) const
inherited

Computes the derivatives of the energy/(baryonic density 1) $^2$.

Parameters
nbar1[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 1 at which the derivatives are to be computed
nbar2[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 2 at which the derivatives are to be computed
x2[input, unit $n_{\rm nuc}^2c^2$] relative velocity $\times$both densities at which the derivative is to be computed
nzetnumber of domains where the derivatives are to be computed.
l_minindex of the innermost domain is which the derivatives are to be computed [default value: 0]; the derivatives are computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.
Returns
derivative $K^{11}$ field (see get_K11 )

Definition at line 750 of file eos_bifluid.C.

References Lorene::Eos_bifluid::calcule(), Lorene::Eos_bifluid::get_K11(), and Lorene::Cmp::get_mp().

◆ get_Knp()

Cmp Lorene::Eos_bifluid::get_Knp ( const Cmp nbar1,
const Cmp nbar2,
const Cmp x2,
int  nzet,
int  l_min = 0 
) const
inherited

Computes the derivatives of the energy with respect to $x^2=n_1n_2\Gamma_\Delta^2$.

Parameters
nbar1[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 1 at which the derivatives are to be computed
nbar2[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 2 at which the derivatives are to be computed
x2[input, unit $n_{\rm nuc}^2c^2$] relative velocity $\times$both densities at which the derivative is to be computed
nzetnumber of domains where the derivatives are to be computed.
l_minindex of the innermost domain is which the derivatives are to be computed [default value: 0]; the derivatives are computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.
Returns
derivative $K^{12}$ field (see get_K12 )

Definition at line 761 of file eos_bifluid.C.

References Lorene::Eos_bifluid::calcule(), Lorene::Eos_bifluid::get_K12(), and Lorene::Cmp::get_mp().

◆ get_Kpp()

Cmp Lorene::Eos_bifluid::get_Kpp ( const Cmp nbar1,
const Cmp nbar2,
const Cmp x2,
int  nzet,
int  l_min = 0 
) const
inherited

Computes the derivatives of the energy/(baryonic density 2) $^2$.

Parameters
nbar1[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 1 at which the derivatives are to be computed
nbar2[input, unit $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$] baryonic density field of fluid 2 at which the derivatives are to be computed
x2[input, unit $n_{\rm nuc}^2c^2$] relative velocity $\times$both densities at which the derivative is to be computed
nzetnumber of domains where the derivatives are to be computed.
l_minindex of the innermost domain is which the derivatives are to be computed [default value: 0]; the derivatives are computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.
Returns
derivative $K^{22}$ field (see get_K12 )

Definition at line 772 of file eos_bifluid.C.

References Lorene::Eos_bifluid::calcule(), Lorene::Eos_bifluid::get_K22(), and Lorene::Cmp::get_mp().

◆ get_m1()

double Lorene::Eos_bifluid::get_m1 ( ) const
inlineinherited

Return the individual particule mass $m_1$.

[unit: $m_B = 1.66\ 10^{-27} \ {\rm kg}$].

Definition at line 263 of file eos_bifluid.h.

References Lorene::Eos_bifluid::m_1.

◆ get_m2()

double Lorene::Eos_bifluid::get_m2 ( ) const
inlineinherited

Return the individual particule mass $m_2$.

[unit: $m_B = 1.66\ 10^{-27} \ {\rm kg}$].

Definition at line 269 of file eos_bifluid.h.

References Lorene::Eos_bifluid::m_2.

◆ get_name()

string Lorene::Eos_bifluid::get_name ( ) const
inlineinherited

Returns the EOS name.

Definition at line 253 of file eos_bifluid.h.

References Lorene::Eos_bifluid::name.

◆ identify()

int Lorene::Eos_bf_poly::identify ( ) const
virtual

Returns a number to identify the sub-classe of Eos_bifluid
the object belongs to.

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 101 of file eos_bf_file.C.

◆ nbar_ent()

void Lorene::Eos_bifluid::nbar_ent ( const Cmp ent1,
const Cmp ent2,
const Cmp delta2,
Cmp nbar1,
Cmp nbar2,
int  nzet,
int  l_min = 0 
) const
inherited

Computes both baryon density fields from the log-enthalpy fields and the relative velocity.

Parameters
ent1[input, unit: $c^2$] log-enthalpy $H_1$
ent2[input, unit: $c^2$] log-enthalpy $H_2$
delta2[input, unit: $c^2$] relative velocity $\Delta^2$
nbar1[output] baryonic density of the first fluid
nbar2[output] baryonic density of the second fluid [unit: $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$]
nzetnumber of domains where the baryon density is to be computed.
l_minindex of the innermost domain is which the baryon density is to be computed [default value: 0]; the baryon density is computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.

Definition at line 416 of file eos_bifluid.C.

References Lorene::Cmp::get_etat().

◆ nbar_ent_p()

bool Lorene::Eos_bf_poly::nbar_ent_p ( const double  ent1,
const double  ent2,
const double  delta2,
double &  nbar1,
double &  nbar2 
) const
virtual

Computes both baryon densities from the log-enthalpies.

Parameters
ent1[input, unit: $c^2$] log-enthalpy $H_1$
ent2[input, unit: $c^2$] log-enthalpy $H_2$
delta2[input, unit: $c^2$] relative velocity $\Delta^2$
nbar1[output] baryonic density of the first fluid
nbar2[output] baryonic density of the second fluid [unit: $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$]

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 510 of file eos_bf_poly.C.

References Lorene::Param::add_double(), Lorene::Param::add_double_mod(), beta, ecart, Lorene::exp(), gam1, gam1m1, gam2, gam2m1, gam3, gam4, gam5, gam6, Lorene::Param::get_double_mod(), kap1, kap2, kap3, Lorene::log(), Lorene::Eos_bifluid::m_1, Lorene::Eos_bifluid::m_2, Lorene::pow(), precis, relax, Lorene::sqrt(), typeos, and Lorene::zerosec_b().

◆ nbar_ent_p1()

double Lorene::Eos_bf_poly::nbar_ent_p1 ( const double  ent1) const
virtual

Computes baryon density out of the log-enthalpy asuming that only fluid 1 is present.

Parameters
ent1[input, unit: $c^2$] log-enthalpy $H_1$
Returns
nbar1 baryonic density of the first fluid

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 881 of file eos_bf_poly.C.

References Lorene::exp(), gam1, gam1m1, kap1, and Lorene::Eos_bifluid::m_1.

◆ nbar_ent_p2()

double Lorene::Eos_bf_poly::nbar_ent_p2 ( const double  ent2) const
virtual

Computes baryon density out of the log-enthalpy assuming that only fluid 2 is present.

Parameters
ent2[input, unit: $c^2$] log-enthalpy $H_1$
Returns
nbar1 baryonic density of the first fluid

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 885 of file eos_bf_poly.C.

References Lorene::exp(), gam2, gam2m1, kap2, and Lorene::Eos_bifluid::m_2.

◆ operator!=()

bool Lorene::Eos_bf_poly::operator!= ( const Eos_bifluid eos_i) const
virtual

Comparison operator (difference)

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 444 of file eos_bf_poly.C.

References operator==().

◆ operator=()

void Lorene::Eos_bf_poly::operator= ( const Eos_bf_poly eosi)

Assignment to another Eos_bf_poly.

Definition at line 297 of file eos_bf_poly.C.

References beta, ecart, gam1, gam2, gam3, kap1, kap2, kap3, Lorene::Eos_bifluid::operator=(), precis, relax, set_auxiliary(), and typeos.

◆ operator==()

bool Lorene::Eos_bf_poly::operator== ( const Eos_bifluid eos_i) const
virtual

Comparison operator (egality)

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 389 of file eos_bf_poly.C.

References beta, gam1, gam2, gam3, gam4, gam5, gam6, Lorene::Eos_bifluid::identify(), identify(), kap1, kap2, kap3, Lorene::Eos_bifluid::m_1, and Lorene::Eos_bifluid::m_2.

◆ operator>>()

ostream & Lorene::Eos_bf_poly::operator>> ( ostream &  ost) const
protectedvirtual

Operator >>

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 475 of file eos_bf_poly.C.

References beta, ecart, gam1, gam2, gam3, gam4, gam5, gam6, kap1, kap2, kap3, precis, relax, and typeos.

◆ press_ent()

Cmp Lorene::Eos_bifluid::press_ent ( const Cmp ent1,
const Cmp ent2,
const Cmp delta2,
int  nzet,
int  l_min = 0 
) const
inherited

Computes the pressure from the log-enthalpy fields and the relative velocity.

Parameters
ent1[input, unit: $c^2$] log-enthalpy $H_1$
ent2[input, unit: $c^2$] log-enthalpy $H_2$
delta2[input, unit: $c^2$] relative velocity $\Delta^2$
nzetnumber of domains where the pressure is to be computed.
l_minindex of the innermost domain is which the pressure is to be computed [default value: 0]; the pressure is computed only in domains whose indices are in [l_min,l_min+nzet-1] . In the other domains, it is set to zero.
Returns
pressure field [unit: $\rho_{\rm nuc} c^2$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$

Definition at line 584 of file eos_bifluid.C.

References Lorene::Cmp::get_etat(), and Lorene::Cmp::get_mp().

◆ press_nbar_p()

double Lorene::Eos_bf_poly::press_nbar_p ( const double  nbar1,
const double  nbar2,
const double  delta2 
) const
virtual

Computes the pressure from the baryonic densities and the relative velocity.

Parameters
nbar1[input] baryonic density of the first fluid
nbar2[input] baryonic density of the second fluid [unit: $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$]
delta2[input, unit: $c^2$] relative velocity $\Delta^2$
Returns
pressure p [unit: $\rho_{\rm nuc} c^2$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 912 of file eos_bf_poly.C.

References beta, gam1, gam1m1, gam2, gam2m1, gam3, gam34m1, gam4, gam5, gam56m1, gam6, kap1, kap2, kap3, and Lorene::pow().

◆ sauve()

void Lorene::Eos_bf_poly::sauve ( FILE *  fich) const
virtual

Save in a file.

Reimplemented from Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 455 of file eos_bf_poly.C.

References beta, ecart, Lorene::fwrite_be(), gam1, gam2, gam3, gam4, gam5, gam6, kap1, kap2, kap3, precis, relax, and Lorene::Eos_bifluid::sauve().

◆ set_auxiliary()

void Lorene::Eos_bf_poly::set_auxiliary ( )
protected

Computes the auxiliary quantities gam1m1 , gam2m1 and gam3m1.

Definition at line 323 of file eos_bf_poly.C.

References gam1, gam1m1, gam2, gam2m1, gam3, gam34m1, gam4, gam5, gam56m1, gam6, kap1, kap2, and kap3.

◆ trans2Eos()

Eos * Lorene::Eos_bf_poly::trans2Eos ( ) const
virtual

Makes a translation from Eos_bifluid to Eos .

This is only useful for the construction of a Et_rot_bifluid star and ought not to be used in other situations.

Implements Lorene::Eos_bifluid.

Reimplemented in Lorene::Eos_bf_poly_newt.

Definition at line 978 of file eos_bf_poly.C.

References gam1, kap1, and Lorene::Eos_bifluid::m_1.

Friends And Related Function Documentation

◆ Eos_bifluid::eos_from_file

The construction functions from a file.

Member Data Documentation

◆ beta

double Lorene::Eos_bf_poly::beta
protected

Coefficient $\beta$ , see Eq.

~eeosbfpolye} [unit: $\rho_{\rm nuc} c^2 / n_{\rm nuc}^\gamma$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$ and $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$.

Definition at line 767 of file eos_bifluid.h.

◆ ecart

double Lorene::Eos_bf_poly::ecart
protected

contains the precision required in the relaxation nbar_ent_p

Definition at line 807 of file eos_bifluid.h.

◆ gam1

double Lorene::Eos_bf_poly::gam1
protected

Adiabatic indexes $\gamma_1$, see Eq.~eeosbfpolye}.

Definition at line 724 of file eos_bifluid.h.

◆ gam1m1

double Lorene::Eos_bf_poly::gam1m1
protected

$\gamma_1-1$

Definition at line 769 of file eos_bifluid.h.

◆ gam2

double Lorene::Eos_bf_poly::gam2
protected

Adiabatic indexes $\gamma_2$, see Eq.~eeosbfpolye}.

Definition at line 727 of file eos_bifluid.h.

◆ gam2m1

double Lorene::Eos_bf_poly::gam2m1
protected

$\gamma_2-1$

Definition at line 770 of file eos_bifluid.h.

◆ gam3

double Lorene::Eos_bf_poly::gam3
protected

Adiabatic indexes $\gamma_3$, see Eq.~eeosbfpolye}.

Definition at line 730 of file eos_bifluid.h.

◆ gam34m1

double Lorene::Eos_bf_poly::gam34m1
protected

$\gamma_3+\gamma_4-1$

Definition at line 771 of file eos_bifluid.h.

◆ gam4

double Lorene::Eos_bf_poly::gam4
protected

Adiabatic indexes $\gamma_4$, see Eq.~eeosbfpolye}.

Definition at line 733 of file eos_bifluid.h.

◆ gam5

double Lorene::Eos_bf_poly::gam5
protected

Adiabatic indexes $\gamma_5$, see Eq.~eeosbfpolye}.

Definition at line 736 of file eos_bifluid.h.

◆ gam56m1

double Lorene::Eos_bf_poly::gam56m1
protected

$\gamma_5+\gamma_6-1$

Definition at line 772 of file eos_bifluid.h.

◆ gam6

double Lorene::Eos_bf_poly::gam6
protected

Adiabatic indexes $\gamma_6$, see Eq.~eeosbfpolye}.

Definition at line 739 of file eos_bifluid.h.

◆ kap1

double Lorene::Eos_bf_poly::kap1
protected

Pressure coefficient $\kappa_1$ , see Eq.

~eeosbfpolye} [unit: $\rho_{\rm nuc} c^2 / n_{\rm nuc}^\gamma$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$ and $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$.

Definition at line 746 of file eos_bifluid.h.

◆ kap2

double Lorene::Eos_bf_poly::kap2
protected

Pressure coefficient $\kappa_2$ , see Eq.

~eeosbfpolye} [unit: $\rho_{\rm nuc} c^2 / n_{\rm nuc}^\gamma$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$ and $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$.

Definition at line 753 of file eos_bifluid.h.

◆ kap3

double Lorene::Eos_bf_poly::kap3
protected

Pressure coefficient $\kappa_3$ , see Eq.

~eeosbfpolye} [unit: $\rho_{\rm nuc} c^2 / n_{\rm nuc}^\gamma$], where $\rho_{\rm nuc} := 1.66\ 10^{17} \ {\rm kg/m}^3$ and $n_{\rm nuc} := 0.1 \ {\rm fm}^{-3}$.

Definition at line 760 of file eos_bifluid.h.

◆ m_1

double Lorene::Eos_bifluid::m_1
protectedinherited

Individual particle mass $m_1$
[unit: $m_B = 1.66\ 10^{-27} \ {\rm kg}$].

Definition at line 191 of file eos_bifluid.h.

◆ m_2

double Lorene::Eos_bifluid::m_2
protectedinherited

Individual particle mass $m_2$
[unit: $m_B = 1.66\ 10^{-27} \ {\rm kg}$].

Definition at line 196 of file eos_bifluid.h.

◆ name

string Lorene::Eos_bifluid::name
protectedinherited

EOS name.

Definition at line 186 of file eos_bifluid.h.

◆ precis

double Lorene::Eos_bf_poly::precis
protected

contains the precision required in zerosec_b

Definition at line 804 of file eos_bifluid.h.

◆ relax

double Lorene::Eos_bf_poly::relax
protected

Parameters needed for some inversions of the EOS.

In particular, it is used for type 4 EOS: contains the relaxation parameter needed in the iteration

Definition at line 802 of file eos_bifluid.h.

◆ typeos

int Lorene::Eos_bf_poly::typeos
protected

The bi-fluid analytical EOS type:

0 - $\gamma_1 = \gamma_2 = 2$ and $\gamma_3 = \gamma_4 = \gamma_5 = \gamma_6 = 1$. In this case, the EOS can be inverted analytically.

1 - $\gamma_3 = \gamma_4 = \gamma_5 = \gamma_6 = 1$, but $\gamma_1 \not= 2$ or $\gamma_2 \not= 2$.

2 - $\gamma_3 = \gamma_5 = 1$, but none of the previous cases.

3 - $\gamma_4 = \gamma_6 = 1$, but none of the previous cases.

4 - None of the previous cases (the most general)

5 - special case of comparison to slow-rotation approximation: this is identical to typeos=0, but using a modified EOS-inversion method, namely we don't switch to a 1-fluid EOS in 1-fluid regions.

Definition at line 796 of file eos_bifluid.h.


The documentation for this class was generated from the following files: