Etoile Class Reference
[Stars and black holes]

Base class for stars *** DEPRECATED : use class Star instead ***. More...

#include <etoile.h>

Inheritance diagram for Etoile:
Etoile_bin Etoile_rot Et_bin_bhns_extr Et_bin_nsbh Et_rot_bifluid Et_rot_diff Et_rot_mag Et_magnetisation

List of all members.

Public Member Functions

 Etoile (Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
 Standard constructor.
 Etoile (const Etoile &)
 Copy constructor.
 Etoile (Map &mp_i, const Eos &eos_i, FILE *fich)
 Constructor from a file (see sauve(FILE*) ).
virtual ~Etoile ()
 Destructor.
void operator= (const Etoile &)
 Assignment to another Etoile.
Mapset_mp ()
 Read/write of the mapping.
void set_enthalpy (const Cmp &)
 Assignment of the enthalpy field.
virtual void equation_of_state ()
 Computes the proper baryon and energy density, as well as pressure from the enthalpy.
virtual void hydro_euler ()
 Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame (nbar , ener and press ).
virtual void equilibrium_spher (double ent_c, double precis=1.e-14, const Tbl *ent_limit=0x0)
 Computes a spherical static configuration.
void equil_spher_regular (double ent_c, double precis=1.e-14)
 Computes a spherical static configuration.
virtual void equil_spher_falloff (double ent_c, double precis=1.e-14)
 Computes a spherical static configuration with the outer boundary condition at a finite radius.
const Mapget_mp () const
 Returns the mapping.
int get_nzet () const
 Returns the number of domains occupied by the star.
bool is_relativistic () const
 Returns true for a relativistic star, false for a Newtonian one.
const Eosget_eos () const
 Returns the equation of state.
const Tenseurget_ent () const
 Returns the enthalpy field.
const Tenseurget_nbar () const
 Returns the proper baryon density.
const Tenseurget_ener () const
 Returns the proper total energy density.
const Tenseurget_press () const
 Returns the fluid pressure.
const Tenseurget_ener_euler () const
 Returns the total energy density with respect to the Eulerian observer.
const Tenseurget_s_euler () const
 Returns the trace of the stress tensor in the Eulerian frame.
const Tenseurget_gam_euler () const
 Returns the Lorentz factor between the fluid and Eulerian observers.
const Tenseurget_u_euler () const
 Returns the fluid 3-velocity with respect to the Eulerian observer.
const Tenseurget_logn_auto () const
 Returns the logarithm of the part of the lapse N generated principaly by the star.
const Tenseurget_logn_auto_regu () const
 Returns the regular part of the logarithm of the part of the lapse N generated principaly by the star.
const Tenseurget_logn_auto_div () const
 Returns the divergent part of the logarithm of the part of the lapse N generated principaly by the star.
const Tenseurget_d_logn_auto_div () const
 Returns the gradient of logn_auto_div.
const Tenseurget_beta_auto () const
 Returns the logarithm of the part of the product AN generated principaly by the star.
const Tenseurget_nnn () const
 Returns the total lapse function N.
const Tenseurget_shift () const
 Returns the total shift vector $N^i$.
const Tenseurget_a_car () const
 Returns the total conformal factor $A^2$.
virtual void sauve (FILE *) const
 Save in a file.
double ray_eq () const
 Coordinate radius at $\phi=0$, $\theta=\pi/2$ [r_unit].
double ray_eq_pis2 () const
 Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$ [r_unit].
double ray_eq_pi () const
 Coordinate radius at $\phi=\pi$, $\theta=\pi/2$ [r_unit].
double ray_eq_3pis2 () const
 Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$ [r_unit].
double ray_pole () const
 Coordinate radius at $\theta=0$ [r_unit].
double ray_eq (int kk) const
 Coordinate radius at $\phi=2k\pi/np$, $\theta=\pi/2$ [r_unit].
virtual const Itbll_surf () const
 Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$.
const Tblxi_surf () const
 Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.
virtual double mass_b () const
 Baryon mass.
virtual double mass_g () const
 Gravitational mass.

Protected Member Functions

virtual void del_deriv () const
 Deletes all the derived quantities.
virtual void set_der_0x0 () const
 Sets to 0x0 all the pointers on derived quantities.
virtual void del_hydro_euler ()
 Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.
virtual ostream & operator>> (ostream &) const
 Operator >> (virtual function called by the operator <<).

Protected Attributes

Mapmp
 Mapping associated with the star.
int nzet
 Number of domains of *mp occupied by the star.
bool relativistic
 Indicator of relativity: true for a relativistic star, false for a Newtonian one.
double unsurc2
 $1/c^2$ : unsurc2=1 for a relativistic star, 0 for a Newtonian one.
int k_div
 Index of regularity of the gravitational potential logn_auto .
const Eoseos
 Equation of state of the stellar matter.
Tenseur ent
 Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Tenseur nbar
 Baryon density in the fluid frame.
Tenseur ener
 Total energy density in the fluid frame.
Tenseur press
 Fluid pressure.
Tenseur ener_euler
 Total energy density in the Eulerian frame.
Tenseur s_euler
 Trace of the stress tensor in the Eulerian frame.
Tenseur gam_euler
 Lorentz factor between the fluid and Eulerian observers.
Tenseur u_euler
 Fluid 3-velocity with respect to the Eulerian observer.
Tenseur logn_auto
 Total of the logarithm of the part of the lapse N generated principaly by the star.
Tenseur logn_auto_regu
 Regular part of the logarithm of the part of the lapse N generated principaly by the star.
Tenseur logn_auto_div
 Divergent part (if k_div!=0 ) of the logarithm of the part of the lapse N generated principaly by the star.
Tenseur d_logn_auto_div
 Gradient of logn_auto_div (if k_div!=0 ).
Tenseur beta_auto
 Logarithm of the part of the product AN generated principaly by by the star.
Tenseur nnn
 Total lapse function.
Tenseur shift
 Total shift vector.
Tenseur a_car
 Total conformal factor $A^2$.
double * p_ray_eq
 Coordinate radius at $\phi=0$, $\theta=\pi/2$.
double * p_ray_eq_pis2
 Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$.
double * p_ray_eq_pi
 Coordinate radius at $\phi=\pi$, $\theta=\pi/2$.
double * p_ray_eq_3pis2
 Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$.
double * p_ray_pole
 Coordinate radius at $\theta=0$.
Itblp_l_surf
 Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$.
Tblp_xi_surf
 Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.
double * p_mass_b
 Baryon mass.
double * p_mass_g
 Gravitational mass.

Friends

ostream & operator<< (ostream &, const Etoile &)
 Display.

Detailed Description

Base class for stars *** DEPRECATED : use class Star instead ***.

()

An Etoile is constructed upon (i) a mapping (derived class of Map ), the center of which defines the center of the star, and (ii) an equation of state (derived class of Eos ). It contains tensor fields (class Tenseur ) which describle the hydrodynamical quantities as well as the gravitational field (spacetime metric).

According to the 3+1 formalism, the spacetime metric is written

\[ \label{eetoilemetrique} ds^2 = - (N^2 - N_i N^i) dt^2 - 2 N_i \, dt\, dx^i + A^2 \, {\tilde h}_{ij} \, dx^i dx^j \]

where ${\tilde h}_{ij}$ is a 3-metric, the exact form of which is specified in the derived classes of Etoile . The base class Etoile by itself provides only storage for the lapse function N (member nnn ), the shift vector $N^i$ (member shift ) and the conformal factor $A^2$ (member a_car ).

The 3+1 formalism introduces two kinds of priviledged observers: the fluid comoving observer and the Eulerian observer, whose 4-velocity is the unit future directed normal to the t = const hypersurfaces. The hydrodynamical quantities measured by the fluid observer correspond to the members ent , nbar , ener , and press . The hydrodynamical quantities measured by the Eulerian observer correspond to the members ener_euler , s_euler , gam_euler , and u_euler .

A star of class Etoile can be either relativistic or Newtonian, depending on the boolean indicator relativistic . For a Newtonian star, the metric coefficients N and A are set to 1, and $N^i$ is set to zero; the only relevant gravitational quantity in this case is logn_auto which represents the (regular part of) the Newtonian gravitational potential generated by the star.

Definition at line 411 of file etoile.h.


Constructor & Destructor Documentation

Etoile::Etoile ( Map mp_i,
int  nzet_i,
bool  relat,
const Eos eos_i 
)

Standard constructor.

Parameters:
mp_i Mapping on which the star will be defined
nzet_i Number of domains occupied by the star
relat should be true for a relativistic star, false for a Newtonian one
eos_i Equation of state of the stellar matter

Definition at line 139 of file etoile.C.

References a_car, beta_auto, d_logn_auto_div, ener, ener_euler, ent, eos, gam_euler, logn_auto, logn_auto_div, logn_auto_regu, nbar, nnn, press, relativistic, s_euler, set_der_0x0(), Tenseur::set_std_base(), shift, u_euler, and unsurc2.

Etoile::Etoile ( const Etoile et  ) 

Copy constructor.

Definition at line 241 of file etoile.C.

References set_der_0x0().

Etoile::Etoile ( Map mp_i,
const Eos eos_i,
FILE *  fich 
)

Constructor from a file (see sauve(FILE*) ).

Parameters:
mp_i Mapping on which the star will be defined
eos_i Equation of state of the stellar matter
fich input file (must have been created by the function sauve )

Definition at line 271 of file etoile.C.

References beta_auto, d_logn_auto_div, ent, eos, Eos::eos_from_file(), fread_be(), Map::get_bvect_spher(), k_div, logn_auto, logn_auto_div, logn_auto_regu, mp, nzet, relativistic, set_der_0x0(), shift, and unsurc2.

Etoile::~Etoile (  )  [virtual]

Destructor.

Definition at line 363 of file etoile.C.

References del_deriv().


Member Function Documentation

void Etoile::del_deriv (  )  const [protected, virtual]

Deletes all the derived quantities.

Reimplemented in Et_rot_bifluid, Et_rot_mag, Etoile_bin, and Etoile_rot.

Definition at line 374 of file etoile.C.

References p_l_surf, p_mass_b, p_mass_g, p_ray_eq, p_ray_eq_3pis2, p_ray_eq_pi, p_ray_eq_pis2, p_ray_pole, p_xi_surf, and set_der_0x0().

void Etoile::del_hydro_euler (  )  [protected, virtual]

Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer.

Reimplemented in Et_rot_bifluid, Et_rot_mag, Etoile_bin, and Etoile_rot.

Definition at line 406 of file etoile.C.

References del_deriv(), ener_euler, gam_euler, s_euler, Tenseur::set_etat_nondef(), and u_euler.

void Etoile::equation_of_state (  )  [virtual]
void Etoile::equil_spher_falloff ( double  ent_c,
double  precis = 1.e-14 
) [virtual]

Computes a spherical static configuration with the outer boundary condition at a finite radius.

Parameters:
ent_c [input] central value of the enthalpy
precis [input] threshold in the relative difference between the enthalpy fields of two consecutive steps to stop the iterative procedure (default value: 1.e-14)

Definition at line 53 of file etoile_eqsph_falloff.C.

References a_car, Tenseur::annule(), beta_auto, diffrel(), Cmp::dsdr(), Map_af::dsdr(), ener, ener_euler, ent, equation_of_state(), exp(), gam_euler, Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Map_af::homothetie(), logn_auto, mass_b(), mass_g(), mp, nbar, nnn, norme(), nzet, press, relativistic, s_euler, Tenseur::set(), Tenseur::set_etat_qcq(), Tenseur::set_std_base(), shift, sqrt(), u_euler, unsurc2, and Map::val_r().

void Etoile::equil_spher_regular ( double  ent_c,
double  precis = 1.e-14 
)

Computes a spherical static configuration.

The sources for Poisson equations are regularized by extracting analytical diverging parts.

Parameters:
ent_c [input] central value of the enthalpy
precis [input] threshold in the relative difference between the enthalpy fields of two consecutive steps to stop the iterative procedure (default value: 1.e-14)

Definition at line 110 of file et_equil_spher_regu.C.

References a_car, Tenseur::annule(), beta_auto, d_logn_auto_div, Eos::der_ener_ent_p(), Eos::der_nbar_ent_p(), diffrel(), Map_af::dsdr(), ener, ener_euler, ent, eos, equation_of_state(), exp(), gam_euler, Map::get_bvect_spher(), Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Tenseur::gradient_spher(), Map_af::homothetie(), k_div, logn_auto, logn_auto_div, logn_auto_regu, mass_b(), mass_g(), mp, nbar, nnn, norme(), nzet, Map_af::poisson(), Map_af::poisson_regular(), press, relativistic, s_euler, Tenseur::set(), Tenseur::set_etat_qcq(), Tenseur::set_std_base(), shift, sqrt(), Cmp::std_base_scal(), u_euler, unsurc2, and Map::val_r().

void Etoile::equilibrium_spher ( double  ent_c,
double  precis = 1.e-14,
const Tbl ent_limit = 0x0 
) [virtual]

Computes a spherical static configuration.

Parameters:
ent_c [input] central value of the enthalpy
precis [input] threshold in the relative difference between the enthalpy fields of two consecutive steps to stop the iterative procedure (default value: 1.e-14)
ent_limit [input] : array of enthalpy values to be set at the boundaries between the domains; if set to 0x0 (default), the initial values will be kept.

Definition at line 83 of file etoile_equil_spher.C.

References a_car, Map_et::adapt(), Param::add_double(), Param::add_int(), Param::add_int_mod(), Param::add_tbl(), Tenseur::annule(), beta_auto, diffrel(), Cmp::dsdr(), Map_af::dsdr(), ener, ener_euler, ent, equation_of_state(), exp(), gam_euler, Map_et::get_alpha(), Map_af::get_alpha(), Map_et::get_beta(), Map_af::get_beta(), get_ent(), Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), get_press(), Map_af::homothetie(), logn_auto, mass_b(), mass_g(), mp, nbar, nnn, norme(), nzet, Map_af::poisson(), press, relativistic, s_euler, Tenseur::set(), Map_af::set_alpha(), Map_af::set_beta(), Tenseur::set_etat_qcq(), Tenseur::set_std_base(), shift, sqrt(), u_euler, unsurc2, and Map::val_r().

