LORENE
Lorene::Etoile Class Reference


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

#include <etoile.h>

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

Public Member Functions

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

Protected Member Functions

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

Protected Attributes

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

Friends

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

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 427 of file etoile.h.

Constructor & Destructor Documentation

◆ Etoile() [1/3]

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

Standard constructor.

Parameters
mp_iMapping on which the star will be defined
nzet_iNumber of domains occupied by the star
relatshould be true for a relativistic star, false for a Newtonian one
eos_iEquation of state of the stellar matter

Definition at line 146 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(), Lorene::Tenseur::set_std_base(), shift, u_euler, and unsurc2.

◆ Etoile() [2/3]

Lorene::Etoile::Etoile ( const Etoile et)

Copy constructor.

Definition at line 248 of file etoile.C.

References set_der_0x0().

◆ Etoile() [3/3]

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

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

Parameters
mp_iMapping on which the star will be defined
eos_iEquation of state of the stellar matter
fichinput file (must have been created by the function sauve )

Definition at line 278 of file etoile.C.

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

◆ ~Etoile()

Lorene::Etoile::~Etoile ( )
virtual

Destructor.

Definition at line 370 of file etoile.C.

References del_deriv().

Member Function Documentation

◆ del_deriv()

void Lorene::Etoile::del_deriv ( ) const
protectedvirtual

Deletes all the derived quantities.

Reimplemented in Lorene::Etoile_rot, Lorene::Etoile_bin, Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.

Definition at line 381 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().

◆ del_hydro_euler()

void Lorene::Etoile::del_hydro_euler ( )
protectedvirtual

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

Reimplemented in Lorene::Etoile_rot, Lorene::Etoile_bin, Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.

Definition at line 413 of file etoile.C.

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

◆ equation_of_state()

◆ equil_spher_falloff()

void Lorene::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 60 of file etoile_eqsph_falloff.C.

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

◆ equil_spher_regular()

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

◆ equilibrium_spher()

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

◆ get_a_car()

const Tenseur& Lorene::Etoile::get_a_car ( ) const
inline

Returns the total conformal factor $A^2$.

Definition at line 736 of file etoile.h.

References a_car.

◆ get_beta_auto()

const Tenseur& Lorene::Etoile::get_beta_auto ( ) const
inline

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

Definition at line 727 of file etoile.h.

References beta_auto.

◆ get_d_logn_auto_div()

const Tenseur& Lorene::Etoile::get_d_logn_auto_div ( ) const
inline

Returns the gradient of logn_auto_div.

Definition at line 722 of file etoile.h.

References d_logn_auto_div.

◆ get_ener()

const Tenseur& Lorene::Etoile::get_ener ( ) const
inline

Returns the proper total energy density.

Definition at line 682 of file etoile.h.

References ener.

◆ get_ener_euler()

const Tenseur& Lorene::Etoile::get_ener_euler ( ) const
inline

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

Definition at line 688 of file etoile.h.

References ener_euler.

◆ get_ent()

const Tenseur& Lorene::Etoile::get_ent ( ) const
inline

Returns the enthalpy field.

Definition at line 676 of file etoile.h.

References ent.

◆ get_eos()

const Eos& Lorene::Etoile::get_eos ( ) const
inline

Returns the equation of state.

Definition at line 673 of file etoile.h.

References eos.

◆ get_gam_euler()

const Tenseur& Lorene::Etoile::get_gam_euler ( ) const
inline

Returns the Lorentz factor between the fluid and Eulerian observers.

Definition at line 694 of file etoile.h.

References gam_euler.

◆ get_logn_auto()

const Tenseur& Lorene::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 704 of file etoile.h.

References logn_auto.

◆ get_logn_auto_div()

const Tenseur& Lorene::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 718 of file etoile.h.

References logn_auto_div.

◆ get_logn_auto_regu()

const Tenseur& Lorene::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 711 of file etoile.h.

References logn_auto_regu.

◆ get_mp()

const Map& Lorene::Etoile::get_mp ( ) const
inline

Returns the mapping.

Definition at line 662 of file etoile.h.

References mp.

◆ get_nbar()

const Tenseur& Lorene::Etoile::get_nbar ( ) const
inline

Returns the proper baryon density.

Definition at line 679 of file etoile.h.

References nbar.

◆ get_nnn()

const Tenseur& Lorene::Etoile::get_nnn ( ) const
inline

