LORENE
|
Class for isolated rotating stars *** DEPRECATED : use class Star_rot
instead ***.
More...
#include <etoile.h>
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 Tenseur & | get_bbb () const |
Returns the metric factor B. More... | |
const Tenseur & | get_b_car () const |
Returns the square of the metric factor B. More... | |
const Tenseur & | get_nphi () const |
Returns the metric coefficient . More... | |
const Tenseur & | get_tnphi () const |
Returns the component of the shift vector. More... | |
const Tenseur & | get_uuu () const |
Returns the norm of u_euler . More... | |
const Tenseur & | get_logn () const |
Returns the metric potential = logn_auto . More... | |
const Tenseur & | get_nuf () const |
Returns the part of the Metric potential = logn generated by the matter terms. More... | |
const Tenseur & | get_nuq () const |
Returns the Part of the Metric potential = logn generated by the quadratic terms. More... | |
const Tenseur & | get_dzeta () const |
Returns the Metric potential = beta_auto . More... | |
const Tenseur & | get_tggg () const |
Returns the Metric potential . More... | |
const Tenseur & | get_w_shift () const |
Returns the vector used in the decomposition of shift , following Shibata's prescription [Prog. More... | |
const Tenseur & | get_khi_shift () const |
Returns the scalar used in the decomposition of shift following Shibata's prescription [Prog. More... | |
const Tenseur_sym & | get_tkij () const |
Returns the tensor related to the extrinsic curvature tensor by . More... | |
const Tenseur & | get_ak_car () const |
Returns the scalar . More... | |
virtual void | sauve (FILE *) const |
Save in a file. More... | |
virtual void | display_poly (ostream &) const |
Display in polytropic units. More... | |
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 . 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... | |
Map & | set_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 Map & | get_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 Eos & | get_eos () const |
Returns the equation of state. More... | |
const Tenseur & | get_ent () const |
Returns the enthalpy field. More... | |
const Tenseur & | get_nbar () const |
Returns the proper baryon density. More... | |
const Tenseur & | get_ener () const |
Returns the proper total energy density. More... | |
const Tenseur & | get_press () const |
Returns the fluid pressure. More... | |
const Tenseur & | get_ener_euler () const |
Returns the total energy density with respect to the Eulerian observer. More... | |
const Tenseur & | get_s_euler () const |
Returns the trace of the stress tensor in the Eulerian frame. More... | |
const Tenseur & | get_gam_euler () const |
Returns the Lorentz factor between the fluid and Eulerian observers. More... | |
const Tenseur & | get_u_euler () const |
Returns the fluid 3-velocity with respect to the Eulerian observer. More... | |
const Tenseur & | get_logn_auto () const |
Returns the logarithm of the part of the lapse N generated principaly by the star. More... | |
const Tenseur & | get_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 Tenseur & | get_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 Tenseur & | get_d_logn_auto_div () const |
Returns the gradient of logn_auto_div . More... | |
const Tenseur & | get_beta_auto () const |
Returns the logarithm of the part of the product AN generated principaly by the star. More... | |
const Tenseur & | get_nnn () const |
Returns the total lapse function N. More... | |
const Tenseur & | get_shift () const |
Returns the total shift vector . More... | |
const Tenseur & | get_a_car () const |
Returns the total conformal factor . More... | |
double | ray_eq () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_eq (int kk) const |
Coordinate radius at , [r_unit]. More... | |
double | ray_eq_pis2 () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_eq_pi () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_eq_3pis2 () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_pole () const |
Coordinate radius at [r_unit]. More... | |
const Tbl & | xi_surf () const |
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 . More... | |
Static Public Member Functions | |
static double | lambda_grv2 (const Cmp &sou_m, const Cmp &sou_q) |
Computes the coefficient 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 . More... | |
Tenseur | tnphi |
Component of the shift vector. More... | |
Tenseur | uuu |
Norm of u_euler . More... | |
Tenseur & | logn |
Metric potential = logn_auto . More... | |
Tenseur | nuf |
Part of the Metric potential = logn generated by the matter terms. More... | |
Tenseur | nuq |
Part of the Metric potential = logn generated by the quadratic terms. More... | |
Tenseur & | dzeta |
Metric potential = beta_auto . More... | |
Tenseur | tggg |
Metric potential . More... | |
Tenseur | w_shift |
Vector used in the decomposition of shift , following Shibata's prescription [Prog. More... | |
Tenseur | khi_shift |
Scalar used in the decomposition of shift , following Shibata's prescription [Prog. More... | |
Tenseur_sym | tkij |
Tensor related to the extrinsic curvature tensor by . More... | |
Tenseur | ak_car |
Scalar . 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 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 . 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... | |
Map & | mp |
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 |
: 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 Eos & | eos |
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 . More... | |
double * | p_ray_eq |
Coordinate radius at , . More... | |
double * | p_ray_eq_pis2 |
Coordinate radius at , . More... | |
double * | p_ray_eq_pi |
Coordinate radius at , . More... | |
double * | p_ray_eq_3pis2 |
Coordinate radius at , . More... | |
double * | p_ray_pole |
Coordinate radius at . More... | |
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 . More... | |
Tbl * | p_xi_surf |
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the surface at the collocation points in . More... | |
double * | p_mass_b |
Baryon mass. More... | |
double * | p_mass_g |
Gravitational mass. More... | |
Standard constructor.
mp_i | Mapping on which the star will be defined |
nzet_i | Number of domains occupied by the star |
relat | should be true for a relativistic star, false for a Newtonian one |
eos_i | Equation of state of the stellar matter |
Definition at line 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.
Lorene::Etoile_rot::Etoile_rot | ( | const Etoile_rot & | et | ) |
Constructor from a file (see sauve(FILE*)
).
mp_i | Mapping on which the star will be defined |
eos_i | Equation of state of the stellar matter |
fich | input file (must have been created by the function sauve ) |
Definition at line 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.
|
virtual |
|
virtual |
Angular momentum.
Reimplemented in Lorene::Et_magnetisation, Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.
Definition at line 197 of file et_rot_global.C.
References Lorene::Etoile::a_car, b_car, Lorene::Etoile::ener_euler, Lorene::Cmp::integrale(), Lorene::Cmp::mult_r(), Lorene::Etoile::nbar, p_angu_mom, Lorene::Etoile::press, Lorene::Etoile::relativistic, Lorene::Cmp::std_base_scal(), uuu, and Lorene::Cmp::va.
|
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().
|
virtual |
Surface area.
Definition at line 476 of file et_rot_global.C.
References Lorene::Valeur::annule_hard(), Lorene::Valeur::c_cf, Lorene::Valeur::dsdt(), Lorene::Etoile::get_a_car(), Lorene::Mg3d::get_angu(), get_bbb(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), l_surf(), Lorene::Etoile::mp, Lorene::Valeur::mult_st(), p_area, Lorene::sqrt(), Lorene::Valeur::std_base_scal(), Lorene::Cmp::va, Lorene::Valeur::val_point_jk(), Lorene::Map_radial::val_r_jk(), and Lorene::Etoile::xi_surf().
|
protectedvirtual |
Deletes all the derived quantities.
Reimplemented from Lorene::Etoile.
Reimplemented in Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.
Definition at line 344 of file etoile_rot.C.
References Lorene::Etoile::del_deriv(), 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 set_der_0x0().
|
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().
|
virtual |
Display in polytropic units.
Reimplemented in Lorene::Et_rot_diff.
Definition at line 708 of file etoile_rot.C.
References angu_mom(), Lorene::Etoile::ener, Lorene::Etoile::eos, Lorene::Eos_poly::get_gam(), Lorene::Eos_poly::get_kap(), get_omega_c(), mass_b(), mass_g(), Lorene::Etoile::nbar, Lorene::pow(), r_circ(), Lorene::Etoile::ray_eq(), and Lorene::sqrt().
|
virtualinherited |
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Reimplemented in Lorene::Et_magnetisation, and Lorene::Et_rot_bifluid.
Definition at line 569 of file etoile.C.
References Lorene::Param::add_int(), Lorene::Cmp::allocate_all(), Lorene::Etoile::del_deriv(), Lorene::Etoile::ener, Lorene::Eos::ener_ent(), Lorene::Etoile::ent, Lorene::Etoile::eos, Lorene::Mg3d::get_grille3d(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Etoile::mp, Lorene::Etoile::nbar, Lorene::Eos::nbar_ent(), Lorene::Etoile::nzet, Lorene::Etoile::press, Lorene::Eos::press_ent(), Lorene::Mtbl::set(), Lorene::Tenseur::set(), Lorene::Cmp::set(), Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), Lorene::Cmp::std_base_scal(), Lorene::Mtbl::t, and Lorene::Grille3d::x.
|
virtualinherited |
Computes a spherical static configuration with the outer boundary condition at a finite radius.
ent_c | [input] central value of the enthalpy |
precis | [input] threshold in the relative difference between the enthalpy fields of two consecutive steps to stop the iterative procedure (default value: 1.e-14) |
Definition at line 60 of file etoile_eqsph_falloff.C.
References Lorene::Etoile::a_car, 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::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), 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::Etoile::press, Lorene::Etoile::relativistic, Lorene::Etoile::s_euler, Lorene::Tenseur::set(), 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().
|
inherited |
Computes a spherical static configuration.
The sources for Poisson equations are regularized by extracting analytical diverging parts.
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().
|
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,
|
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 :
|
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.
|
virtualinherited |
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 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().
|
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().
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.
|
virtual |
Computation of frequency of eccentric orbits.
ecc | eccentricity of the orbit |
periasrt | periastron of the orbit |
ost | output stream to give details of the computation; if set to 0x0 [default value], no details will be given. |
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().
|
virtual |
Orbital frequency at the equator.
Definition at line 322 of file et_rot_isco.C.
|
virtual |
Orbital frequency at the innermost stable circular orbit (ISCO).
Definition at line 270 of file et_rot_isco.C.
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.
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)] :
Definition at line 766 of file etoile_rot.C.
References Lorene::Tenseur::get_etat(), Lorene::Tenseur::gradient(), and khi_shift.
|
inlineinherited |
Returns the total conformal factor .
Definition at line 736 of file etoile.h.
References Lorene::Etoile::a_car.
|
inline |
|
inline |
|
inline |
|
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.
|
inlineinherited |
Returns the gradient of logn_auto_div
.
Definition at line 722 of file etoile.h.
References Lorene::Etoile::d_logn_auto_div.
|
inline |
|
inlineinherited |
Returns the proper total energy density.
Definition at line 682 of file etoile.h.
References Lorene::Etoile::ener.
|
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.
|
inlineinherited |
Returns the enthalpy field.
Definition at line 676 of file etoile.h.
References Lorene::Etoile::ent.
|
inlineinherited |
Returns the equation of state.
Definition at line 673 of file etoile.h.
References Lorene::Etoile::eos.
|
inlineinherited |
Returns the Lorentz factor between the fluid and Eulerian observers.
Definition at line 694 of file etoile.h.
References Lorene::Etoile::gam_euler.
|
inline |
|
inline |
|
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 ).
Definition at line 704 of file etoile.h.
References Lorene::Etoile::logn_auto.
|
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 ).
Definition at line 718 of file etoile.h.
References Lorene::Etoile::logn_auto_div.
|
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 ).
Definition at line 711 of file etoile.h.
References Lorene::Etoile::logn_auto_regu.
|
inlineinherited |
|
inlineinherited |
Returns the proper baryon density.
Definition at line 679 of file etoile.h.
References Lorene::Etoile::nbar.
|
inlineinherited |
Returns the total lapse function N.
Definition at line 730 of file etoile.h.
References Lorene::Etoile::nnn.
|
inline |
|
inline |
|
inline |
|
inlineinherited |
Returns the number of domains occupied by the star.
Definition at line 665 of file etoile.h.
References Lorene::Etoile::nzet.
|
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.
|
inlineinherited |
Returns the fluid pressure.
Definition at line 685 of file etoile.h.
References Lorene::Etoile::press.
|
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.
|
inlineinherited |
Returns the total shift vector .
Definition at line 733 of file etoile.h.
References Lorene::Etoile::shift.
|
inline |
|
inline |
|
inline |
|
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.
|
inline |
|
inline |
|
virtual |
Error on the virial identity GRV2.
This indicator is only valid for relativistic computations.
Reimplemented in Lorene::Et_magnetisation, Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.
Definition at line 263 of file et_rot_global.C.
References Lorene::Etoile::a_car, ak_car, Lorene::Etoile::ener_euler, Lorene::flat_scalar_prod(), Lorene::Tenseur::gradient_spher(), lambda_grv2(), logn, Lorene::Etoile::mp, Lorene::Etoile::nbar, p_grv2, Lorene::Etoile::press, Lorene::Etoile::relativistic, and uuu.
|
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. |
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().
|
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.
|
inlineinherited |
Returns true
for a relativistic star, false
for a Newtonian one.
Definition at line 670 of file etoile.h.
References Lorene::Etoile::relativistic.
|
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 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.
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 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().
|
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().
|
virtual |
Baryon mass.
Reimplemented from Lorene::Etoile.
Reimplemented in Lorene::Et_rot_bifluid.
Definition at line 134 of file et_rot_global.C.
References Lorene::Etoile::a_car, bbb, Lorene::Etoile::gam_euler, Lorene::Tenseur::get_etat(), Lorene::Cmp::integrale(), Lorene::Etoile::nbar, Lorene::Etoile::p_mass_b, Lorene::Etoile::relativistic, and Lorene::Cmp::std_base_scal().
|
virtual |
Gravitational mass.
Reimplemented from Lorene::Etoile.
Reimplemented in Lorene::Et_magnetisation, Lorene::Et_rot_bifluid, and Lorene::Et_rot_mag.
Definition at line 166 of file et_rot_global.C.
References Lorene::Etoile::a_car, bbb, Lorene::Etoile::ener_euler, mass_b(), Lorene::Etoile::nnn, Lorene::Etoile::p_mass_g, Lorene::Etoile::press, Lorene::Etoile::relativistic, Lorene::Etoile::s_euler, Lorene::Tenseur::set_std_base(), tnphi, and uuu.
|
virtual |
|
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, , where 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.
|
virtual |
Part of the quadrupole moment.
is defined as , 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().
|
virtual |
Part of the quadrupole moment.
This term is defined by Laarakkers and Poisson, Astrophys. J. 512 , 282 (1999). Note that 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().
void Lorene::Etoile_rot::operator= | ( | const Etoile_rot & | et | ) |
Assignment to another Etoile_rot
.
Definition at line 415 of file etoile_rot.C.
References ak_car, b_car, bbb, del_deriv(), khi_shift, nphi, nuf, nuq, omega, Lorene::Etoile::operator=(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, tggg, tkij, tnphi, uuu, and w_shift.
|
protectedvirtual |
Operator >> (virtual function called by the operator <<).
Reimplemented from Lorene::Etoile.
Reimplemented in Lorene::Et_magnetisation, Lorene::Et_rot_bifluid, Lorene::Et_rot_mag, and Lorene::Et_rot_diff.
Definition at line 477 of file etoile_rot.C.
References Lorene::Etoile::a_car, angu_mom(), aplat(), area(), b_car, Lorene::diffrel(), dzeta, f_isco(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), get_omega_c(), grv2(), grv3(), khi_shift, logn, mass_g(), mean_radius(), mom_quad(), Lorene::Etoile::mp, nphi, Lorene::Etoile::nzet, omega, Lorene::Etoile::operator>>(), Lorene::pow(), r_circ(), r_circ_merid(), r_isco(), Lorene::Etoile::ray_eq(), Lorene::sqrt(), tsw(), uuu, w_shift, z_eqb(), z_eqf(), and z_pole().
|
protectedvirtual |
Printing of some informations, excluding all global quantities.
Reimplemented in Lorene::Et_rot_bifluid.
Definition at line 645 of file etoile_rot.C.
References Lorene::Etoile::a_car, aplat(), b_car, dzeta, Lorene::Etoile::ener, Lorene::Etoile::ent, Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), get_omega_c(), logn, Lorene::Etoile::mp, Lorene::Etoile::nbar, Lorene::Etoile::nnn, nphi, Lorene::Etoile::nzet, Lorene::Etoile::press, r_circ(), r_circ_merid(), Lorene::Etoile::ray_eq(), and uuu.
|
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.
|
virtual |
Circumferential meridional radius.
Definition at line 425 of file et_rot_global.C.
References Lorene::Valeur::annule_hard(), Lorene::Valeur::dsdt(), Lorene::Etoile::get_a_car(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), l_surf(), Lorene::Etoile::mp, p_r_circ_merid, Lorene::sqrt(), Lorene::Valeur::std_base_scal(), Lorene::Cmp::std_base_scal(), Lorene::Cmp::va, Lorene::Valeur::val_point_jk(), Lorene::Map_radial::val_r_jk(), and Lorene::Etoile::xi_surf().
|
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 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().
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 123 of file etoile_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Etoile::mp, and Lorene::Etoile::p_ray_eq.
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 443 of file etoile_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), and Lorene::Etoile::mp.
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 338 of file etoile_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Etoile::mp, and Lorene::Etoile::p_ray_eq_3pis2.
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 259 of file etoile_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Etoile::mp, and Lorene::Etoile::p_ray_eq_pi.
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 172 of file etoile_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Etoile::mp, and Lorene::Etoile::p_ray_eq_pis2.
|
inherited |
Coordinate radius at [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.
|
virtual |
Save in a file.
Reimplemented from Lorene::Etoile.
Reimplemented in Lorene::Et_magnetisation, Lorene::Et_rot_bifluid, Lorene::Et_rot_mag, and Lorene::Et_rot_diff.
Definition at line 452 of file etoile_rot.C.
References Lorene::fwrite_be(), khi_shift, nuf, nuq, omega, Lorene::Etoile::sauve(), Lorene::Tenseur::sauve(), Lorene::Cmp::sauve(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, tggg, and w_shift.
|
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().
|
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().
|
inlineinherited |
|
virtual |
Ratio T/W.
Reimplemented in Lorene::Et_rot_mag, and Lorene::Et_rot_diff.
Definition at line 229 of file et_rot_global.C.
References Lorene::Etoile::a_car, angu_mom(), bbb, Lorene::Etoile::ener, Lorene::Etoile::gam_euler, Lorene::Cmp::integrale(), logn, mass_g(), Lorene::Etoile::nbar, omega, p_tsw, Lorene::Etoile::relativistic, and Lorene::Cmp::std_base_scal().
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.
|
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 104 of file etoile_global.C.
References Lorene::Etoile::l_surf(), Lorene::Etoile::p_l_surf, and Lorene::Etoile::p_xi_surf.
|
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.
|
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.
|
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().
|
protectedinherited |
|
protected |
|
protected |
|
protected |
|
protectedinherited |
|
protectedinherited |
|
protected |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protected |
|
protected |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protected |
|
protected |
|
protected |
|
protectedinherited |
|
protected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotectedinherited |
|
mutableprotected |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotected |
|
mutableprotectedinherited |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protected |
|
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
|
protected |
Effective source at the previous step for the resolution of the Poisson equation for nuf
by means of Map_et::poisson
.
|
protected |
Effective source at the previous step for the resolution of the Poisson equation for nuq
by means of Map_et::poisson
.
|
protected |
|
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
)
|
protected |
|
protected |
|
protected |
|
protectedinherited |
|
protectedinherited |
|
protected |
|
protected |