const Tenseur& Etoile::get_a_car (  )  const [inline]

Returns the total conformal factor $A^2$.

Definition at line 718 of file etoile.h.

References a_car.

const Tenseur& Etoile::get_beta_auto (  )  const [inline]

Returns the logarithm of the part of the product AN generated principaly by the star.

Definition at line 709 of file etoile.h.

References beta_auto.

const Tenseur& Etoile::get_d_logn_auto_div (  )  const [inline]

Returns the gradient of logn_auto_div.

Definition at line 704 of file etoile.h.

References d_logn_auto_div.

const Tenseur& Etoile::get_ener (  )  const [inline]

Returns the proper total energy density.

Definition at line 664 of file etoile.h.

References ener.

const Tenseur& Etoile::get_ener_euler (  )  const [inline]

Returns the total energy density with respect to the Eulerian observer.

Definition at line 670 of file etoile.h.

References ener_euler.

const Tenseur& Etoile::get_ent (  )  const [inline]

Returns the enthalpy field.

Definition at line 658 of file etoile.h.

References ent.

const Eos& Etoile::get_eos (  )  const [inline]

Returns the equation of state.

Reimplemented in Et_rot_bifluid.

Definition at line 655 of file etoile.h.

References eos.

const Tenseur& Etoile::get_gam_euler (  )  const [inline]

