LORENE
Lorene::Etoile_rot Class Reference

Class for isolated rotating stars *** DEPRECATED : use class Star_rot instead ***. More...

#include <etoile.h>

Inheritance diagram for Lorene::Etoile_rot:
Lorene::Etoile Lorene::Et_rot_bifluid Lorene::Et_rot_diff Lorene::Et_rot_mag Lorene::Et_magnetisation

Public Member Functions

 Etoile_rot (Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
 Standard constructor. More...
 
 Etoile_rot (const Etoile_rot &)
 Copy constructor. More...
 
 Etoile_rot (Map &mp_i, const Eos &eos_i, FILE *fich)
 Constructor from a file (see sauve(FILE*) ). More...
 
virtual ~Etoile_rot ()
 Destructor. More...
 
void operator= (const Etoile_rot &)
 Assignment to another Etoile_rot. More...
 
virtual double get_omega_c () const
 Returns the central value of the rotation angular velocity ([f_unit] ) More...
 
const Tenseurget_bbb () const
 Returns the metric factor B. More...
 
const Tenseurget_b_car () const
 Returns the square of the metric factor B. More...
 
const Tenseurget_nphi () const
 Returns the metric coefficient $N^\varphi$. More...
 
const Tenseurget_tnphi () const
 Returns the component $\tilde N^\varphi = N^\varphi r\sin\theta$ of the shift vector. More...
 
const Tenseurget_uuu () const
 Returns the norm of u_euler. More...
 
const Tenseurget_logn () const
 Returns the metric potential $\nu = \ln N$ = logn_auto. More...
 
const Tenseurget_nuf () const
 Returns the part of the Metric potential $\nu = \ln N$ = logn generated by the matter terms. More...
 
const Tenseurget_nuq () const
 Returns the Part of the Metric potential $\nu = \ln N$ = logn generated by the quadratic terms. More...
 
const Tenseurget_dzeta () const
 Returns the Metric potential $\zeta = \ln(AN)$ = beta_auto. More...
 
const Tenseurget_tggg () const
 Returns the Metric potential $\tilde G = (NB-1) r\sin\theta$. More...
 
const Tenseurget_w_shift () const
 Returns the vector $W^i$ used in the decomposition of shift , following Shibata's prescription [Prog. More...
 
const Tenseurget_khi_shift () const
 Returns the scalar $\chi$ used in the decomposition of shift
following Shibata's prescription [Prog. More...
 
const Tenseur_symget_tkij () const
 Returns the tensor ${\tilde K_{ij}}$ related to the extrinsic curvature tensor by ${\tilde K_{ij}} = B^{-2} K_{ij}$. More...
 
const Tenseurget_ak_car () const
 Returns the scalar $A^2 K_{ij} K^{ij}$. More...
 
virtual void sauve (FILE *) const
 Save in a file. More...
 
virtual void display_poly (ostream &) const
 Display in polytropic units. 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...
 
virtual double mass_b () const
 Baryon mass. More...
 
virtual double mass_g () const
 Gravitational mass. More...
 
virtual double angu_mom () const
 Angular momentum. More...
 
virtual double tsw () const
 Ratio T/W. More...
 
virtual double grv2 () const
 Error on the virial identity GRV2. More...
 
virtual double grv3 (ostream *ost=0x0) const
 Error on the virial identity GRV3. More...
 
virtual double r_circ () const
 Circumferential equatorial radius. More...
 
virtual double r_circ_merid () const
 Circumferential meridional radius. More...
 
virtual double area () const
 Surface area. More...
 
virtual double mean_radius () const
 Mean radius. More...
 
virtual double aplat () const
 Flatening r_pole/r_eq. More...
 
virtual double z_eqf () const
 Forward redshift factor at equator. More...
 
virtual double z_eqb () const
 Backward redshift factor at equator. More...
 
virtual double z_pole () const
 Redshift factor at North pole. More...
 
virtual double mom_quad () const
 Quadrupole moment. More...
 
virtual double mom_quad_old () const
 Part of the quadrupole moment. More...
 
virtual double mom_quad_Bo () const
 Part of the quadrupole moment. More...
 
virtual double r_isco (ostream *ost=0x0) const
 Circumferential radius of the innermost stable circular orbit (ISCO). More...
 
virtual double f_isco () const
 Orbital frequency at the innermost stable circular orbit (ISCO). More...
 
virtual double espec_isco () const
 Energy of a particle on the ISCO. More...
 
virtual double lspec_isco () const
 Angular momentum of a particle on the ISCO. More...
 
virtual double f_eccentric (double ecc, double periast, ostream *ost=0x0) const
 Computation of frequency of eccentric orbits. More...
 
virtual double f_eq () const
 Orbital frequency at the equator. More...
 
virtual void hydro_euler ()
 Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame. More...
 
void update_metric ()
 Computes metric coefficients from known potentials. More...
 
void fait_shift ()
 Computes shift from w_shift and khi_shift according to Shibata's prescription [Prog. More...
 
void fait_nphi ()
 Computes tnphi and nphi from the Cartesian components of the shift, stored in shift . More...
 
void extrinsic_curvature ()
 Computes tkij and ak_car from shift , nnn and b_car . More...
 
virtual void equilibrium (double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff, Param *=0x0)
 Computes an equilibrium configuration. 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 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...
 
double ray_eq () const
 Coordinate radius at $\phi=0$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_eq (int kk) const
 Coordinate radius at $\phi=2k\pi/np$, $\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...
 
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...
 

Static Public Member Functions

static double lambda_grv2 (const Cmp &sou_m, const Cmp &sou_q)
 Computes the coefficient $\lambda$ which ensures that the GRV2 virial identity is satisfied. 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...
 
virtual void partial_display (ostream &) const
 Printing of some informations, excluding all global quantities. More...
 

Protected Attributes

double omega
 Rotation angular velocity ([f_unit] ) More...
 
Tenseur bbb
 Metric factor B. More...
 
Tenseur b_car
 Square of the metric factor B. More...
 
Tenseur nphi
 Metric coefficient $N^\varphi$. More...
 
Tenseur tnphi
 Component $\tilde N^\varphi = N^\varphi r\sin\theta$ of the shift vector. More...
 
Tenseur uuu
 Norm of u_euler. More...
 
Tenseurlogn
 Metric potential $\nu = \ln N$ = logn_auto. More...
 
Tenseur nuf
 Part of the Metric potential $\nu = \ln N$ = logn generated by the matter terms. More...
 
Tenseur nuq
 Part of the Metric potential $\nu = \ln N$ = logn generated by the quadratic terms. More...
 
Tenseurdzeta
 Metric potential $\zeta = \ln(AN)$ = beta_auto. More...
 
Tenseur tggg
 Metric potential $\tilde G = (NB-1) r\sin\theta$. More...
 
Tenseur w_shift
 Vector $W^i$ used in the decomposition of shift , following Shibata's prescription [Prog. More...
 
Tenseur khi_shift
 Scalar $\chi$ used in the decomposition of shift , following Shibata's prescription [Prog. More...
 
Tenseur_sym tkij
 Tensor ${\tilde K_{ij}}$ related to the extrinsic curvature tensor by ${\tilde K_{ij}} = B^{-2} K_{ij}$. More...
 
Tenseur ak_car
 Scalar $A^2 K_{ij} K^{ij}$. More...
 
Cmp ssjm1_nuf
 Effective source at the previous step for the resolution of the Poisson equation for nuf by means of Map_et::poisson . More...
 
Cmp ssjm1_nuq
 Effective source at the previous step for the resolution of the Poisson equation for nuq by means of Map_et::poisson . More...
 
Cmp ssjm1_dzeta
 Effective source at the previous step for the resolution of the Poisson equation for dzeta . More...
 
Cmp ssjm1_tggg
 Effective source at the previous step for the resolution of the Poisson equation for tggg . More...
 
Cmp ssjm1_khi
 Effective source at the previous step for the resolution of the Poisson equation for the scalar $\chi$ by means of Map_et::poisson . More...
 
Tenseur ssjm1_wshift
 Effective source at the previous step for the resolution of the vector Poisson equation for $W^i$. More...
 
double * p_angu_mom
 Angular momentum. More...
 
double * p_tsw
 Ratio T/W. More...
 
double * p_grv2
 Error on the virial identity GRV2. More...
 
double * p_grv3
 Error on the virial identity GRV3. More...
 
double * p_r_circ
 Circumferential equatorial radius. More...
 
double * p_r_circ_merid
 Circumferential meridional radius. More...
 
double * p_area
 Surface area. More...
 
double * p_aplat
 Flatening r_pole/r_eq. More...
 
double * p_z_eqf
 Forward redshift factor at equator. More...
 
double * p_z_eqb
 Backward redshift factor at equator. More...
 
double * p_z_pole
 Redshift factor at North pole. More...
 
double * p_mom_quad
 Quadrupole moment. More...
 
double * p_mom_quad_old
 Part of the quadrupole moment. More...
 
double * p_mom_quad_Bo
 Part of the quadrupole moment. More...
 
double * p_r_isco
 Circumferential radius of the ISCO. More...
 
double * p_f_isco
 Orbital frequency of the ISCO. More...
 
double * p_espec_isco
 Specific energy of a particle on the ISCO. More...
 
double * p_lspec_isco
 Specific angular momentum of a particle on the ISCO. More...
 
double * p_f_eq
 Orbital frequency at the equator. More...
 
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...
 

Detailed Description

Class for isolated rotating stars *** DEPRECATED : use class Star_rot instead ***.

()

The metric is

\[ ds^2 = - N^2 dt^2 + A^2 (dr^2 + r^2 d\theta^2) + B^2 r^2 \sin^2\theta (d\varphi - N^\varphi dt)^2 \]

Definition at line 1499 of file etoile.h.

Constructor & Destructor Documentation

◆ Etoile_rot() [1/3]

Lorene::Etoile_rot::Etoile_rot ( 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_rot.C.

References ak_car, b_car, bbb, khi_shift, nphi, nuf, nuq, omega, Lorene::Tenseur::set(), set_der_0x0(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_etat_zero(), Lorene::Tenseur::set_std_base(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, tggg, tkij, tnphi, uuu, and w_shift.

◆ Etoile_rot() [2/3]

Lorene::Etoile_rot::Etoile_rot ( const Etoile_rot et)

Copy constructor.

Definition at line 216 of file etoile_rot.C.

References omega, and set_der_0x0().

◆ Etoile_rot() [3/3]

Lorene::Etoile_rot::Etoile_rot ( 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 248 of file etoile_rot.C.

References b_car, bbb, fait_nphi(), fait_shift(), Lorene::fread_be(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), khi_shift, Lorene::Etoile::mp, nuf, nuq, omega, set_der_0x0(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, tggg, uuu, and w_shift.

◆ ~Etoile_rot()

Lorene::Etoile_rot::~Etoile_rot ( )
virtual

Destructor.

Definition at line 334 of file etoile_rot.C.

References del_deriv().

Member Function Documentation

◆ angu_mom()

◆ aplat()

double Lorene::Etoile_rot::aplat ( ) const
virtual

Flatening r_pole/r_eq.

Definition at line 554 of file et_rot_global.C.

References p_aplat, Lorene::Etoile::ray_eq(), and Lorene::Etoile::ray_pole().

◆ area()

◆ del_deriv()

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

◆ del_hydro_euler()

void Lorene::Etoile_rot::del_hydro_euler ( )
protectedvirtual

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

Reimplemented from Lorene::Etoile.

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

Definition at line 400 of file etoile_rot.C.

References del_deriv(), and Lorene::Etoile::del_hydro_euler().

◆ display_poly()

void Lorene::Etoile_rot::display_poly ( ostream &  ost) const
virtual

◆ equation_of_state()

◆ equil_spher_falloff()

◆ equil_spher_regular()

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

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 118 of file et_equil_spher_regu.C.

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

◆ equilibrium()

void Lorene::Etoile_rot::equilibrium ( double  ent_c,
double  omega0,
double  fact_omega,
int  nzadapt,
const Tbl ent_limit,
const Itbl icontrol,
const Tbl control,
double  mbar_wanted,
double  aexp_mass,
Tbl diff,
Param = 0x0 
)
virtual

Computes an equilibrium configuration.

Parameters
ent_c[input] Central enthalpy
omega0[input] Requested angular velocity (if fact_omega=1. )
fact_omega[input] 1.01 = search for the Keplerian frequency,
  1. = otherwise.
nzadapt[input] Number of (inner) domains where the mapping adaptation to an iso-enthalpy surface should be performed
ent_limit[input] 1-D Tbl of dimension nzet which defines the enthalpy at the outer boundary of each domain
icontrol[input] Set of integer parameters (stored as a 1-D Itbl of size 8) to control the iteration:
  • icontrol(0) = mer_max : maximum number of steps
  • icontrol(1) = mer_rot : step at which the rotation is switched on
  • icontrol(2) = mer_change_omega : step at which the rotation velocity is changed to reach the final one
  • icontrol(3) = mer_fix_omega : step at which the final rotation velocity must have been reached
  • icontrol(4) = mer_mass : the absolute value of mer_mass is the step from which the baryon mass is forced to converge, by varying the central enthalpy (mer_mass>0 ) or the angular velocity (mer_mass<0 )
  • icontrol(5) = mermax_poisson : maximum number of steps in Map_et::poisson
  • icontrol(6) = mer_triax : step at which the 3-D perturbation is switched on
  • icontrol(7) = delta_mer_kep : number of steps after mer_fix_omega when omega starts to be increased by fact_omega to search for the Keplerian velocity
control[input] Set of parameters (stored as a 1-D Tbl of size 7) to control the iteration:
  • control(0) = precis : threshold on the enthalpy relative change for ending the computation
  • control(1) = omega_ini : initial angular velocity, switched on only if mer_rot<0 , otherwise 0 is used
  • control(2) = relax : relaxation factor in the main iteration
  • control(3) = relax_poisson : relaxation factor in Map_et::poisson
  • control(4) = thres_adapt : threshold on dH/dr for freezing the adaptation of the mapping
  • control(5) = ampli_triax : relative amplitude of the 3-D perturbation
  • control(6) = precis_adapt : precision for Map_et::adapt
mbar_wanted[input] Requested baryon mass (effective only if mer_mass > mer_max )
aexp_mass[input] Exponent for the increase factor of the central enthalpy to converge to the requested baryon mass
diff[output] 1-D Tbl of size 7 for the storage of some error indicators :
  • diff(0) : Relative change in the enthalpy field between two successive steps
  • diff(1) : Relative error in the resolution of the Poisson equation for nuf
  • diff(2) : Relative error in the resolution of the Poisson equation for nuq
  • diff(3) : Relative error in the resolution of the Poisson equation for dzeta
  • diff(4) : Relative error in the resolution of the Poisson equation for tggg
  • diff(5) : Relative error in the resolution of the equation for shift (x comp.)
  • diff(6) : Relative error in the resolution of the equation for shift (y comp.)

Reimplemented in Lorene::Et_rot_diff.

Definition at line 150 of file et_rot_equilibrium.C.

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

◆ equilibrium_spher()

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

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 90 of file etoile_equil_spher.C.

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

◆ espec_isco()

double Lorene::Etoile_rot::espec_isco ( ) const
virtual

Energy of a particle on the ISCO.

Definition at line 304 of file et_rot_isco.C.

References p_espec_isco, and r_isco().

◆ extrinsic_curvature()

void Lorene::Etoile_rot::extrinsic_curvature ( )

Computes tkij and ak_car from shift , nnn and b_car .

Definition at line 60 of file et_rot_extr_curv.C.

References Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Etoile::mp, nphi, Lorene::Tenseur::set_etat_zero(), and tkij.

◆ f_eccentric()

double Lorene::Etoile_rot::f_eccentric ( double  ecc,
double  periast,
ostream *  ost = 0x0 
) const
virtual

Computation of frequency of eccentric orbits.

Parameters
ecceccentricity of the orbit
periasrtperiastron of the orbit
ostoutput stream to give details of the computation; if set to 0x0 [default value], no details will be given.
Returns
orbital frequency

Definition at line 81 of file et_rot_f_eccentric.C.

References Lorene::Param::add_cmp(), Lorene::Param::add_int(), Lorene::Cmp::annule(), bbb, Lorene::Cmp::dsdr(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Etoile::mp, Lorene::Etoile::nnn, nphi, Lorene::Etoile::nzet, p_f_isco, p_r_isco, Lorene::Map::r, Lorene::Etoile::ray_eq(), Lorene::sqrt(), Lorene::Cmp::std_base_scal(), Lorene::Cmp::va, Lorene::Valeur::val_point(), and Lorene::Map::val_r().

◆ f_eq()

double Lorene::Etoile_rot::f_eq ( ) const
virtual

Orbital frequency at the equator.

Definition at line 322 of file et_rot_isco.C.

References p_f_eq, and r_isco().

◆ f_isco()

double Lorene::Etoile_rot::f_isco ( ) const
virtual

Orbital frequency at the innermost stable circular orbit (ISCO).

Definition at line 270 of file et_rot_isco.C.

References p_f_isco, and r_isco().

◆ fait_nphi()

void Lorene::Etoile_rot::fait_nphi ( )

Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .

Definition at line 803 of file etoile_rot.C.

References Lorene::Tenseur::get_etat(), and Lorene::Etoile::shift.

◆ fait_shift()

void Lorene::Etoile_rot::fait_shift ( )

Computes shift from w_shift and khi_shift according to Shibata's prescription [Prog.

Theor. Phys. 101 , 1199 (1999)] :

\[ N^i = {7\over 8} W^i - {1\over 8} \left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

Definition at line 766 of file etoile_rot.C.

References Lorene::Tenseur::get_etat(), Lorene::Tenseur::gradient(), and khi_shift.

◆ get_a_car()

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

Returns the total conformal factor $A^2$.

Definition at line 736 of file etoile.h.

References Lorene::Etoile::a_car.

◆ get_ak_car()

const Tenseur& Lorene::Etoile_rot::get_ak_car ( ) const
inline

Returns the scalar $A^2 K_{ij} K^{ij}$.

For axisymmetric stars, this quantity is related to the derivatives of $N^\varphi$ by

\[ A^2 K_{ij} K^{ij} = {B^2 \over 2 N^2} \, r^2\sin^2\theta \, \left[ \left( {\partial N^\varphi \over \partial r} \right) ^2 + {1\over r^2} \left( {\partial N^\varphi \over \partial \theta} \right) ^2 \right] \ . \]

In particular it is related to the quantities $k_1$ and $k_2$ introduced by Eqs.~(3.7) and (3.8) of Bonazzola et al. Astron. Astrophys. 278 , 421 (1993) by

\[ A^2 K_{ij} K^{ij} = 2 A^2 (k_1^2 + k_2^2) \ . \]

Definition at line 1803 of file etoile.h.

References ak_car.

◆ get_b_car()

const Tenseur& Lorene::Etoile_rot::get_b_car ( ) const
inline

Returns the square of the metric factor B.

Definition at line 1719 of file etoile.h.

References b_car.

◆ get_bbb()

const Tenseur& Lorene::Etoile_rot::get_bbb ( ) const
inline

Returns the metric factor B.

Definition at line 1716 of file etoile.h.

References bbb.

◆ get_beta_auto()

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

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

Definition at line 727 of file etoile.h.

References Lorene::Etoile::beta_auto.

◆ get_d_logn_auto_div()

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

Returns the gradient of logn_auto_div.

Definition at line 722 of file etoile.h.

References Lorene::Etoile::d_logn_auto_div.

◆ get_dzeta()

const Tenseur& Lorene::Etoile_rot::get_dzeta ( ) const
inline

Returns the Metric potential $\zeta = \ln(AN)$ = beta_auto.

Definition at line 1746 of file etoile.h.

References dzeta.

◆ get_ener()

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

Returns the proper total energy density.

Definition at line 682 of file etoile.h.

References Lorene::Etoile::ener.

◆ get_ener_euler()

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

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

Definition at line 688 of file etoile.h.

References Lorene::Etoile::ener_euler.

◆ get_ent()

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

Returns the enthalpy field.

Definition at line 676 of file etoile.h.

References Lorene::Etoile::ent.

◆ get_eos()

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

Returns the equation of state.

Definition at line 673 of file etoile.h.

References Lorene::Etoile::eos.

◆ get_gam_euler()

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

Returns the Lorentz factor between the fluid and Eulerian observers.

Definition at line 694 of file etoile.h.

References Lorene::Etoile::gam_euler.

◆ get_khi_shift()

const Tenseur& Lorene::Etoile_rot::get_khi_shift ( ) const
inline

Returns the scalar $\chi$ used in the decomposition of shift
following Shibata's prescription [Prog.

Theor. Phys. 101 , 1199 (1999)] :

\[ N^i = {7\over 8} W^i - {1\over 8} \left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

NB: w_shift contains the components of $W^i$ with respect to the Cartesian triad associated with the mapping mp .

Definition at line 1777 of file etoile.h.

References khi_shift.

◆ get_logn()

const Tenseur& Lorene::Etoile_rot::get_logn ( ) const
inline

Returns the metric potential $\nu = \ln N$ = logn_auto.

Definition at line 1733 of file etoile.h.

References logn.

◆ get_logn_auto()

const Tenseur& Lorene::Etoile::get_logn_auto ( ) const
inlineinherited

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 Lorene::Etoile::logn_auto.

◆ get_logn_auto_div()

const Tenseur& Lorene::Etoile::get_logn_auto_div ( ) const
inlineinherited

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 Lorene::Etoile::logn_auto_div.

◆ get_logn_auto_regu()

const Tenseur& Lorene::Etoile::get_logn_auto_regu ( ) const
inlineinherited

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 Lorene::Etoile::logn_auto_regu.

◆ get_mp()

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

Returns the mapping.

Definition at line 662 of file etoile.h.

References Lorene::Etoile::mp.

◆ get_nbar()

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

Returns the proper baryon density.

Definition at line 679 of file etoile.h.

References Lorene::Etoile::nbar.

◆ get_nnn()

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

Returns the total lapse function N.

Definition at line 730 of file etoile.h.

References Lorene::Etoile::nnn.

◆ get_nphi()

const Tenseur& Lorene::Etoile_rot::get_nphi ( ) const
inline

Returns the metric coefficient $N^\varphi$.

Definition at line 1722 of file etoile.h.

References nphi.

◆ get_nuf()

const Tenseur& Lorene::Etoile_rot::get_nuf ( ) const
inline

Returns the part of the Metric potential $\nu = \ln N$ = logn generated by the matter terms.

Definition at line 1738 of file etoile.h.

References nuf.

◆ get_nuq()

const Tenseur& Lorene::Etoile_rot::get_nuq ( ) const
inline

Returns the Part of the Metric potential $\nu = \ln N$ = logn generated by the quadratic terms.

Definition at line 1743 of file etoile.h.

References nuq.

◆ get_nzet()

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

Returns the number of domains occupied by the star.

Definition at line 665 of file etoile.h.

References Lorene::Etoile::nzet.

◆ get_omega_c()

double Lorene::Etoile_rot::get_omega_c ( ) const
virtual

Returns the central value of the rotation angular velocity ([f_unit] )

Reimplemented in Lorene::Et_rot_diff.

Definition at line 698 of file etoile_rot.C.

References omega.

◆ get_press()

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

Returns the fluid pressure.

Definition at line 685 of file etoile.h.

References Lorene::Etoile::press.

◆ get_s_euler()

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

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

Definition at line 691 of file etoile.h.

References Lorene::Etoile::s_euler.

◆ get_shift()

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

Returns the total shift vector $N^i$.

Definition at line 733 of file etoile.h.

References Lorene::Etoile::shift.

◆ get_tggg()

const Tenseur& Lorene::Etoile_rot::get_tggg ( ) const
inline

Returns the Metric potential $\tilde G = (NB-1) r\sin\theta$.

Definition at line 1749 of file etoile.h.

References tggg.

◆ get_tkij()

const Tenseur_sym& Lorene::Etoile_rot::get_tkij ( ) const
inline

Returns the tensor ${\tilde K_{ij}}$ related to the extrinsic curvature tensor by ${\tilde K_{ij}} = B^{-2} K_{ij}$.

tkij contains the Cartesian components of ${\tilde K_{ij}}$.

Definition at line 1784 of file etoile.h.

References tkij.

◆ get_tnphi()

const Tenseur& Lorene::Etoile_rot::get_tnphi ( ) const
inline

Returns the component $\tilde N^\varphi = N^\varphi r\sin\theta$ of the shift vector.

Definition at line 1727 of file etoile.h.

References tnphi.

◆ get_u_euler()

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

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

Definition at line 697 of file etoile.h.

References Lorene::Etoile::u_euler.

◆ get_uuu()

const Tenseur& Lorene::Etoile_rot::get_uuu ( ) const
inline

Returns the norm of u_euler.

Definition at line 1730 of file etoile.h.

References uuu.

◆ get_w_shift()

const Tenseur& Lorene::Etoile_rot::get_w_shift ( ) const
inline

Returns the vector $W^i$ used in the decomposition of shift , following Shibata's prescription [Prog.

Theor. Phys. 101 , 1199 (1999)] :

\[ N^i = {7\over 8} W^i - {1\over 8} \left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

NB: w_shift contains the components of $W^i$ with respect to the Cartesian triad associated with the mapping mp .

Definition at line 1763 of file etoile.h.

References w_shift.

◆ grv2()

double Lorene::Etoile_rot::grv2 ( ) const
virtual

◆ grv3()

double Lorene::Etoile_rot::grv3 ( ostream *  ost = 0x0) const
virtual

Error on the virial identity GRV3.

The error is computed as the integral defined by Eq. (43) of [Gourgoulhon and Bonazzola, Class. Quantum Grav. 11, 443 (1994)] divided by the integral of the matter terms.

Parameters
ostoutput stream to give details of the computation; if set to 0x0 [default value], no details will be given.

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

Definition at line 294 of file et_rot_global.C.

References ak_car, bbb, dzeta, Lorene::flat_scalar_prod(), Lorene::Cmp::get_etat(), Lorene::Tenseur::gradient_spher(), Lorene::log(), logn, Lorene::Etoile::mp, p_grv3, Lorene::Etoile::relativistic, Lorene::Tenseur::set_std_base(), and Lorene::Cmp::srdsdt().

◆ hydro_euler()

void Lorene::Etoile_rot::hydro_euler ( )
virtual

Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame.

The calculation is performed starting from the quantities ent , ener , press , and a_car ,
which are supposed to be up to date.
From these, the following fields are updated: gam_euler , u_euler , ener_euler , s_euler .

Reimplemented from Lorene::Etoile.

Reimplemented in Lorene::Et_rot_bifluid, and Lorene::Et_rot_diff.

Definition at line 86 of file et_rot_hydro.C.

References Lorene::Tenseur::annule(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Etoile::mp, Lorene::Etoile::nnn, omega, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), Lorene::Tenseur::set_triad(), Lorene::Etoile::shift, Lorene::Etoile::u_euler, Lorene::Map::x, and Lorene::Map::y.

◆ is_relativistic()

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

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

Definition at line 670 of file etoile.h.

References Lorene::Etoile::relativistic.

◆ l_surf()

const Itbl & Lorene::Etoile_rot::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 from Lorene::Etoile.

Reimplemented in Lorene::Et_rot_bifluid.

Definition at line 104 of file et_rot_global.C.

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

◆ lambda_grv2()

double Lorene::Etoile_rot::lambda_grv2 ( const Cmp sou_m,
const Cmp sou_q 
)
static

Computes the coefficient $\lambda$ which ensures that the GRV2 virial identity is satisfied.

$\lambda$ is the coefficient by which one must multiply the quadratic source term $\sigma_q$ of the 2-D Poisson equation

\[ \Delta_2 u = \sigma_m + \sigma_q \]

in order that the total source does not contain any monopolar term, i.e. in order that

\[ \int_0^{2\pi} \int_0^{+\infty} \sigma(r, \theta) \, r \, dr \, d\theta = 0 \ , \]

where $\sigma = \sigma_m + \sigma_q$. $\lambda$ is computed according to the formula

\[ \lambda = - { \int_0^{2\pi} \int_0^{+\infty} \sigma_m(r, \theta) \, r \, dr \, d\theta \over \int_0^{2\pi} \int_0^{+\infty} \sigma_q(r, \theta) \, r \, dr \, d\theta } \ . \]

Then, by construction, the new source $\sigma' = \sigma_m + \lambda \sigma_q$ has a vanishing monopolar term.

Parameters
sou_m[input] matter source term $\sigma_m$
sou_q[input] quadratic source term $\sigma_q$
Returns
value of $\lambda$

Definition at line 82 of file et_rot_lambda_grv2.C.

References Lorene::Cmp::check_dzpuis(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::val_r().

◆ lspec_isco()

double Lorene::Etoile_rot::lspec_isco ( ) const
virtual

Angular momentum of a particle on the ISCO.

Definition at line 287 of file et_rot_isco.C.

References p_lspec_isco, and r_isco().

◆ mass_b()

◆ mass_g()

◆ mean_radius()

double Lorene::Etoile_rot::mean_radius ( ) const
virtual

Mean radius.

Definition at line 530 of file et_rot_global.C.

References area(), and Lorene::sqrt().

◆ mom_quad()

double Lorene::Etoile_rot::mom_quad ( ) const
virtual

Quadrupole moment.

The quadrupole moment Q is defined according to Eq. (11) of [Pappas and Apostolatos, Physical Review Letters 108, 231104 (2012)]. This is a corrected version of the quadrupole moment defined by [Salgado, Bonazzola, Gourgoulhon and Haensel, Astron. Astrophys. 291 , 155 (1994)]. Following this definition, $Q = {\bar Q } - 4/3 (1/4 + b) M^3 $, where ${\bar Q }$ is defined as the negative of the (wrong) quadrupole moment defined in Eq. (7) of [Salgado, Bonazzola, Gourgoulhon and Haensel, Astron. Astrophys. 291 , 155 (1994)], b is defined by Eq. (3.37) of [Friedman and Stergioulas, Rotating Relativistic Stars, Cambridge Monograph on mathematical physics] and M is the gravitational mass of the star.

Reimplemented in Lorene::Et_rot_bifluid.

Definition at line 658 of file et_rot_global.C.

References mass_g(), mom_quad_Bo(), mom_quad_old(), p_mom_quad, Lorene::pow(), and Lorene::Etoile::relativistic.

◆ mom_quad_Bo()

double Lorene::Etoile_rot::mom_quad_Bo ( ) const
virtual

Part of the quadrupole moment.

$B_o$ is defined as $bM^2$, where b is given by Eq. (3.37) of [Friedman and Stergioulas, Rotating Relativistic Stars, Cambridge Monograph on mathematical physics] and M is the the gravitational mass of the star.

Reimplemented in Lorene::Et_magnetisation, and Lorene::Et_rot_bifluid.

Definition at line 676 of file et_rot_global.C.

References Lorene::Etoile::a_car, bbb, Lorene::Cmp::integrale(), Lorene::Etoile::mp, Lorene::Cmp::mult_rsint(), Lorene::Etoile::nnn, p_mom_quad_Bo, Lorene::Etoile::press, and Lorene::Cmp::std_base_scal().

◆ mom_quad_old()

double Lorene::Etoile_rot::mom_quad_old ( ) const
virtual

Part of the quadrupole moment.

This term ${\bar Q }$ is defined by Laarakkers and Poisson, Astrophys. J. 512 , 282 (1999). Note that ${\bar Q }$ is the negative of the (wrong) quadrupole moment defined in Eq. (7) of [Salgado, Bonazzola, Gourgoulhon and Haensel, Astron. Astrophys. 291 , 155 (1994)].

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

Definition at line 699 of file et_rot_global.C.

References Lorene::Etoile::a_car, ak_car, bbb, Lorene::Cmp::check_dzpuis(), Lorene::Etoile::ener_euler, Lorene::flat_scalar_prod(), Lorene::Cmp::get_etat(), Lorene::Tenseur::gradient_spher(), Lorene::Cmp::inc2_dzpuis(), Lorene::log(), logn, Lorene::Etoile::mp, Lorene::Cmp::mult_r(), Lorene::Etoile::nbar, p_mom_quad_old, Lorene::Etoile::relativistic, Lorene::Etoile::s_euler, Lorene::Tenseur::set(), and Lorene::Tenseur::set_std_base().

◆ operator=()

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

◆ operator>>()

◆ partial_display()

void Lorene::Etoile_rot::partial_display ( ostream &  ost) const
protectedvirtual

◆ r_circ()

double Lorene::Etoile_rot::r_circ ( ) const
virtual

Circumferential equatorial radius.

Definition at line 399 of file et_rot_global.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_type_t(), Lorene::Etoile::mp, Lorene::Etoile::nzet, and p_r_circ.

◆ r_circ_merid()

◆ r_isco()

double Lorene::Etoile_rot::r_isco ( ostream *  ost = 0x0) const
virtual

Circumferential radius of the innermost stable circular orbit (ISCO).

Parameters
ostoutput stream to give details of the computation; if set to 0x0 [default value], no details will be given.

Definition at line 87 of file et_rot_isco.C.

References Lorene::Param::add_cmp(), Lorene::Param::add_int(), Lorene::Cmp::annule(), bbb, Lorene::Cmp::dsdr(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Etoile::mp, Lorene::Etoile::nnn, nphi, Lorene::Etoile::nzet, p_espec_isco, p_f_eq, p_f_isco, p_lspec_isco, p_r_isco, Lorene::Map::r, Lorene::Etoile::ray_eq(), Lorene::sqrt(), Lorene::Cmp::std_base_scal(), Lorene::Cmp::va, Lorene::Valeur::val_point(), Lorene::Map::val_r(), and Lorene::zerosec().

◆ ray_eq() [1/2]

double Lorene::Etoile::ray_eq ( ) const
inherited

◆ ray_eq() [2/2]

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

◆ ray_eq_3pis2()

double Lorene::Etoile::ray_eq_3pis2 ( ) const
inherited

◆ ray_eq_pi()

double Lorene::Etoile::ray_eq_pi ( ) const
inherited

◆ ray_eq_pis2()

double Lorene::Etoile::ray_eq_pis2 ( ) const
inherited

◆ ray_pole()

double Lorene::Etoile::ray_pole ( ) const
inherited

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(), Lorene::Etoile::mp, and Lorene::Etoile::p_ray_pole.

◆ sauve()

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

◆ set_der_0x0()

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

Sets to 0x0 all the pointers on derived quantities.

Reimplemented from Lorene::Etoile.

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

Definition at line 374 of file etoile_rot.C.

References p_angu_mom, p_aplat, p_area, p_espec_isco, p_f_eq, p_f_isco, p_grv2, p_grv3, p_lspec_isco, p_mom_quad, p_mom_quad_Bo, p_mom_quad_old, p_r_circ, p_r_circ_merid, p_r_isco, p_tsw, p_z_eqb, p_z_eqf, p_z_pole, and Lorene::Etoile::set_der_0x0().

◆ set_enthalpy()

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

Assignment of the enthalpy field.

Definition at line 468 of file etoile.C.

References Lorene::Etoile::del_deriv(), Lorene::Etoile::ent, and Lorene::Etoile::equation_of_state().

◆ set_mp()

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

Read/write of the mapping.

Definition at line 604 of file etoile.h.

References Lorene::Etoile::mp.

◆ tsw()

◆ update_metric()

void Lorene::Etoile_rot::update_metric ( )

Computes metric coefficients from known potentials.

The calculation is performed starting from the quantities logn , dzeta , tggg and shift , which are supposed to be up to date.
From these, the following fields are updated: nnn , a_car , bbb and b_car .

Definition at line 72 of file et_rot_upmetr.C.

References Lorene::Etoile::a_car, b_car, bbb, del_deriv(), dzeta, Lorene::exp(), extrinsic_curvature(), logn, Lorene::Etoile::nnn, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), tggg, and Lorene::Etoile::unsurc2.

◆ xi_surf()

const Tbl & Lorene::Etoile::xi_surf ( ) const
inherited

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 Lorene::Etoile::l_surf(), Lorene::Etoile::p_l_surf, and Lorene::Etoile::p_xi_surf.

◆ z_eqb()

double Lorene::Etoile_rot::z_eqb ( ) const
virtual

Backward redshift factor at equator.

Definition at line 603 of file et_rot_global.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_type_t(), Lorene::Etoile::mp, Lorene::Etoile::nzet, and p_z_eqb.

◆ z_eqf()

double Lorene::Etoile_rot::z_eqf ( ) const
virtual

Forward redshift factor at equator.

Definition at line 571 of file et_rot_global.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_type_t(), Lorene::Etoile::mp, Lorene::Etoile::nzet, and p_z_eqf.

◆ z_pole()

double Lorene::Etoile_rot::z_pole ( ) const
virtual

Redshift factor at North pole.

Definition at line 638 of file et_rot_global.C.

References Lorene::Etoile::nnn, p_z_pole, and Lorene::Etoile::ray_pole().

Member Data Documentation

◆ a_car

Tenseur Lorene::Etoile::a_car
protectedinherited

Total conformal factor $A^2$.

Definition at line 518 of file etoile.h.

◆ ak_car

Tenseur Lorene::Etoile_rot::ak_car
protected

Scalar $A^2 K_{ij} K^{ij}$.

For axisymmetric stars, this quantity is related to the derivatives of $N^\varphi$ by

\[ A^2 K_{ij} K^{ij} = {B^2 \over 2 N^2} \, r^2\sin^2\theta \, \left[ \left( {\partial N^\varphi \over \partial r} \right) ^2 + {1\over r^2} \left( {\partial N^\varphi \over \partial \theta} \right) ^2 \right] \ . \]

In particular it is related to the quantities $k_1$ and $k_2$ introduced by Eqs.~(3.7) and (3.8) of Bonazzola et al. Astron. Astrophys. 278 , 421 (1993) by

\[ A^2 K_{ij} K^{ij} = 2 A^2 (k_1^2 + k_2^2) \ . \]

Definition at line 1589 of file etoile.h.

◆ b_car

Tenseur Lorene::Etoile_rot::b_car
protected

Square of the metric factor B.

Definition at line 1510 of file etoile.h.

◆ bbb

Tenseur Lorene::Etoile_rot::bbb
protected

Metric factor B.

Definition at line 1507 of file etoile.h.

◆ beta_auto

Tenseur Lorene::Etoile::beta_auto
protectedinherited

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
protectedinherited

Gradient of logn_auto_div (if k_div!=0 )

Definition at line 504 of file etoile.h.

◆ dzeta

Tenseur& Lorene::Etoile_rot::dzeta
protected

Metric potential $\zeta = \ln(AN)$ = beta_auto.

Definition at line 1537 of file etoile.h.

◆ ener

Tenseur Lorene::Etoile::ener
protectedinherited

Total energy density in the fluid frame.

Definition at line 463 of file etoile.h.

◆ ener_euler

Tenseur Lorene::Etoile::ener_euler
protectedinherited

Total energy density in the Eulerian frame.

Definition at line 468 of file etoile.h.

◆ ent

Tenseur Lorene::Etoile::ent
protectedinherited

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

Definition at line 460 of file etoile.h.

◆ eos

const Eos& Lorene::Etoile::eos
protectedinherited

Equation of state of the stellar matter.

Definition at line 454 of file etoile.h.

◆ gam_euler

Tenseur Lorene::Etoile::gam_euler
protectedinherited

Lorentz factor between the fluid and Eulerian observers.

Definition at line 474 of file etoile.h.

◆ k_div

int Lorene::Etoile::k_div
protectedinherited

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.

◆ khi_shift

Tenseur Lorene::Etoile_rot::khi_shift
protected

Scalar $\chi$ used in the decomposition of shift , following Shibata's prescription [Prog.

Theor. Phys. 101 , 1199 (1999)] :

\[ N^i = {7\over 8} W^i - {1\over 8} \left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

Definition at line 1563 of file etoile.h.

◆ logn

Tenseur& Lorene::Etoile_rot::logn
protected

Metric potential $\nu = \ln N$ = logn_auto.

Definition at line 1524 of file etoile.h.

◆ logn_auto

Tenseur Lorene::Etoile::logn_auto
protectedinherited

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
protectedinherited

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
protectedinherited

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
protectedinherited

Mapping associated with the star.

Definition at line 432 of file etoile.h.

◆ nbar

Tenseur Lorene::Etoile::nbar
protectedinherited

Baryon density in the fluid frame.

Definition at line 462 of file etoile.h.

◆ nnn

Tenseur Lorene::Etoile::nnn
protectedinherited

Total lapse function.

Definition at line 512 of file etoile.h.

◆ nphi

Tenseur Lorene::Etoile_rot::nphi
protected

Metric coefficient $N^\varphi$.

Definition at line 1513 of file etoile.h.

◆ nuf

Tenseur Lorene::Etoile_rot::nuf
protected

Part of the Metric potential $\nu = \ln N$ = logn generated by the matter terms.

Definition at line 1529 of file etoile.h.

◆ nuq

Tenseur Lorene::Etoile_rot::nuq
protected

Part of the Metric potential $\nu = \ln N$ = logn generated by the quadratic terms.

Definition at line 1534 of file etoile.h.

◆ nzet

int Lorene::Etoile::nzet
protectedinherited

Number of domains of *mp occupied by the star.

Definition at line 435 of file etoile.h.

◆ omega

double Lorene::Etoile_rot::omega
protected

Rotation angular velocity ([f_unit] )

Definition at line 1504 of file etoile.h.

◆ p_angu_mom

double* Lorene::Etoile_rot::p_angu_mom
mutableprotected

Angular momentum.

Definition at line 1634 of file etoile.h.

◆ p_aplat

double* Lorene::Etoile_rot::p_aplat
mutableprotected

Flatening r_pole/r_eq.

Definition at line 1641 of file etoile.h.

◆ p_area

double* Lorene::Etoile_rot::p_area
mutableprotected

Surface area.

Definition at line 1640 of file etoile.h.

◆ p_espec_isco

double* Lorene::Etoile_rot::p_espec_isco
mutableprotected

Specific energy of a particle on the ISCO.

Definition at line 1651 of file etoile.h.

◆ p_f_eq

double* Lorene::Etoile_rot::p_f_eq
mutableprotected

Orbital frequency at the equator.

Definition at line 1654 of file etoile.h.

◆ p_f_isco

double* Lorene::Etoile_rot::p_f_isco
mutableprotected

Orbital frequency of the ISCO.

Definition at line 1649 of file etoile.h.

◆ p_grv2

double* Lorene::Etoile_rot::p_grv2
mutableprotected

Error on the virial identity GRV2.

Definition at line 1636 of file etoile.h.

◆ p_grv3

double* Lorene::Etoile_rot::p_grv3
mutableprotected

Error on the virial identity GRV3.

Definition at line 1637 of file etoile.h.

◆ p_l_surf

Itbl* Lorene::Etoile::p_l_surf
mutableprotectedinherited

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_lspec_isco

double* Lorene::Etoile_rot::p_lspec_isco
mutableprotected

Specific angular momentum of a particle on the ISCO.

Definition at line 1653 of file etoile.h.

◆ p_mass_b

double* Lorene::Etoile::p_mass_b
mutableprotectedinherited

Baryon mass.

Definition at line 550 of file etoile.h.

◆ p_mass_g

double* Lorene::Etoile::p_mass_g
mutableprotectedinherited

Gravitational mass.

Definition at line 551 of file etoile.h.

◆ p_mom_quad

double* Lorene::Etoile_rot::p_mom_quad
mutableprotected

Quadrupole moment.

Definition at line 1645 of file etoile.h.

◆ p_mom_quad_Bo

double* Lorene::Etoile_rot::p_mom_quad_Bo
mutableprotected

Part of the quadrupole moment.

Definition at line 1647 of file etoile.h.

◆ p_mom_quad_old

double* Lorene::Etoile_rot::p_mom_quad_old
mutableprotected

Part of the quadrupole moment.

Definition at line 1646 of file etoile.h.

◆ p_r_circ

double* Lorene::Etoile_rot::p_r_circ
mutableprotected

Circumferential equatorial radius.

Definition at line 1638 of file etoile.h.

◆ p_r_circ_merid

double* Lorene::Etoile_rot::p_r_circ_merid
mutableprotected

Circumferential meridional radius.

Definition at line 1639 of file etoile.h.

◆ p_r_isco

double* Lorene::Etoile_rot::p_r_isco
mutableprotected

Circumferential radius of the ISCO.

Definition at line 1648 of file etoile.h.

◆ p_ray_eq

double* Lorene::Etoile::p_ray_eq
mutableprotectedinherited

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
mutableprotectedinherited

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
mutableprotectedinherited

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
mutableprotectedinherited

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
mutableprotectedinherited

Coordinate radius at $\theta=0$.

Definition at line 536 of file etoile.h.

◆ p_tsw

double* Lorene::Etoile_rot::p_tsw
mutableprotected

Ratio T/W.

Definition at line 1635 of file etoile.h.

◆ p_xi_surf

Tbl* Lorene::Etoile::p_xi_surf
mutableprotectedinherited

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.

◆ p_z_eqb

double* Lorene::Etoile_rot::p_z_eqb
mutableprotected

Backward redshift factor at equator.

Definition at line 1643 of file etoile.h.

◆ p_z_eqf

double* Lorene::Etoile_rot::p_z_eqf
mutableprotected

Forward redshift factor at equator.

Definition at line 1642 of file etoile.h.

◆ p_z_pole

double* Lorene::Etoile_rot::p_z_pole
mutableprotected

Redshift factor at North pole.

Definition at line 1644 of file etoile.h.

◆ press

Tenseur Lorene::Etoile::press
protectedinherited

Fluid pressure.

Definition at line 464 of file etoile.h.

◆ relativistic

bool Lorene::Etoile::relativistic
protectedinherited

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
protectedinherited

Trace of the stress tensor in the Eulerian frame.

Definition at line 471 of file etoile.h.

◆ shift

Tenseur Lorene::Etoile::shift
protectedinherited

Total shift vector.

Definition at line 515 of file etoile.h.

◆ ssjm1_dzeta

Cmp Lorene::Etoile_rot::ssjm1_dzeta
protected

Effective source at the previous step for the resolution of the Poisson equation for dzeta .

Definition at line 1606 of file etoile.h.

◆ ssjm1_khi

Cmp Lorene::Etoile_rot::ssjm1_khi
protected

Effective source at the previous step for the resolution of the Poisson equation for the scalar $\chi$ by means of Map_et::poisson .

$\chi$ is an intermediate quantity for the resolution of the elliptic equation for the shift vector $N^i$

Definition at line 1619 of file etoile.h.

◆ ssjm1_nuf

Cmp Lorene::Etoile_rot::ssjm1_nuf
protected

Effective source at the previous step for the resolution of the Poisson equation for nuf by means of Map_et::poisson .

Definition at line 1595 of file etoile.h.

◆ ssjm1_nuq

Cmp Lorene::Etoile_rot::ssjm1_nuq
protected

Effective source at the previous step for the resolution of the Poisson equation for nuq by means of Map_et::poisson .

Definition at line 1601 of file etoile.h.

◆ ssjm1_tggg

Cmp Lorene::Etoile_rot::ssjm1_tggg
protected

Effective source at the previous step for the resolution of the Poisson equation for tggg .

Definition at line 1611 of file etoile.h.

◆ ssjm1_wshift

Tenseur Lorene::Etoile_rot::ssjm1_wshift
protected

Effective source at the previous step for the resolution of the vector Poisson equation for $W^i$.

$W^i$ is an intermediate quantity for the resolution of the elliptic equation for the shift vector $N^i$ (Components with respect to the Cartesian triad associated with the mapping mp )

Definition at line 1628 of file etoile.h.

◆ tggg

Tenseur Lorene::Etoile_rot::tggg
protected

Metric potential $\tilde G = (NB-1) r\sin\theta$.

Definition at line 1540 of file etoile.h.

◆ tkij

Tenseur_sym Lorene::Etoile_rot::tkij
protected

Tensor ${\tilde K_{ij}}$ related to the extrinsic curvature tensor by ${\tilde K_{ij}} = B^{-2} K_{ij}$.

tkij contains the Cartesian components of ${\tilde K_{ij}}$.

Definition at line 1570 of file etoile.h.

◆ tnphi

Tenseur Lorene::Etoile_rot::tnphi
protected

Component $\tilde N^\varphi = N^\varphi r\sin\theta$ of the shift vector.

Definition at line 1518 of file etoile.h.

◆ u_euler

Tenseur Lorene::Etoile::u_euler
protectedinherited

Fluid 3-velocity with respect to the Eulerian observer.

Definition at line 477 of file etoile.h.

◆ unsurc2

double Lorene::Etoile::unsurc2
protectedinherited

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

Definition at line 445 of file etoile.h.

◆ uuu

Tenseur Lorene::Etoile_rot::uuu
protected

Norm of u_euler.

Definition at line 1521 of file etoile.h.

◆ w_shift

Tenseur Lorene::Etoile_rot::w_shift
protected

Vector $W^i$ used in the decomposition of shift , following Shibata's prescription [Prog.

Theor. Phys. 101 , 1199 (1999)] :

\[ N^i = {7\over 8} W^i - {1\over 8} \left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

NB: w_shift contains the components of $W^i$ with respect to the Cartesian triad associated with the mapping mp .

Definition at line 1553 of file etoile.h.


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