Returns the total lapse function N.

Definition at line 730 of file etoile.h.

References nnn.

◆ get_nzet()

int Lorene::Etoile::get_nzet ( ) const
inline

Returns the number of domains occupied by the star.

Definition at line 665 of file etoile.h.

References nzet.

◆ get_press()

const Tenseur& Lorene::Etoile::get_press ( ) const
inline

Returns the fluid pressure.

Definition at line 685 of file etoile.h.

References press.

◆ get_s_euler()

const Tenseur& Lorene::Etoile::get_s_euler ( ) const
inline

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

Definition at line 691 of file etoile.h.

References s_euler.

◆ get_shift()

const Tenseur& Lorene::Etoile::get_shift ( ) const
inline

Returns the total shift vector $N^i$.

Definition at line 733 of file etoile.h.

References shift.

◆ get_u_euler()

const Tenseur& Lorene::Etoile::get_u_euler ( ) const
inline

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

Definition at line 697 of file etoile.h.

References u_euler.

◆ hydro_euler()

void Lorene::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 Lorene::Etoile_rot, Lorene::Etoile_bin, Lorene::Et_rot_bifluid, and Lorene::Et_rot_diff.

Definition at line 687 of file etoile.C.

◆ is_relativistic()

bool Lorene::Etoile::is_relativistic ( ) const
inline

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

Definition at line 670 of file etoile.h.

References relativistic.

◆ l_surf()

const Itbl & Lorene::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 Lorene::Etoile_rot, and Lorene::Et_rot_bifluid.

Definition at line 78 of file etoile_global.C.

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

◆ mass_b()

double Lorene::Etoile::mass_b ( ) const
virtual

Baryon mass.

Reimplemented in Lorene::Etoile_rot, Lorene::Etoile_bin, and Lorene::Et_rot_bifluid.

Definition at line 528 of file etoile_global.C.

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

◆ mass_g()

double Lorene::Etoile::mass_g ( ) const
virtual

Gravitational mass.

Reimplemented in Lorene::Etoile_rot, Lorene::Etoile_bin, Lorene::Et_magnetisation, Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.

Definition at line 554 of file etoile_global.C.

References mass_b(), p_mass_g, and relativistic.

◆ operator=()

void Lorene::Etoile::operator= ( const Etoile et)

◆ operator>>()

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

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

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

Definition at line 514 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.

◆ ray_eq() [1/2]

double Lorene::Etoile::ray_eq ( ) const

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

Definition at line 123 of file etoile_global.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), mp, and p_ray_eq.

◆ ray_eq() [2/2]

double Lorene::Etoile::ray_eq ( int  kk) const

◆ ray_eq_3pis2()

double Lorene::Etoile::ray_eq_3pis2 ( ) const

◆ ray_eq_pi()

double Lorene::Etoile::ray_eq_pi ( ) const

◆ ray_eq_pis2()

double Lorene::Etoile::ray_eq_pis2 ( ) const

◆ ray_pole()

double Lorene::Etoile::ray_pole ( ) const

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

Definition at line 418 of file etoile_global.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_type_t(), mp, and p_ray_pole.

◆ sauve()

◆ set_der_0x0()

void Lorene::Etoile::set_der_0x0 ( ) const
protectedvirtual

Sets to 0x0 all the pointers on derived quantities.

Reimplemented in Lorene::Etoile_rot, Lorene::Etoile_bin, Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.

Definition at line 399 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.

◆ set_enthalpy()

void Lorene::Etoile::set_enthalpy ( const Cmp ent_i)

Assignment of the enthalpy field.

Definition at line 468 of file etoile.C.

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

◆ set_mp()

Map& Lorene::Etoile::set_mp ( )
inline

Read/write of the mapping.

Definition at line 604 of file etoile.h.

References mp.

◆ xi_surf()

const Tbl & Lorene::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 104 of file etoile_global.C.

References l_surf(), p_l_surf, and p_xi_surf.

Friends And Related Function Documentation

◆ operator<<

ostream& operator<< ( ostream &  ost,
const Etoile et 
)
friend

Display.

Definition at line 509 of file etoile.C.

Member Data Documentation

◆ a_car

Tenseur Lorene::Etoile::a_car
protected

Total conformal factor $A^2$.

Definition at line 518 of file etoile.h.

◆ beta_auto