Returns the Lorentz factor between the fluid and Eulerian observers.

Definition at line 676 of file etoile.h.

References gam_euler.

const Tenseur& Etoile::get_logn_auto (  )  const [inline]

Returns the logarithm of the part of the lapse N generated principaly by the star.

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 686 of file etoile.h.

References logn_auto.

const Tenseur& Etoile::get_logn_auto_div (  )  const [inline]

Returns the divergent part of the logarithm of the part of the lapse N generated principaly by the star.

In the Newtonian case, this is the diverging part of the Newtonian gravitational potential (in units of $c^2$).

Definition at line 700 of file etoile.h.

References logn_auto_div.

const Tenseur& Etoile::get_logn_auto_regu (  )  const [inline]

Returns the regular part of the logarithm of the part of the lapse N generated principaly by the star.

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 693 of file etoile.h.

References logn_auto_regu.

const Map& Etoile::get_mp (  )  const [inline]

Returns the mapping.

Definition at line 644 of file etoile.h.

References mp.

const Tenseur& Etoile::get_nbar (  )  const [inline]

Returns the proper baryon density.

Definition at line 661 of file etoile.h.

References nbar.

const Tenseur& Etoile::get_nnn (  )  const [inline]

Returns the total lapse function N.

Definition at line 712 of file etoile.h.

