Class for isolated rotating stars. More...
#include <star_rot.h>
Public Member Functions | |
Star_rot (Map &mp_i, int nzet_i, bool relat, const Eos &eos_i) | |
Standard constructor. | |
Star_rot (const Star_rot &) | |
Copy constructor. | |
Star_rot (Map &mp_i, const Eos &eos_i, FILE *fich) | |
Constructor from a file (see sauve(FILE*) ). | |
virtual | ~Star_rot () |
Destructor. | |
void | operator= (const Star_rot &) |
Assignment to another Star_rot . | |
bool | is_relativistic () const |
Returns true for a relativistic star, false for a Newtonian one. | |
virtual double | get_omega_c () const |
Returns the central value of the rotation angular velocity ( [f_unit] ). | |
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_tnphi () const |
Returns the component ![]() | |
const Scalar & | get_uuu () const |
Returns the norm of u_euler . | |
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 ![]() shift , following Shibata's prescription [Prog. | |
const Scalar & | get_khi_shift () const |
Returns the scalar ![]() shift following Shibata's prescription [Prog. | |
const Sym_tensor & | get_tkij () const |
Returns the tensor ![]() ![]() | |
const Scalar & | get_ak_car () const |
Returns the scalar ![]() | |
virtual void | sauve (FILE *) const |
Save in a file. | |
virtual void | display_poly (ostream &) const |
Display in polytropic units. | |
virtual const Itbl & | l_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 ![]() | |
virtual double | mass_b () const |
Baryon mass. | |
virtual double | mass_g () const |
Gravitational mass. | |
virtual double | angu_mom () const |
Angular momentum. | |
virtual double | tsw () const |
Ratio T/W. | |
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 | r_circ () const |
Circumferential radius. | |
virtual double | aplat () const |
Flatening r_pole/r_eq. | |
virtual double | z_eqf () const |
Forward redshift factor at equator. | |
virtual double | z_eqb () const |
Backward redshift factor at equator. | |
virtual double | z_pole () const |
Redshift factor at North pole. | |
virtual double | mom_quad () const |
Quadrupole moment. | |
virtual double | r_isco (ostream *ost=0x0) const |
Circumferential radius of the innermost stable circular orbit (ISCO). | |
virtual double | f_isco () const |
Orbital frequency at the innermost stable circular orbit (ISCO). | |
virtual double | espec_isco () const |
Energy of a particle on the ISCO. | |
virtual double | lspec_isco () const |
Angular momentum of a particle on the ISCO. | |
virtual double | f_eq () const |
Orbital frequency at the equator. | |
virtual void | hydro_euler () |
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame. | |
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 . | |
void | extrinsic_curvature () |
Computes tkij and ak_car from shift , nnn and b_car . | |
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. | |
Map & | set_mp () |
Read/write of the mapping. | |
void | set_enthalpy (const Scalar &) |
Assignment of the enthalpy field. | |
void | equation_of_state () |
Computes the proper baryon and energy density, as well as pressure from the enthalpy. | |
virtual void | equilibrium_spher (double ent_c, double precis=1.e-14, const Tbl *pent_limit=0x0) |
Computes a spherical static configuration. | |
const Map & | get_mp () const |
Returns the mapping. | |
int | get_nzet () const |
Returns the number of domains occupied by the star. | |
const Eos & | get_eos () const |
Returns the equation of state. | |
const Scalar & | get_ent () const |
Returns the enthalpy field. | |
const Scalar & | get_nbar () const |
Returns the proper baryon density. | |
const Scalar & | get_ener () const |
Returns the proper total energy density. | |
const Scalar & | get_press () const |
Returns the fluid pressure. | |
const Scalar & | get_ener_euler () const |
Returns the total energy density with respect to the Eulerian observer. | |
const Scalar & | get_s_euler () const |
Returns the trace of the stress tensor in the Eulerian frame. | |
const Scalar & | get_gam_euler () const |
Returns the Lorentz factor between the fluid and Eulerian observers. | |
const Vector & | get_u_euler () const |
Returns the fluid 3-velocity with respect to the Eulerian observer. | |
const Tensor & | get_stress_euler () const |
Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer. | |
const Scalar & | get_logn () const |
Returns the logarithm of the lapse N. | |
const Scalar & | get_nn () const |
Returns the lapse function N. | |
const Vector & | get_beta () const |
Returns the shift vector ![]() | |
const Scalar & | get_lnq () const |
const Metric & | get_gamma () const |
Returns the 3-metric ![]() | |
double | ray_eq () const |
Coordinate radius at ![]() ![]() | |
double | ray_eq_pis2 () const |
Coordinate radius at ![]() ![]() | |
double | ray_eq_pi () const |
Coordinate radius at ![]() ![]() | |
double | ray_eq_3pis2 () const |
Coordinate radius at ![]() ![]() | |
double | ray_pole () const |
Coordinate radius at ![]() | |
const Tbl & | xi_surf () const |
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate ![]() ![]() | |
Static Public Member Functions | |
static double | lambda_grv2 (const Scalar &sou_m, const Scalar &sou_q) |
Computes the coefficient ![]() | |
Protected Member Functions | |
virtual void | del_deriv () const |
Deletes all the derived quantities. | |
virtual void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. | |
virtual void | del_hydro_euler () |
Sets to ETATNONDEF (undefined state) the hydrodynamical quantities relative to the Eulerian observer. | |
virtual ostream & | operator>> (ostream &) const |
Operator >> (virtual function called by the operator <<). | |
virtual void | partial_display (ostream &) const |
Printing of some informations, excluding all global quantities. | |
Protected Attributes | |
bool | relativistic |
Indicator of relativity: true for a relativistic star, false for a Newtonian one. | |
double | unsurc2 |
![]() unsurc2=1 for a relativistic star, 0 for a Newtonian one. | |
double | omega |
Rotation angular velocity ( [f_unit] ). | |
Scalar | a_car |
Square of the metric factor A. | |
Scalar | bbb |
Metric factor B. | |
Scalar | b_car |
Square of the metric factor B. | |
Scalar | nphi |
Metric coefficient ![]() | |
Scalar | tnphi |
Component ![]() | |
Scalar | uuu |
Norm of u_euler . | |
Scalar | nuf |
Part of the Metric potential ![]() logn generated by the matter terms. | |
Scalar | nuq |
Part of the Metric potential ![]() logn generated by the quadratic terms. | |
Scalar | dzeta |
Metric potential ![]() | |
Scalar | tggg |
Metric potential ![]() | |
Vector | w_shift |
Vector ![]() shift , following Shibata's prescription [Prog. | |
Scalar | khi_shift |
Scalar ![]() shift , following Shibata's prescription [Prog. | |
Sym_tensor | tkij |
Tensor ![]() ![]() | |
Scalar | ak_car |
Scalar ![]() | |
Scalar | ssjm1_nuf |
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of Map_et::poisson . | |
Scalar | ssjm1_nuq |
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of Map_et::poisson . | |
Scalar | ssjm1_dzeta |
Effective source at the previous step for the resolution of the Poisson equation for dzeta . | |
Scalar | ssjm1_tggg |
Effective source at the previous step for the resolution of the Poisson equation for tggg . | |
Scalar | ssjm1_khi |
Effective source at the previous step for the resolution of the Poisson equation for the scalar ![]() Map_et::poisson . | |
Vector | ssjm1_wshift |
Effective source at the previous step for the resolution of the vector Poisson equation for ![]() | |
double * | p_angu_mom |
Angular momentum. | |
double * | p_tsw |
Ratio T/W. | |
double * | p_grv2 |
Error on the virial identity GRV2. | |
double * | p_grv3 |
Error on the virial identity GRV3. | |
double * | p_r_circ |
Circumferential radius. | |
double * | p_aplat |
Flatening r_pole/r_eq. | |
double * | p_z_eqf |
Forward redshift factor at equator. | |
double * | p_z_eqb |
Backward redshift factor at equator. | |
double * | p_z_pole |
Redshift factor at North pole. | |
double * | p_mom_quad |
Quadrupole moment. | |
double * | p_r_isco |
Circumferential radius of the ISCO. | |
double * | p_f_isco |
Orbital frequency of the ISCO. | |
double * | p_espec_isco |
Specific energy of a particle on the ISCO. | |
double * | p_lspec_isco |
Specific angular momentum of a particle on the ISCO. | |
double * | p_f_eq |
Orbital frequency at the equator. | |
Map & | mp |
Mapping associated with the star. | |
int | nzet |
Number of domains of *mp occupied by the star. | |
const Eos & | eos |
Equation of state of the stellar matter. | |
Scalar | ent |
Log-enthalpy. | |
Scalar | nbar |
Baryon density in the fluid frame. | |
Scalar | ener |
Total energy density in the fluid frame. | |
Scalar | press |
Fluid pressure. | |
Scalar | ener_euler |
Total energy density in the Eulerian frame. | |
Scalar | s_euler |
Trace of the stress scalar in the Eulerian frame. | |
Scalar | gam_euler |
Lorentz factor between the fluid and Eulerian observers. | |
Vector | u_euler |
Fluid 3-velocity with respect to the Eulerian observer. | |
Sym_tensor | stress_euler |
Spatial part of the stress-energy tensor with respect to the Eulerian observer. | |
Scalar | logn |
Logarithm of the lapse N . | |
Scalar | nn |
Lapse function N . | |
Vector | beta |
Shift vector. | |
Scalar | lnq |
Metric | gamma |
3-metric | |
double * | p_ray_eq |
Coordinate radius at ![]() ![]() | |
double * | p_ray_eq_pis2 |
Coordinate radius at ![]() ![]() | |
double * | p_ray_eq_pi |
Coordinate radius at ![]() ![]() | |
double * | p_ray_eq_3pis2 |
Coordinate radius at ![]() ![]() | |
double * | p_ray_pole |
Coordinate radius at ![]() | |
Itbl * | p_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 ![]() | |
Tbl * | p_xi_surf |
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate ![]() ![]() | |
double * | p_mass_b |
Baryon mass. | |
double * | p_mass_g |
Gravitational mass. | |
Friends | |
ostream & | operator<< (ostream &, const Star &) |
Display. |
Class for isolated rotating stars.
()
The metric is
A star of class Star_rot
can be either relativistic or Newtonian, depending on the boolean indicator relativistic
. For a Newtonian star, the metric coefficients N, A, and B are set to 1, and is set to zero; the only relevant gravitational quantity in this case is
logn
which represents the Newtonian gravitational potential generated by the star.
Definition at line 79 of file star_rot.h.
Standard constructor.
mp_i | Mapping on which the star is contructed | |
nzet_i | Number of domains occupied by the star | |
relat | true for a relativistic star, false for a Newtonian one | |
eos_i | Equation of state of the stellar matter |
Definition at line 79 of file star_rot.C.
References a_car, ak_car, b_car, bbb, Star::beta, dzeta, Map::get_bvect_cart(), khi_shift, Star::mp, nphi, nuf, nuq, omega, relativistic, set_der_0x0(), Tensor::set_etat_zero(), Tensor::set_triad(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, Scalar::std_spectral_base(), tggg, tkij, tnphi, unsurc2, uuu, and w_shift.
Star_rot::Star_rot | ( | const Star_rot & | et | ) |
Constructor from a file (see sauve(FILE*)
).
mp_i | Mapping on which the star is constructed | |
eos_i | Equation of state of the stellar matter | |
fich | input file (must have been created by the function Star_rot::sauve ) |
Definition at line 182 of file star_rot.C.
References a_car, ak_car, b_car, bbb, dzeta, fait_nphi(), fait_shift(), fread_be(), Map::get_bvect_cart(), Map::get_mg(), khi_shift, Star::mp, nuf, nuq, omega, relativistic, set_der_0x0(), Tensor::set_etat_zero(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, tggg, tkij, unsurc2, uuu, and w_shift.
Star_rot::~Star_rot | ( | ) | [virtual] |
double Star_rot::angu_mom | ( | ) | const [virtual] |
Angular momentum.
Definition at line 150 of file star_rot_global.C.
References a_car, b_car, Star::ener_euler, Scalar::get_spectral_va(), Scalar::integrale(), Scalar::mult_r(), Star::nbar, p_angu_mom, Star::press, relativistic, Scalar::set_spectral_va(), and uuu.
double Star_rot::aplat | ( | ) | const [virtual] |
Flatening r_pole/r_eq.
Definition at line 373 of file star_rot_global.C.
References p_aplat, Star::ray_eq(), and Star::ray_pole().
void Star_rot::del_deriv | ( | ) | const [protected, virtual] |
Deletes all the derived quantities.
Reimplemented from Star.
Definition at line 290 of file star_rot.C.
References p_angu_mom, p_aplat, p_espec_isco, p_f_eq, p_f_isco, p_grv2, p_grv3, p_lspec_isco, p_mom_quad, p_r_circ, p_r_isco, p_tsw, p_z_eqb, p_z_eqf, p_z_pole, and set_der_0x0().
void Star_rot::del_hydro_euler | ( | ) | [protected, virtual] |
Sets to ETATNONDEF
(undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Reimplemented from Star.
Definition at line 336 of file star_rot.C.
References del_deriv().
void Star_rot::display_poly | ( | ostream & | ost | ) | const [virtual] |
Display in polytropic units.
Definition at line 636 of file star_rot.C.
References angu_mom(), Star::ener, Star::eos, Eos_poly::get_gam(), Eos_poly::get_kap(), get_omega_c(), mass_b(), mass_g(), Star::nbar, pow(), r_circ(), Star::ray_eq(), sqrt(), and Scalar::val_grid_point().
void Star::equation_of_state | ( | ) | [inherited] |
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Reimplemented in Gravastar.
Definition at line 458 of file star.C.
References Param::add_int(), Scalar::allocate_all(), Star::del_deriv(), Star::ener, Eos::ener_ent(), Star::ent, Star::eos, Mg3d::get_grille3d(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Star::mp, Star::nbar, Eos::nbar_ent(), Star::nzet, Star::press, Eos::press_ent(), Mtbl::set(), Scalar::set_domain(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Mtbl::set_etat_qcq(), Scalar::std_spectral_base(), Mtbl::t, and Grille3d::x.
void Star_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.
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:
| |
control | [input] Set of parameters (stored as a 1-D Tbl of size 7) to control the iteration:
| |
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 :
|
Definition at line 67 of file star_rot_equil.C.
References a_car, abs(), Map::adapt(), Param::add_cmp_mod(), Param::add_double(), Param::add_double_mod(), Param::add_int(), Param::add_int_mod(), Param::add_tbl(), Param::add_tenseur_mod(), ak_car, Tensor::annule_domain(), bbb, Star::beta, Valeur::c_cf, Vector::change_triad(), Map::cmp_zero(), Valeur::coef(), contract(), cos(), Scalar::derive_con(), Scalar::derive_cov(), diffrel(), Scalar::dsdr(), dzeta, Star::ener_euler, Star::ent, Star::equation_of_state(), fait_nphi(), Map::flat_met_cart(), Map::flat_met_spher(), Star::gam_euler, Map::get_bvect_cart(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_t(), grv2(), Map_et::homothetie(), hydro_euler(), khi_shift, log(), log10(), Star::logn, mass_b(), mass_g(), Star::mp, Scalar::mult_rsint(), Star::nbar, Star::nn, nphi, nuf, nuq, Star::nzet, omega, partial_display(), Map::phi, Scalar::poisson(), Map::poisson2d(), pow(), Star::press, Star::ray_eq(), Star::ray_pole(), Map::reevaluate(), relativistic, Star::s_euler, Vector::set(), Tenseur::set(), Tbl::set(), Scalar::set_dzpuis(), Tenseur::set_etat_qcq(), Tbl::set_etat_qcq(), Map::sint, sqrt(), ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, Scalar::std_spectral_base(), Scalar::test_poisson(), tggg, tkij, Star::u_euler, Tensor::up(), update_metric(), uuu, Scalar::val_grid_point(), and w_shift.
void Star::equilibrium_spher | ( | double | ent_c, | |
double | precis = 1.e-14 , |
|||
const Tbl * | pent_limit = 0x0 | |||
) | [virtual, inherited] |
Computes a spherical static configuration.
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 91 of file star_equil_spher.C.
References Map_et::adapt(), Param::add_double(), Param::add_int(), Param::add_int_mod(), Param::add_tbl(), Scalar::annule(), diffrel(), Scalar::dsdr(), Map_af::dsdr(), Star::ener, Star::ener_euler, Star::ent, Star::equation_of_state(), exp(), Star::gam_euler, Star::gamma, Map_et::get_alpha(), Map_af::get_alpha(), Map_et::get_beta(), Map_af::get_beta(), Map::get_bvect_cart(), Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Map_af::homothetie(), Scalar::integrale(), Star::logn, Star::mass_b(), Star::mass_g(), Star::mp, Star::nn, norme(), Star::nzet, Map_af::poisson(), Star::press, Star::s_euler, Vector::set(), Map_af::set_alpha(), Map_af::set_beta(), Scalar::set_dzpuis(), Cmp::set_etat_qcq(), Scalar::set_etat_zero(), sqrt(), Scalar::std_spectral_base(), Star::u_euler, Scalar::val_grid_point(), and Map::val_r().
double Star_rot::espec_isco | ( | ) | const [virtual] |
Energy of a particle on the ISCO.
Definition at line 285 of file star_rot_isco.C.
References p_espec_isco, and r_isco().
void Star_rot::extrinsic_curvature | ( | ) |
Computes tkij
and ak_car
from shift
, nnn
and b_car
.
Definition at line 106 of file star_rot_upmetr.C.
References ak_car, b_car, Star::beta, contract(), Tensor::derive_cov(), Scalar::dsdr(), Map::flat_met_cart(), Scalar::get_etat(), Map::get_mg(), Scalar::get_spectral_va(), Star::mp, Scalar::mult_rsint(), Star::nn, nphi, Tensor::set(), Tensor::set_etat_qcq(), Tensor::set_etat_zero(), Scalar::set_spectral_va(), Scalar::srdsdt(), and tkij.
double Star_rot::f_eq | ( | ) | const [virtual] |
Orbital frequency at the equator.
Definition at line 303 of file star_rot_isco.C.
double Star_rot::f_isco | ( | ) | const [virtual] |
Orbital frequency at the innermost stable circular orbit (ISCO).
Definition at line 251 of file star_rot_isco.C.
void Star_rot::fait_nphi | ( | ) |
Computes tnphi
and nphi
from the Cartesian components of the shift, stored in shift
.
Definition at line 735 of file star_rot.C.
References Star::beta, Map::comp_p_from_cartesian(), Scalar::div_rsint(), Star::mp, nphi, and tnphi.
void Star_rot::fait_shift | ( | ) |
Computes shift
from w_shift
and khi_shift
according to Shibata's prescription [Prog.
Theor. Phys. 101 , 1199 (1999)] :
Definition at line 690 of file star_rot.C.
References Star::beta, contract(), Map::cosp, Map::cost, Tensor::dec_dzpuis(), Tensor::derive_con(), Scalar::derive_con(), Scalar::domain(), Map::flat_met_cart(), Map::get_bvect_cart(), Map::get_mg(), Mg3d::get_nzone(), khi_shift, Star::mp, Tbl::set(), Tensor::set(), Scalar::set_domain(), Map::sinp, Map::sint, w_shift, Map::x, Map::y, and Map::z.
const Scalar& Star_rot::get_a_car | ( | ) | const [inline] |
Returns the square of the metric factor A.
Definition at line 311 of file star_rot.h.
References a_car.
const Scalar& Star_rot::get_ak_car | ( | ) | const [inline] |
Returns the scalar .
For axisymmetric stars, this quantity is related to the derivatives of by
In particular it is related to the quantities and
introduced by Eqs. (3.7) and (3.8) of Bonazzola et al. Astron. Astrophys. 278 , 421 (1993) by
Definition at line 395 of file star_rot.h.
References ak_car.
const Scalar& Star_rot::get_b_car | ( | ) | const [inline] |
Returns the square of the metric factor B.
Definition at line 314 of file star_rot.h.
References b_car.
const Scalar& Star_rot::get_bbb | ( | ) | const [inline] |
const Vector& Star::get_beta | ( | ) | const [inline, inherited] |
const Scalar& Star_rot::get_dzeta | ( | ) | const [inline] |
const Scalar& Star::get_ener | ( | ) | const [inline, inherited] |
Returns the proper total energy density.
Definition at line 366 of file star.h.
References Star::ener.
const Scalar& Star::get_ener_euler | ( | ) | const [inline, inherited] |
Returns the total energy density with respect to the Eulerian observer.
Definition at line 372 of file star.h.
References Star::ener_euler.
const Scalar& Star::get_ent | ( | ) | const [inline, inherited] |
const Eos& Star::get_eos | ( | ) | const [inline, inherited] |
const Scalar& Star::get_gam_euler | ( | ) | const [inline, inherited] |
Returns the Lorentz factor between the fluid and Eulerian observers.
Definition at line 378 of file star.h.
References Star::gam_euler.
const Metric& Star::get_gamma | ( | ) | const [inline, inherited] |
const Scalar& Star_rot::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 369 of file star_rot.h.
References khi_shift.
const Scalar& Star::get_logn | ( | ) | const [inline, inherited] |
Returns the logarithm of the lapse N.
In the Newtonian case, this is the Newtonian gravitational potential (in units of ).
Definition at line 392 of file star.h.
References Star::logn.
const Map& Star::get_mp | ( | ) | const [inline, inherited] |
const Scalar& Star::get_nbar | ( | ) | const [inline, inherited] |
const Scalar& Star::get_nn | ( | ) | const [inline, inherited] |
const Scalar& Star_rot::get_nphi | ( | ) | const [inline] |
const Scalar& Star_rot::get_nuf | ( | ) | const [inline] |
Returns the part of the Metric potential =
logn
generated by the matter terms.
Definition at line 330 of file star_rot.h.
References nuf.
const Scalar& Star_rot::get_nuq | ( | ) | const [inline] |
Returns the Part of the Metric potential =
logn
generated by the quadratic terms.
Definition at line 335 of file star_rot.h.
References nuq.
int Star::get_nzet | ( | ) | const [inline, inherited] |
Returns the number of domains occupied by the star.
Definition at line 354 of file star.h.
References Star::nzet.
double Star_rot::get_omega_c | ( | ) | const [virtual] |
Returns the central value of the rotation angular velocity ([f_unit] ).
Definition at line 627 of file star_rot.C.
References omega.
const Scalar& Star::get_press | ( | ) | const [inline, inherited] |
const Scalar& Star::get_s_euler | ( | ) | const [inline, inherited] |
Returns the trace of the stress tensor in the Eulerian frame.
Definition at line 375 of file star.h.
References Star::s_euler.
const Tensor& Star::get_stress_euler | ( | ) | const [inline, inherited] |
Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer.
Definition at line 386 of file star.h.
References Star::stress_euler.
const Scalar& Star_rot::get_tggg | ( | ) | const [inline] |
const Sym_tensor& Star_rot::get_tkij | ( | ) | const [inline] |
Returns the tensor related to the extrinsic curvature tensor by
.
tkij
contains the Cartesian components of .
Definition at line 376 of file star_rot.h.
References tkij.
const Scalar& Star_rot::get_tnphi | ( | ) | const [inline] |
Returns the component of the shift vector.
Definition at line 322 of file star_rot.h.
References tnphi.
const Vector& Star::get_u_euler | ( | ) | const [inline, inherited] |
Returns the fluid 3-velocity with respect to the Eulerian observer.
Definition at line 381 of file star.h.
References Star::u_euler.
const Scalar& Star_rot::get_uuu | ( | ) | const [inline] |
const Vector& Star_rot::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 355 of file star_rot.h.
References w_shift.
double Star_rot::grv2 | ( | ) | const [virtual] |
Error on the virial identity GRV2.
This indicator is only valid for relativistic computations.
Definition at line 210 of file star_rot_global.C.
References a_car, ak_car, Scalar::derive_cov(), Star::ener_euler, Map::flat_met_spher(), lambda_grv2(), Star::logn, Star::mp, Star::nbar, p_grv2, Star::press, relativistic, and uuu.
double Star_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.
ost | output stream to give details of the computation; if set to 0x0 [default value], no details will be given. |
Definition at line 241 of file star_rot_global.C.
References a_car, ak_car, bbb, Scalar::derive_cov(), Scalar::dsdr(), dzeta, Map::flat_met_spher(), Scalar::get_dzpuis(), Scalar::get_etat(), Scalar::integrale(), log(), Star::logn, Star::mp, Valeur::mult_ct(), Star::nbar, p_grv3, Star::press, relativistic, Star::s_euler, Scalar::set_dzpuis(), Scalar::set_spectral_va(), Scalar::srdsdt(), Valeur::ssint(), Scalar::std_spectral_base(), Valeur::sx(), uuu, and Map_radial::xsr.
void Star_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 Star.
Definition at line 53 of file star_rot_hydro.C.
References a_car, Tensor::annule_domain(), b_car, Star::beta, Vector::change_triad(), del_deriv(), Star::ener, Star::ener_euler, Star::gam_euler, Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_nzone(), Star::mp, Star::nn, omega, Star::press, Star::s_euler, Vector::set(), Tensor::set_etat_qcq(), Tensor::set_etat_zero(), Scalar::set_spectral_va(), Tensor::set_triad(), sqrt(), Scalar::std_spectral_base(), Vector::std_spectral_base(), Star::u_euler, unsurc2, uuu, Map::x, and Map::y.
bool Star_rot::is_relativistic | ( | ) | const [inline] |
Returns true
for a relativistic star, false
for a Newtonian one.
Definition at line 300 of file star_rot.h.
References relativistic.
const Itbl & Star_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 .
The stellar surface is defined as the location where the enthalpy (member ent
) vanishes.
Reimplemented from Star.
Definition at line 61 of file star_rot_global.C.
References Star::ent, Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Scalar::get_spectral_va(), Star::mp, Star::nzet, Star::p_l_surf, and Star::p_xi_surf.
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.
sou_m | [input] matter source term ![]() | |
sou_q | [input] quadratic source term ![]() |
Definition at line 58 of file star_rot_lambda_grv2.C.
References Valeur::c, Scalar::check_dzpuis(), Valeur::coef_i(), Map_radial::dxdr, Map_af::get_alpha(), Map_af::get_beta(), Valeur::get_etat(), Scalar::get_etat(), Tbl::get_etat(), Mg3d::get_grille3d(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Map_af::set_alpha(), Map_af::set_beta(), Tbl::t, Mtbl::t, Map::val_r(), Grille3d::x, and Map_radial::xsr.
double Star_rot::lspec_isco | ( | ) | const [virtual] |
Angular momentum of a particle on the ISCO.
Definition at line 268 of file star_rot_isco.C.
References p_lspec_isco, and r_isco().
double Star_rot::mass_b | ( | ) | const [virtual] |
Baryon mass.
Implements Star.
Definition at line 91 of file star_rot_global.C.
References a_car, bbb, Star::gam_euler, Scalar::get_etat(), Scalar::integrale(), Star::nbar, Star::p_mass_b, and relativistic.
double Star_rot::mass_g | ( | ) | const [virtual] |
Gravitational mass.
Implements Star.
Definition at line 120 of file star_rot_global.C.
References a_car, bbb, Star::ener_euler, mass_b(), Star::nn, Star::p_mass_g, Star::press, relativistic, Star::s_euler, Scalar::std_spectral_base(), tnphi, and uuu.
double Star_rot::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 468 of file star_rot_global.C.
References a_car, ak_car, bbb, Scalar::check_dzpuis(), Scalar::derive_cov(), Star::ener_euler, Map::flat_met_spher(), Scalar::get_etat(), Scalar::inc_dzpuis(), Scalar::integrale(), log(), Star::logn, Star::mp, Valeur::mult_ct(), Scalar::mult_r(), Star::nbar, p_mom_quad, relativistic, Star::s_euler, Scalar::set_spectral_va(), and Scalar::std_spectral_base().
void Star_rot::operator= | ( | const Star_rot & | et | ) |
Assignment to another Star_rot
.
Reimplemented from Star.
Reimplemented in Gravastar.
Definition at line 350 of file star_rot.C.
References a_car, ak_car, b_car, bbb, del_deriv(), dzeta, khi_shift, nphi, nuf, nuq, omega, relativistic, ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, tggg, tkij, tnphi, unsurc2, uuu, and w_shift.
ostream & Star_rot::operator>> | ( | ostream & | ost | ) | const [protected, virtual] |
Operator >> (virtual function called by the operator <<).
Reimplemented from Star.
Reimplemented in Gravastar.
Definition at line 419 of file star_rot.C.
References a_car, angu_mom(), aplat(), b_car, diffrel(), dzeta, f_isco(), Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), get_omega_c(), Tbl::get_taille(), grv2(), grv3(), khi_shift, Star::logn, mass_g(), mom_quad(), Star::mp, nphi, Star::nzet, omega, pow(), r_circ(), r_isco(), Star::ray_eq(), relativistic, sqrt(), tsw(), uuu, Scalar::val_grid_point(), w_shift, z_eqb(), z_eqf(), and z_pole().
void Star_rot::partial_display | ( | ostream & | ost | ) | const [protected, virtual] |
Printing of some informations, excluding all global quantities.
Definition at line 577 of file star_rot.C.
References a_car, aplat(), b_car, dzeta, Star::ener, Star::ent, Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), get_omega_c(), Star::logn, Star::mp, Star::nbar, Star::nn, nphi, Star::nzet, Star::press, Star::ray_eq(), uuu, and Scalar::val_grid_point().
double Star_rot::r_circ | ( | ) | const [virtual] |
Circumferential radius.
Definition at line 348 of file star_rot_global.C.
References bbb, Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_type_t(), Star::mp, Star::nzet, p_r_circ, Star::ray_eq(), and Scalar::val_grid_point().
double Star_rot::r_isco | ( | ostream * | ost = 0x0 |
) | const [virtual] |
Circumferential radius of the innermost stable circular orbit (ISCO).
ost | output stream to give details of the computation; if set to 0x0 [default value], no details will be given. |
Definition at line 67 of file star_rot_isco.C.
References Param::add_int(), Param::add_scalar(), Tensor::annule_domain(), bbb, Scalar::dsdr(), Map::get_mg(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Star::mp, Star::nn, nphi, Star::nzet, p_espec_isco, p_f_eq, p_f_isco, p_lspec_isco, p_r_isco, Map::r, Star::ray_eq(), sqrt(), Scalar::std_spectral_base(), Valeur::val_point(), Map::val_r(), and zerosec().
double Star::ray_eq | ( | ) | const [inherited] |
Coordinate radius at ,
[r_unit].
Definition at line 104 of file star_global.C.
References Map::get_mg(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), Star::l_surf(), Star::mp, Star::p_ray_eq, Map::val_r(), and Star::xi_surf().
double Star::ray_eq_3pis2 | ( | ) | const [inherited] |
Coordinate radius at ,
[r_unit].
Definition at line 229 of file star_global.C.
References Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), Star::l_surf(), Star::mp, Star::p_ray_eq_3pis2, Star::ray_eq_pis2(), Map::val_r(), and Star::xi_surf().
double Star::ray_eq_pi | ( | ) | const [inherited] |
Coordinate radius at ,
[r_unit].
Definition at line 182 of file star_global.C.
References Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), Star::l_surf(), Star::mp, Star::p_ray_eq_pi, Star::ray_eq(), Map::val_r(), and Star::xi_surf().
double Star::ray_eq_pis2 | ( | ) | const [inherited] |
Coordinate radius at ,
[r_unit].
Definition at line 134 of file star_global.C.
References Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_type_p(), Mg3d::get_type_t(), Star::l_surf(), Star::mp, Star::p_ray_eq_pis2, Map::val_r(), and Star::xi_surf().
double Star::ray_pole | ( | ) | const [inherited] |
Coordinate radius at [r_unit].
Definition at line 274 of file star_global.C.
References Map::get_mg(), Mg3d::get_type_t(), Star::l_surf(), Star::mp, Star::p_ray_pole, Map::val_r(), and Star::xi_surf().
void Star_rot::sauve | ( | FILE * | fich | ) | const [virtual] |
Save in a file.
Reimplemented from Star.
Definition at line 392 of file star_rot.C.
References dzeta, fwrite_be(), khi_shift, nuf, nuq, omega, relativistic, Tensor::sauve(), Scalar::sauve(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, tggg, and w_shift.
void Star_rot::set_der_0x0 | ( | ) | const [protected, virtual] |
Sets to 0x0
all the pointers on derived quantities.
Reimplemented from Star.
Definition at line 314 of file star_rot.C.
References p_angu_mom, p_aplat, p_espec_isco, p_f_eq, p_f_isco, p_grv2, p_grv3, p_lspec_isco, p_mom_quad, p_r_circ, p_r_isco, p_tsw, p_z_eqb, p_z_eqf, and p_z_pole.
void Star::set_enthalpy | ( | const Scalar & | ent_i | ) | [inherited] |
Assignment of the enthalpy field.
Definition at line 375 of file star.C.
References Star::del_deriv(), Star::ent, and Star::equation_of_state().
Map& Star::set_mp | ( | ) | [inline, inherited] |
double Star_rot::tsw | ( | ) | const [virtual] |
Ratio T/W.
Definition at line 179 of file star_rot_global.C.
References a_car, angu_mom(), bbb, Star::ener, Star::gam_euler, Scalar::integrale(), Star::logn, mass_g(), Star::nbar, omega, p_tsw, and relativistic.
void Star_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
, as well as the 3-metric gamma
.
Definition at line 47 of file star_rot_upmetr.C.
References a_car, b_car, bbb, del_deriv(), Scalar::div_rsint(), dzeta, exp(), extrinsic_curvature(), Star::gamma, Map::get_bvect_spher(), Star::logn, Star::mp, Star::nn, Tensor::set(), Scalar::std_spectral_base(), tggg, and unsurc2.
const Tbl & Star::xi_surf | ( | ) | const [inherited] |
Description of the stellar surface: returns a 2-D Tbl
containing the values of the radial coordinate on the surface at the collocation points in
.
The stellar surface is defined as the location where the enthalpy (member ent
) vanishes.
Definition at line 85 of file star_global.C.
References Star::l_surf(), Star::p_l_surf, and Star::p_xi_surf.
double Star_rot::z_eqb | ( | ) | const [virtual] |
Backward redshift factor at equator.
Definition at line 418 of file star_rot_global.C.
References Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_type_t(), Star::mp, Star::nn, nphi, Star::nzet, p_z_eqb, r_circ(), sqrt(), uuu, and Scalar::val_grid_point().
double Star_rot::z_eqf | ( | ) | const [virtual] |
Forward redshift factor at equator.
Definition at line 390 of file star_rot_global.C.
References Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_type_t(), Star::mp, Star::nn, nphi, Star::nzet, p_z_eqf, r_circ(), sqrt(), uuu, and Scalar::val_grid_point().
double Star_rot::z_pole | ( | ) | const [virtual] |
Redshift factor at North pole.
Definition at line 449 of file star_rot_global.C.
References Star::nn, p_z_pole, Star::ray_pole(), and Scalar::val_point().
ostream& operator<< | ( | ostream & | , | |
const Star & | ||||
) | [friend, inherited] |
Display.
Scalar Star_rot::a_car [protected] |
Square of the metric factor A.
Definition at line 97 of file star_rot.h.
Scalar Star_rot::ak_car [protected] |
Scalar .
For axisymmetric stars, this quantity is related to the derivatives of by
In particular it is related to the quantities and
introduced by Eqs.~(3.7) and (3.8) of Bonazzola et al. Astron. Astrophys. 278 , 421 (1993) by
Definition at line 179 of file star_rot.h.
Scalar Star_rot::b_car [protected] |
Square of the metric factor B.
Definition at line 103 of file star_rot.h.
Scalar Star_rot::bbb [protected] |
Metric factor B.
Definition at line 100 of file star_rot.h.
Vector Star::beta [protected, inherited] |
Scalar Star_rot::dzeta [protected] |
Metric potential .
Definition at line 127 of file star_rot.h.
Scalar Star::ener [protected, inherited] |
Scalar Star::ener_euler [protected, inherited] |
Scalar Star::gam_euler [protected, inherited] |
Metric Star::gamma [protected, inherited] |
Scalar Star_rot::khi_shift [protected] |
Scalar used in the decomposition of
shift
, following Shibata's prescription [Prog.
Theor. Phys. 101 , 1199 (1999)] :
Definition at line 153 of file star_rot.h.
Scalar Star::logn [protected, inherited] |
Scalar Star::nbar [protected, inherited] |
Scalar Star_rot::nphi [protected] |
Metric coefficient .
Definition at line 106 of file star_rot.h.
Scalar Star_rot::nuf [protected] |
Part of the Metric potential =
logn
generated by the matter terms.
Definition at line 119 of file star_rot.h.
Scalar Star_rot::nuq [protected] |
Part of the Metric potential =
logn
generated by the quadratic terms.
Definition at line 124 of file star_rot.h.
int Star::nzet [protected, inherited] |
double Star_rot::omega [protected] |
Rotation angular velocity ([f_unit] ).
Definition at line 94 of file star_rot.h.
double* Star_rot::p_angu_mom [mutable, protected] |
Angular momentum.
Definition at line 225 of file star_rot.h.
double* Star_rot::p_aplat [mutable, protected] |
Flatening r_pole/r_eq.
Definition at line 230 of file star_rot.h.
double* Star_rot::p_espec_isco [mutable, protected] |
Specific energy of a particle on the ISCO.
Definition at line 238 of file star_rot.h.
double* Star_rot::p_f_eq [mutable, protected] |
Orbital frequency at the equator.
Definition at line 241 of file star_rot.h.
double* Star_rot::p_f_isco [mutable, protected] |
Orbital frequency of the ISCO.
Definition at line 236 of file star_rot.h.
double* Star_rot::p_grv2 [mutable, protected] |
Error on the virial identity GRV2.
Definition at line 227 of file star_rot.h.
double* Star_rot::p_grv3 [mutable, protected] |
Error on the virial identity GRV3.
Definition at line 228 of file star_rot.h.
Itbl* Star::p_l_surf [mutable, protected, inherited] |
double* Star_rot::p_lspec_isco [mutable, protected] |
Specific angular momentum of a particle on the ISCO.
Definition at line 240 of file star_rot.h.
double* Star::p_mass_b [mutable, protected, inherited] |
double* Star::p_mass_g [mutable, protected, inherited] |
double* Star_rot::p_mom_quad [mutable, protected] |
Quadrupole moment.
Definition at line 234 of file star_rot.h.
double* Star_rot::p_r_circ [mutable, protected] |
Circumferential radius.
Definition at line 229 of file star_rot.h.
double* Star_rot::p_r_isco [mutable, protected] |
Circumferential radius of the ISCO.
Definition at line 235 of file star_rot.h.
double* Star::p_ray_eq [mutable, protected, inherited] |
double* Star::p_ray_eq_3pis2 [mutable, protected, inherited] |
double* Star::p_ray_eq_pi [mutable, protected, inherited] |
double* Star::p_ray_eq_pis2 [mutable, protected, inherited] |
double* Star::p_ray_pole [mutable, protected, inherited] |
double* Star_rot::p_tsw [mutable, protected] |
Ratio T/W.
Definition at line 226 of file star_rot.h.
Tbl* Star::p_xi_surf [mutable, protected, inherited] |
double* Star_rot::p_z_eqb [mutable, protected] |
Backward redshift factor at equator.
Definition at line 232 of file star_rot.h.
double* Star_rot::p_z_eqf [mutable, protected] |
Forward redshift factor at equator.
Definition at line 231 of file star_rot.h.
double* Star_rot::p_z_pole [mutable, protected] |
Redshift factor at North pole.
Definition at line 233 of file star_rot.h.
Scalar Star::press [protected, inherited] |
bool Star_rot::relativistic [protected] |
Indicator of relativity: true
for a relativistic star, false
for a Newtonian one.
Definition at line 87 of file star_rot.h.
Scalar Star::s_euler [protected, inherited] |
Scalar Star_rot::ssjm1_dzeta [protected] |
Effective source at the previous step for the resolution of the Poisson equation for dzeta
.
Definition at line 196 of file star_rot.h.
Scalar Star_rot::ssjm1_khi [protected] |
Effective source at the previous step for the resolution of the Poisson equation for the scalar by means of
Map_et::poisson
.
is an intermediate quantity for the resolution of the elliptic equation for the shift vector
Definition at line 209 of file star_rot.h.
Scalar Star_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 185 of file star_rot.h.
Scalar Star_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 191 of file star_rot.h.
Scalar Star_rot::ssjm1_tggg [protected] |
Effective source at the previous step for the resolution of the Poisson equation for tggg
.
Definition at line 201 of file star_rot.h.
Vector Star_rot::ssjm1_wshift [protected] |
Effective source at the previous step for the resolution of the vector Poisson equation for .
is an intermediate quantity for the resolution of the elliptic equation for the shift vector
(Components with respect to the Cartesian triad associated with the mapping
mp
)
Definition at line 218 of file star_rot.h.
Sym_tensor Star::stress_euler [protected, inherited] |
Scalar Star_rot::tggg [protected] |
Metric potential .
Definition at line 130 of file star_rot.h.
Sym_tensor Star_rot::tkij [protected] |
Tensor related to the extrinsic curvature tensor by
.
tkij
contains the Cartesian components of .
Definition at line 160 of file star_rot.h.
Scalar Star_rot::tnphi [protected] |
Component of the shift vector.
Definition at line 111 of file star_rot.h.
Vector Star::u_euler [protected, inherited] |
double Star_rot::unsurc2 [protected] |
:
unsurc2=1
for a relativistic star, 0 for a Newtonian one.
Definition at line 92 of file star_rot.h.
Scalar Star_rot::uuu [protected] |
Norm of u_euler
.
Definition at line 114 of file star_rot.h.
Vector Star_rot::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 143 of file star_rot.h.