Base class for axisymmetric stationary compact stars in Quasi-Isotropic coordinates (under development ).
More...
Star_QI (Map &mp_i)
Standard constructor.
Star_QI (const Star_QI &)
Copy constructor.
Star_QI (Map &mp_i, FILE *fich)
Constructor from a file (see sauve(FILE*) ).
virtual ~Star_QI ()
Destructor.
void operator= (const Star_QI &)
Assignment to another Star_QI .
const Scalar & get_logn () const
Returns the logarithm of the lapse N .
const Scalar & get_tnphi () const
Returns the component of the shift vector.
const Scalar & get_nuf () const
Returns the part of the Metric potential = logn generated by the matter terms.
const Scalar & get_nuq () const
Returns the Part of the Metric potential = logn generated by the quadratic terms.
const Scalar & get_dzeta () const
Returns the Metric potential .
const Scalar & get_tggg () const
Returns the Metric potential .
const Vector & get_w_shift () const
Returns the vector used in the decomposition of shift , following Shibata's prescription [Prog .
const Scalar & get_khi_shift () const
Returns the scalar used in the decomposition of shift
following Shibata's prescription [Prog .
virtual void sauve (FILE *) const
Save in a file.
virtual double mass_g () const
Gravitational mass.
virtual double angu_mom () const
Angular momentum.
virtual double grv2 () const
Error on the virial identity GRV2.
virtual double grv3 (ostream *ost=0x0) const
Error on the virial identity GRV3.
virtual double mom_quad () const
Quadrupole moment.
void update_metric ()
Computes metric coefficients from known potentials.
void fait_shift ()
Computes shift from w_shift and khi_shift according to Shibata's prescription [Prog .
void fait_nphi ()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
const Scalar & get_bbb () const
Returns the metric factor B .
const Scalar & get_a_car () const
Returns the square of the metric factor A .
const Scalar & get_b_car () const
Returns the square of the metric factor B .
const Scalar & get_nphi () const
Returns the metric coefficient .
const Scalar & get_ak_car () const
Returns the scalar .
void gyoto_data (const char *file_name) const
Save in a file for GYOTO.
virtual double r_isco (int lmin, ostream *ost=0x0) const
Coordinate r of the innermost stable circular orbit (ISCO).
virtual double f_isco (int lmin) const
Orbital frequency at the innermost stable circular orbit (ISCO).
virtual double espec_isco (int lmin) const
Energy of a particle at the ISCO.
virtual double lspec_isco (int lmin) const
Angular momentum of a particle at the ISCO.
virtual double r_mb (int lmin, ostream *ost=0x0) const
Coordinate r of the marginally bound circular orbit (R_mb).
virtual void extrinsic_curvature ()
Computes the extrinsic curvature and ak_car from nphi , nn and b_car .
Map & set_mp ()
Read/write of the mapping.
const Map & get_mp () const
Returns the mapping.
const Scalar & get_nn () const
Returns the lapse function N .
const Vector & get_beta () const
Returns the shift vector .
const Metric & get_gamma () const
Returns the 3-metric .
const Scalar & get_ener_euler () const
Returns the total energy density E in the Eulerian frame.
const Vector & get_mom_euler () const
Returns the total 3-momentum density in the Eulerian frame.
const Sym_tensor & get_stress_euler () const
Returns the stress tensor with respect to the Eulerian observer.
const Sym_tensor & get_kk () const
Returns the extrinsic curvature tensor .
virtual double adm_mass () const
ADM mass (computed as a surface integral at spatial infinity).
Base class for axisymmetric stationary compact stars in Quasi-Isotropic coordinates (under development ).
()
The time slice has the topology of and the metric is expressed in Quasi-Isotropic (QI) coordinates :
Definition at line 487 of file compobj.h .
Lorene::Star_QI::Star_QI
(
Map & mp_i ,
FILE * fich )
Constructor from a file (see sauve(FILE*) ).
Parameters
mp_i Mapping on which the star is constructed
fich input file (must have been created by the function Star_QI::sauve )
Definition at line 137 of file star_QI.C .
References Lorene::Compobj_QI::Compobj_QI() , dzeta , fait_nphi() , fait_shift() , Lorene::get_bvect_cart() , khi_shift , logn , Lorene::Map() , Lorene::Compobj::mp , nuf , nuq , set_der_0x0() , ssjm1_dzeta , ssjm1_khi , ssjm1_nuf , ssjm1_nuq , ssjm1_tggg , ssjm1_wshift , tggg , tnphi , update_metric() , and w_shift .
const Scalar & Lorene::Star_QI::get_khi_shift
(
)
const
inline
Returns the scalar used in the decomposition of shift
following Shibata's prescription [Prog .
Theor . Phys . 101 , 1199 (1999)] :
NB: w_shift contains the components of with respect to the Cartesian triad associated with the mapping mp .
Definition at line 686 of file compobj.h .
References khi_shift .
const Vector & Lorene::Star_QI::get_w_shift
(
)
const
inline
Returns the vector used in the decomposition of shift , following Shibata's prescription [Prog .
Theor . Phys . 101 , 1199 (1999)] :
NB: w_shift contains the components of with respect to the Cartesian triad associated with the mapping mp .
Definition at line 672 of file compobj.h .
References w_shift .
double Lorene::Star_QI::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
ost output stream to give details of the computation; if set to 0x0 [default value], no details will be given.
Definition at line 148 of file star_QI_global.C .
References Lorene::Compobj_QI::a_car , Lorene::Compobj_QI::ak_car , Lorene::Compobj_QI::bbb , Lorene::Scalar::derive_cov() , Lorene::Scalar::dsdr() , dzeta , Lorene::Compobj::gamma , Lorene::Scalar::get_dzpuis() , Lorene::Scalar::get_etat() , Lorene::Scalar::integrale() , Lorene::log() , logn , Lorene::Compobj::mp , Lorene::Valeur::mult_ct() , p_grv3 , Lorene::Scalar::set_dzpuis() , Lorene::Scalar::set_spectral_va() , Lorene::Scalar::srdsdt() , Lorene::Valeur::ssint() , Lorene::Scalar::std_spectral_base() , Lorene::Compobj::stress_euler , Lorene::Valeur::sx() , and Lorene::Map_radial::xsr .
double Lorene::Star_QI::lambda_grv2
(
const Scalar & sou_m ,
const Scalar & sou_q )
static
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
is the coefficient by which one must multiply the quadratic source term of the 2-D Poisson equation
in order that the total source does not contain any monopolar term, i.e. in order that
where . is computed according to the formula
Then, by construction, the new source has a vanishing monopolar term.
Parameters
sou_m [input] matter source term
sou_q [input] quadratic source term
Returns value of
Definition at line 321 of file star_QI_global.C .
References Lorene::Scalar::check_dzpuis() , Lorene::Valeur::coef_i() , Lorene::Map_radial::dxdr , Lorene::Map_af::get_alpha() , Lorene::Map_af::get_beta() , Lorene::Scalar::get_etat() , Lorene::Tbl::get_etat() , Lorene::Valeur::get_etat() , Lorene::Mg3d::get_grille3d() , Lorene::Tensor::get_mp() , Lorene::Mg3d::get_np() , Lorene::Mg3d::get_nr() , Lorene::Mg3d::get_nt() , Lorene::Mg3d::get_nzone() , Lorene::Scalar::get_spectral_va() , Lorene::Mg3d::get_type_r() , Lorene::Map_af::set_alpha() , Lorene::Map_af::set_beta() , Lorene::Mtbl::t , Lorene::Tbl::t , Lorene::Grille3d::x , and Lorene::Map_radial::xsr .
double Lorene::Star_QI::mom_quad
(
)
const
virtual
Quadrupole moment.
The quadrupole moment Q is defined according to Eq. (7) of [Salgado, Bonazzola, Gourgoulhon and Haensel, Astron . Astrophys . 291 , 155 (1994)]. At the Newtonian limit it is related to the component of the MTW (1973) reduced quadrupole moment by: . Note that Q is the negative of the quadrupole moment defined by Laarakkers and Poisson, Astrophys . J . 512 , 282 (1999).
Definition at line 250 of file star_QI_global.C .
References Lorene::Compobj_QI::a_car , Lorene::Compobj_QI::ak_car , Lorene::Compobj_QI::bbb , Lorene::Scalar::check_dzpuis() , Lorene::Scalar::derive_cov() , Lorene::Compobj::ener_euler , Lorene::Compobj::gamma , Lorene::Scalar::get_etat() , Lorene::Scalar::inc_dzpuis() , Lorene::Scalar::integrale() , Lorene::log() , logn , Lorene::Compobj::mp , Lorene::Valeur::mult_ct() , Lorene::Scalar::mult_r() , p_mom_quad , Lorene::Scalar::set_spectral_va() , Lorene::Scalar::std_spectral_base() , and Lorene::Compobj::stress_euler .
ostream & Lorene::Star_QI::operator>>
(
ostream & ost )
const
protected virtual
Operator >> (virtual function called by the operator <<).
Reimplemented from Lorene::Compobj_QI .
Reimplemented in Lorene::Boson_star .
Definition at line 300 of file star_QI.C .
References angu_mom() , dzeta , grv2() , grv3() , logn , mass_g() , mom_quad() , nuf , nuq , Lorene::Compobj_QI::operator>>() , and Lorene::pow() .
double Lorene::Compobj_QI::r_isco
(
int lmin ,
ostream * ost = 0x0 ) const
virtual inherited
Coordinate r of the innermost stable circular orbit (ISCO).
Parameters
lmin index of the innermost domain in which the ISCO is searched: the ISCO is searched inwards from the last but one domain to the domain of index lmin.
ost output stream to give details of the computation; if set to 0x0 [default value], no details will be given.
Definition at line 108 of file compobj_QI_global.C .
References Lorene::Param::add_int() , Lorene::Param::add_scalar() , Lorene::Tensor::annule_domain() , bbb , Lorene::Scalar::dsdr() , Lorene::Scalar::get_spectral_va() , Lorene::Compobj::mp , Lorene::Compobj::nn , nphi , p_espec_isco , p_f_isco , p_lspec_isco , p_r_isco , Lorene::r , Lorene::sqrt() , Lorene::Scalar::std_spectral_base() , Lorene::Valeur::val_point() , and Lorene::zerosec() .
void Lorene::Star_QI::update_metric
(
)
virtual
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, as well as the 3-metric gamma.
Reimplemented from Lorene::Compobj_QI .
Definition at line 410 of file star_QI.C .
References Lorene::Compobj_QI::a_car , Lorene::Compobj_QI::b_car , Lorene::Compobj_QI::bbb , del_deriv() , Lorene::Scalar::div_rsint() , dzeta , Lorene::exp() , logn , Lorene::Compobj::nn , tggg , and Lorene::Compobj_QI::update_metric() .
Vector Lorene::Star_QI::w_shift
protected
Vector used in the decomposition of shift , following Shibata's prescription [Prog .
Theor . Phys . 101 , 1199 (1999)] :
NB: w_shift contains the components of with respect to the Cartesian triad associated with the mapping mp .
Definition at line 529 of file compobj.h .