References nnn.

int Etoile::get_nzet (  )  const [inline]

Returns the number of domains occupied by the star.

Definition at line 647 of file etoile.h.

References nzet.

const Tenseur& Etoile::get_press (  )  const [inline]

Returns the fluid pressure.

Definition at line 667 of file etoile.h.

References press.

const Tenseur& Etoile::get_s_euler (  )  const [inline]

Returns the trace of the stress tensor in the Eulerian frame.

Definition at line 673 of file etoile.h.

References s_euler.

const Tenseur& Etoile::get_shift (  )  const [inline]

Returns the total shift vector $N^i$.

Definition at line 715 of file etoile.h.

References shift.

const Tenseur& Etoile::get_u_euler (  )  const [inline]

Returns the fluid 3-velocity with respect to the Eulerian observer.

Definition at line 679 of file etoile.h.

References u_euler.

void Etoile::hydro_euler (  )  [virtual]

Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame (nbar , ener and press ).

Reimplemented in Et_rot_bifluid, Et_rot_diff, Etoile_bin, and Etoile_rot.

Definition at line 680 of file etoile.C.

bool Etoile::is_relativistic (  )  const [inline]

Returns true for a relativistic star, false for a Newtonian one.

Definition at line 652 of file etoile.h.

References relativistic.