Tenseur Lorene::Etoile::beta_auto
protected

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

Definition at line 509 of file etoile.h.

◆ d_logn_auto_div

Tenseur Lorene::Etoile::d_logn_auto_div
protected

Gradient of logn_auto_div (if k_div!=0 )

Definition at line 504 of file etoile.h.

◆ ener

Tenseur Lorene::Etoile::ener
protected

Total energy density in the fluid frame.

Definition at line 463 of file etoile.h.

◆ ener_euler

Tenseur Lorene::Etoile::ener_euler
protected

Total energy density in the Eulerian frame.

Definition at line 468 of file etoile.h.

◆ ent

Tenseur Lorene::Etoile::ent
protected

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

Definition at line 460 of file etoile.h.

◆ eos

const Eos& Lorene::Etoile::eos
protected

Equation of state of the stellar matter.

Definition at line 454 of file etoile.h.

◆ gam_euler

Tenseur Lorene::Etoile::gam_euler
protected

Lorentz factor between the fluid and Eulerian observers.

Definition at line 474 of file etoile.h.

◆ k_div

int Lorene::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 452 of file etoile.h.

◆ logn_auto

Tenseur Lorene::Etoile::logn_auto
protected

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 487 of file etoile.h.

◆ logn_auto_div

Tenseur Lorene::Etoile::logn_auto_div
protected

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

Definition at line 500 of file etoile.h.

◆ logn_auto_regu

Tenseur Lorene::Etoile::logn_auto_regu
protected

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 494 of file etoile.h.

◆ mp

Map& Lorene::Etoile::mp
protected

Mapping associated with the star.

Definition at line 432 of file etoile.h.

◆ nbar

Tenseur Lorene::Etoile::nbar
protected

Baryon density in the fluid frame.

Definition at line 462 of file etoile.h.

◆ nnn

Tenseur Lorene::Etoile::nnn
protected

Total lapse function.

Definition at line 512 of file etoile.h.

◆ nzet

int Lorene::Etoile::nzet
protected

Number of domains of *mp occupied by the star.

Definition at line 435 of file etoile.h.

◆ p_l_surf

Itbl* Lorene::Etoile::p_l_surf
mutableprotected

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 542 of file etoile.h.

◆ p_mass_b

double* Lorene::Etoile::p_mass_b
mutableprotected

Baryon mass.

Definition at line 550 of file etoile.h.

◆ p_mass_g

double* Lorene::Etoile::p_mass_g
mutableprotected

Gravitational mass.

Definition at line 551 of file etoile.h.

◆ p_ray_eq

double* Lorene::Etoile::p_ray_eq
mutableprotected

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

Definition at line 524 of file etoile.h.

◆ p_ray_eq_3pis2

double* Lorene::Etoile::p_ray_eq_3pis2
mutableprotected

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

Definition at line 533 of file etoile.h.

◆ p_ray_eq_pi

double* Lorene::Etoile::p_ray_eq_pi
mutableprotected

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

Definition at line 530 of file etoile.h.

◆ p_ray_eq_pis2

double* Lorene::Etoile::p_ray_eq_pis2
mutableprotected

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

Definition at line 527 of file etoile.h.

◆ p_ray_pole

double* Lorene::Etoile::p_ray_pole
mutableprotected

Coordinate radius at $\theta=0$.

Definition at line 536 of file etoile.h.

◆ p_xi_surf

Tbl* Lorene::Etoile::p_xi_surf
mutableprotected

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 548 of file etoile.h.

◆ press

Tenseur Lorene::Etoile::press
protected

Fluid pressure.

Definition at line 464 of file etoile.h.

◆ relativistic

bool Lorene::Etoile::relativistic
protected

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

Definition at line 440 of file etoile.h.

◆ s_euler

Tenseur Lorene::Etoile::s_euler
protected

Trace of the stress tensor in the Eulerian frame.

Definition at line 471 of file etoile.h.

◆ shift

Tenseur Lorene::Etoile::shift
protected

Total shift vector.

Definition at line 515 of file etoile.h.

◆ u_euler

Tenseur Lorene::Etoile::u_euler
protected

Fluid 3-velocity with respect to the Eulerian observer.

Definition at line 477 of file etoile.h.

◆ unsurc2

double Lorene::Etoile::unsurc2
protected

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

Definition at line 445 of file etoile.h.


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