const Itbl & Etoile::l_surf (  )  const [virtual]

Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$.

The stellar surface is defined as the location where the enthalpy (member ent ) vanishes.

Reimplemented in Et_rot_bifluid, and Etoile_rot.

Definition at line 71 of file etoile_global.C.

References ent, Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), mp, nzet, p_l_surf, and p_xi_surf.

double Etoile::mass_b (  )  const [virtual]

Baryon mass.

Reimplemented in Et_rot_bifluid, Etoile_bin, and Etoile_rot.

Definition at line 521 of file etoile_global.C.

References Tenseur::get_etat(), nbar, p_mass_b, and relativistic.

double Etoile::mass_g (  )  const [virtual]

Gravitational mass.

Reimplemented in Et_rot_bifluid, Et_rot_mag, Etoile_bin, and Etoile_rot.

Definition at line 547 of file etoile_global.C.

References mass_b(), p_mass_g, and relativistic.

void Etoile::operator= ( const Etoile et  ) 
ostream & Etoile::operator>> ( ostream &  ost  )  const [protected, virtual]

Operator >> (virtual function called by the operator <<).

Reimplemented in Et_bin_nsbh, Et_rot_bifluid, Et_rot_diff, Et_rot_mag, Et_magnetisation, Etoile_bin, and Etoile_rot.

Definition at line 507 of file etoile.C.

References a_car, ener, ent, eos, k_div, mass_b(), mass_g(), nbar, nnn, nzet, press, ray_eq(), ray_eq_pi(), ray_eq_pis2(), ray_pole(), and relativistic.

double Etoile::ray_eq ( int  kk  )  const

Coordinate radius at $\phi=2k\pi/np$, $\theta=\pi/2$ [r_unit].

Definition at line 436 of file etoile_global.C.

References Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), l_surf(), mp, Map::val_r(), and xi_surf().

double Etoile::ray_eq (  )  const

Coordinate radius at $\phi=0$, $\theta=\pi/2$ [r_unit].

Definition at line 116 of file etoile_global.C.

References Map::get_mg(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), l_surf(), mp, p_ray_eq, Map::val_r(), and xi_surf().

double Etoile::ray_eq_3pis2 (  )  const
double Etoile::ray_eq_pi (  )  const

Coordinate radius at $\phi=\pi$, $\theta=\pi/2$ [r_unit].

Definition at line 252 of file etoile_global.C.

References Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), l_surf(), mp, p_ray_eq_pi, ray_eq(), Map::val_r(), and xi_surf().

double Etoile::ray_eq_pis2 (  )  const

Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$ [r_unit].

Definition at line 165 of file etoile_global.C.

References Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), l_surf(), mp, p_ray_eq_pis2, Map::val_r(), and xi_surf().

double Etoile::ray_pole (  )  const

Coordinate radius at $\theta=0$ [r_unit].

Definition at line 411 of file etoile_global.C.

References Map::get_mg(), Mg3d::get_type_t(), l_surf(), mp, p_ray_pole, Map::val_r(), and xi_surf().

void Etoile::sauve ( FILE *  fich  )  const [virtual]
void Etoile::set_der_0x0 (  )  const [protected, virtual]

Sets to 0x0 all the pointers on derived quantities.

Reimplemented in Et_rot_bifluid, Et_rot_mag, Etoile_bin, and Etoile_rot.

Definition at line 392 of file etoile.C.

References p_l_surf, p_mass_b, p_mass_g, p_ray_eq, p_ray_eq_3pis2, p_ray_eq_pi, p_ray_eq_pis2, p_ray_pole, and p_xi_surf.

void Etoile::set_enthalpy ( const Cmp ent_i  ) 

Assignment of the enthalpy field.

Definition at line 461 of file etoile.C.

References del_deriv(), ent, and equation_of_state().

Map& Etoile::set_mp (  )  [inline]

Read/write of the mapping.

Definition at line 591 of file etoile.h.

References mp.

const Tbl & Etoile::xi_surf (  )  const

Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.

The stellar surface is defined as the location where the enthalpy (member ent ) vanishes.

Definition at line 97 of file etoile_global.C.

References l_surf(), p_l_surf, and p_xi_surf.


Friends And Related Function Documentation

ostream& operator<< ( ostream &  ,
const Etoile  
) [friend]

Display.


Member Data Documentation

Tenseur Etoile::a_car [protected]

Total conformal factor $A^2$.

Definition at line 503 of file etoile.h.

Logarithm of the part of the product AN generated principaly by by the star.

Definition at line 494 of file etoile.h.

Gradient of logn_auto_div (if k_div!=0 ).

Definition at line 489 of file etoile.h.

Tenseur Etoile::ener [protected]

Total energy density in the fluid frame.

Definition at line 448 of file etoile.h.

Total energy density in the Eulerian frame.

Definition at line 453 of file etoile.h.

Tenseur Etoile::ent [protected]

Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).

Definition at line 445 of file etoile.h.

const Eos& Etoile::eos [protected]

Equation of state of the stellar matter.

Reimplemented in Et_rot_bifluid.

Definition at line 440 of file etoile.h.

Lorentz factor between the fluid and Eulerian observers.

Definition at line 459 of file etoile.h.

int Etoile::k_div [protected]

Index of regularity of the gravitational potential logn_auto .

If k_div=0 , logn_auto contains the total potential generated principaly by the star, otherwise it should be supplemented by logn_auto_div .

Definition at line 438 of file etoile.h.

Total of the logarithm of the part of the lapse N generated principaly by the star.

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 472 of file etoile.h.

Divergent part (if k_div!=0 ) of the logarithm of the part of the lapse N generated principaly by the star.

Definition at line 485 of file etoile.h.

Regular part of the logarithm of the part of the lapse N generated principaly by the star.

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 479 of file etoile.h.

Map& Etoile::mp [protected]

Mapping associated with the star.

Definition at line 416 of file etoile.h.

Tenseur Etoile::nbar [protected]

Baryon density in the fluid frame.

Definition at line 447 of file etoile.h.

Tenseur Etoile::nnn [protected]

Total lapse function.

Definition at line 497 of file etoile.h.

int Etoile::nzet [protected]

Number of domains of *mp occupied by the star.

Definition at line 421 of file etoile.h.

Itbl* Etoile::p_l_surf [mutable, protected]

Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$.

Definition at line 527 of file etoile.h.

double* Etoile::p_mass_b [mutable, protected]

Baryon mass.

Definition at line 535 of file etoile.h.

double* Etoile::p_mass_g [mutable, protected]

Gravitational mass.

Definition at line 536 of file etoile.h.

double* Etoile::p_ray_eq [mutable, protected]

Coordinate radius at $\phi=0$, $\theta=\pi/2$.

Definition at line 509 of file etoile.h.

double* Etoile::p_ray_eq_3pis2 [mutable, protected]

Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$.

Definition at line 518 of file etoile.h.

double* Etoile::p_ray_eq_pi [mutable, protected]

Coordinate radius at $\phi=\pi$, $\theta=\pi/2$.

Definition at line 515 of file etoile.h.

double* Etoile::p_ray_eq_pis2 [mutable, protected]

Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$.

Definition at line 512 of file etoile.h.

double* Etoile::p_ray_pole [mutable, protected]

Coordinate radius at $\theta=0$.

Definition at line 521 of file etoile.h.

Tbl* Etoile::p_xi_surf [mutable, protected]

Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.

Definition at line 533 of file etoile.h.

Tenseur Etoile::press [protected]

Fluid pressure.

Definition at line 449 of file etoile.h.

bool Etoile::relativistic [protected]

Indicator of relativity: true for a relativistic star, false for a Newtonian one.

Definition at line 426 of file etoile.h.

Tenseur Etoile::s_euler [protected]

Trace of the stress tensor in the Eulerian frame.

Definition at line 456 of file etoile.h.

Tenseur Etoile::shift [protected]

Total shift vector.

Definition at line 500 of file etoile.h.

Tenseur Etoile::u_euler [protected]

Fluid 3-velocity with respect to the Eulerian observer.

Definition at line 462 of file etoile.h.

double Etoile::unsurc2 [protected]

$1/c^2$ : unsurc2=1 for a relativistic star, 0 for a Newtonian one.

Definition at line 431 of file etoile.h.


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

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1