LORENE
Lorene::Star_rot Class Reference

Class for isolated rotating stars. More...

#include <star_rot.h>

Inheritance diagram for Lorene::Star_rot:
Lorene::Star Lorene::Gravastar

Public Member Functions

 Star_rot (Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
 Standard constructor. More...
 
 Star_rot (const Star_rot &)
 Copy constructor. More...
 
 Star_rot (Map &mp_i, const Eos &eos_i, FILE *fich)
 Constructor from a file (see sauve(FILE*) ). More...
 
virtual ~Star_rot ()
 Destructor. More...
 
void operator= (const Star_rot &)
 Assignment to another Star_rot. More...
 
bool is_relativistic () const
 Returns true for a relativistic star, false for a Newtonian one. More...
 
virtual double get_omega_c () const
 Returns the central value of the rotation angular velocity ([f_unit] ) More...
 
const Scalarget_bbb () const
 Returns the metric factor B. More...
 
const Scalarget_a_car () const
 Returns the square of the metric factor A. More...
 
const Scalarget_b_car () const
 Returns the square of the metric factor B. More...
 
const Scalarget_nphi () const
 Returns the metric coefficient $N^\varphi$. More...
 
const Scalarget_tnphi () const
 Returns the component $\tilde N^\varphi = N^\varphi r\sin\theta$ of the shift vector. More...
 
const Scalarget_uuu () const
 Returns the norm of u_euler. More...
 
const Scalarget_nuf () const
 Returns the part of the Metric potential $\nu = \ln N$ = logn generated by the matter terms. More...
 
const Scalarget_nuq () const
 Returns the Part of the Metric potential $\nu = \ln N$ = logn generated by the quadratic terms. More...
 
const Scalarget_dzeta () const
 Returns the Metric potential $\zeta = \ln(AN)$. More...
 
const Scalarget_tggg () const
 Returns the Metric potential $\tilde G = (NB-1) r\sin\theta$. More...
 
const Vectorget_w_shift () const
 Returns the vector $W^i$ used in the decomposition of shift , following Shibata's prescription [Prog. More...
 
const Scalarget_khi_shift () const
 Returns the scalar $\chi$ used in the decomposition of shift
following Shibata's prescription [Prog. More...
 
const Sym_tensorget_tkij () const
 Returns the tensor ${\tilde K_{ij}}$ related to the extrinsic curvature tensor by ${\tilde K_{ij}} = B^{-2} K_{ij}$. More...
 
const Scalarget_ak_car () const
 Returns the scalar $A^2 K_{ij} K^{ij}$. More...
 
virtual void sauve (FILE *) const
 Save in a file. More...
 
virtual void display_poly (ostream &) const
 Display in polytropic units. More...
 
virtual const Itbll_surf () const
 Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$. More...
 
virtual double mass_b () const
 Baryon mass. More...
 
virtual double mass_g () const
 Gravitational mass. More...
 
virtual double angu_mom () const
 Angular momentum. More...
 
virtual double tsw () const
 Ratio T/W. More...
 
virtual double grv2 () const
 Error on the virial identity GRV2. More...
 
virtual double grv3 (ostream *ost=0x0) const
 Error on the virial identity GRV3. More...
 
virtual double r_circ () const
 Circumferential equatorial radius. More...
 
virtual double r_circ_merid () const
 Crcumferential meridional radius. More...
 
virtual double aplat () const
 Flatening r_pole/r_eq. More...
 
virtual double area () const
 Integrated surface area in ${\rm km}^2$. More...
 
virtual double mean_radius () const
 Mean star radius from the area $ r_{\rm mean} = \sqrt{\mathcal{A}} / 4\pi$. 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 const Tblsurf_grav () const
 Surface gravity (table along the theta direction) 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_eq () const
 Orbital frequency at the equator. More...
 
virtual void hydro_euler ()
 Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame. More...
 
void update_metric ()
 Computes metric coefficients from known potentials. More...
 
void fait_shift ()
 Computes shift from w_shift and khi_shift according to Shibata's prescription [Prog. More...
 
void fait_nphi ()
 Computes tnphi and nphi from the Cartesian components of the shift, stored in shift . More...
 
void extrinsic_curvature ()
 Computes tkij and ak_car from shift , nnn and b_car . More...
 
virtual void equilibrium (double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff, Param *=0x0)
 Computes an equilibrium configuration. More...
 
Mapset_mp ()
 Read/write of the mapping. More...
 
void set_enthalpy (const Scalar &)
 Assignment of the enthalpy field. More...
 
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 *pent_limit=0x0)
 Computes a spherical static configuration. More...
 
const Mapget_mp () const
 Returns the mapping. More...
 
int get_nzet () const
 Returns the number of domains occupied by the star. More...
 
const Eosget_eos () const
 Returns the equation of state. More...
 
const Scalarget_ent () const
 Returns the enthalpy field. More...
 
const Scalarget_nbar () const
 Returns the proper baryon density. More...
 
const Scalarget_ener () const
 Returns the proper total energy density. More...
 
const Scalarget_press () const
 Returns the fluid pressure. More...
 
const Scalarget_ener_euler () const
 Returns the total energy density with respect to the Eulerian observer. More...
 
const Scalarget_s_euler () const
 Returns the trace of the stress tensor in the Eulerian frame. More...
 
const Scalarget_gam_euler () const
 Returns the Lorentz factor between the fluid and Eulerian observers. More...
 
const Vectorget_u_euler () const
 Returns the fluid 3-velocity with respect to the Eulerian observer. More...
 
const Tensorget_stress_euler () const
 Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer. More...
 
const Scalarget_logn () const
 Returns the logarithm of the lapse N. More...
 
const Scalarget_nn () const
 Returns the lapse function N. More...
 
const Vectorget_beta () const
 Returns the shift vector $\beta^i$. More...
 
const Scalarget_lnq () const
 
const Metricget_gamma () const
 Returns the 3-metric $\gamma$. More...
 
double ray_eq () const
 Coordinate radius at $\phi=0$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_eq_pis2 () const
 Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_eq_pi () const
 Coordinate radius at $\phi=\pi$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_eq_3pis2 () const
 Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$ [r_unit]. More...
 
double ray_pole () const
 Coordinate radius at $\theta=0$ [r_unit]. More...
 
const Tblxi_surf () const
 Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$. More...
 

Static Public Member Functions

static double lambda_grv2 (const Scalar &sou_m, const Scalar &sou_q)
 Computes the coefficient $\lambda$ which ensures that the GRV2 virial identity is satisfied. More...
 

Protected Member Functions

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

Protected Attributes

bool relativistic
 Indicator of relativity: true for a relativistic star, false for a Newtonian one. More...
 
double unsurc2
 $1/c^2$ : unsurc2=1 for a relativistic star, 0 for a Newtonian one. More...
 
double omega
 Rotation angular velocity ([f_unit] ) More...
 
Scalar a_car
 Square of the metric factor A. More...
 
Scalar bbb
 Metric factor B. More...
 
Scalar b_car
 Square of the metric factor B. More...
 
Scalar nphi
 Metric coefficient $N^\varphi$. More...
 
Scalar tnphi
 Component $\tilde N^\varphi = N^\varphi r\sin\theta$ of the shift vector. More...
 
Scalar uuu
 Norm of u_euler. More...
 
Scalar nuf
 Part of the Metric potential $\nu = \ln N$ = logn generated by the matter terms. More...
 
Scalar nuq
 Part of the Metric potential $\nu = \ln N$ = logn generated by the quadratic terms. More...
 
Scalar dzeta
 Metric potential $\zeta = \ln(AN)$. More...
 
Scalar tggg
 Metric potential $\tilde G = (NB-1) r\sin\theta$. More...
 
Vector w_shift
 Vector $W^i$ used in the decomposition of shift , following Shibata's prescription [Prog. More...
 
Scalar khi_shift
 Scalar $\chi$ used in the decomposition of shift , following Shibata's prescription [Prog. More...
 
Sym_tensor tkij
 Tensor ${\tilde K_{ij}}$ related to the extrinsic curvature tensor by ${\tilde K_{ij}} = B^{-2} K_{ij}$. More...
 
Scalar ak_car
 Scalar $A^2 K_{ij} K^{ij}$. More...
 
Scalar ssjm1_nuf
 Effective source at the previous step for the resolution of the Poisson equation for nuf by means of Map_et::poisson . More...
 
Scalar ssjm1_nuq
 Effective source at the previous step for the resolution of the Poisson equation for nuq by means of Map_et::poisson . More...
 
Scalar ssjm1_dzeta
 Effective source at the previous step for the resolution of the Poisson equation for dzeta . More...
 
Scalar ssjm1_tggg
 Effective source at the previous step for the resolution of the Poisson equation for tggg . More...
 
Scalar ssjm1_khi
 Effective source at the previous step for the resolution of the Poisson equation for the scalar $\chi$ by means of Map_et::poisson . More...
 
Vector ssjm1_wshift
 Effective source at the previous step for the resolution of the vector Poisson equation for $W^i$. More...
 
double * p_angu_mom
 Angular momentum. More...
 
double * p_tsw
 Ratio T/W. More...
 
double * p_grv2
 Error on the virial identity GRV2. More...
 
double * p_grv3
 Error on the virial identity GRV3. More...
 
double * p_r_circ
 Circumferential radius (equator) More...
 
double * p_r_circ_merid
 Circumferential radius (meridian) More...
 
double * p_aplat
 Flatening r_pole/r_eq. More...
 
double * p_area
 Integrated surface area. 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_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...
 
Tblp_surf_grav
 Surface gravity (along the theta direction) More...
 
Mapmp
 Mapping associated with the star. More...
 
int nzet
 Number of domains of *mp occupied by the star. More...
 
const Eoseos
 Equation of state of the stellar matter. More...
 
Scalar ent
 Log-enthalpy. More...
 
Scalar nbar
 Baryon density in the fluid frame. More...
 
Scalar ener
 Total energy density in the fluid frame. More...
 
Scalar press
 Fluid pressure. More...
 
Scalar ener_euler
 Total energy density in the Eulerian frame. More...
 
Scalar s_euler
 Trace of the stress scalar in the Eulerian frame. More...
 
Scalar gam_euler
 Lorentz factor between the fluid and Eulerian observers. More...
 
Vector u_euler
 Fluid 3-velocity with respect to the Eulerian observer. More...
 
Sym_tensor stress_euler
 Spatial part of the stress-energy tensor with respect to the Eulerian observer. More...
 
Scalar logn
 Logarithm of the lapse N . More...
 
Scalar nn
 Lapse function N . More...
 
Vector beta
 Shift vector. More...
 
Scalar lnq
 
Metric gamma
 3-metric More...
 
double * p_ray_eq
 Coordinate radius at $\phi=0$, $\theta=\pi/2$. More...
 
double * p_ray_eq_pis2
 Coordinate radius at $\phi=\pi/2$, $\theta=\pi/2$. More...
 
double * p_ray_eq_pi
 Coordinate radius at $\phi=\pi$, $\theta=\pi/2$. More...
 
double * p_ray_eq_3pis2
 Coordinate radius at $\phi=3\pi/2$, $\theta=\pi/2$. More...
 
double * p_ray_pole
 Coordinate radius at $\theta=0$. More...
 
Itblp_l_surf
 Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$. More...
 
Tblp_xi_surf
 Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$. More...
 
double * p_mass_b
 Baryon mass. More...
 
double * p_mass_g
 Gravitational mass. More...
 

Friends

class Star
 

Detailed Description

Class for isolated rotating stars.

()

The metric is

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

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 $N^\varphi$ is set to zero; the only relevant gravitational quantity in this case is logn which represents the Newtonian gravitational potential generated by the star.

Version
Id
star_rot.h,v 1.10 2023/07/04 09:01:09 j_novak Exp

Definition at line 98 of file star_rot.h.

Constructor & Destructor Documentation

◆ Star_rot() [1/3]

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

Standard constructor.

Parameters
mp_iMapping on which the star is contructed
nzet_iNumber of domains occupied by the star
relattrue for a relativistic star, false for a Newtonian one
eos_iEquation of state of the stellar matter

Definition at line 97 of file star_rot.C.

References a_car, ak_car, b_car, bbb, Lorene::Star::beta, dzeta, Lorene::Map::get_bvect_cart(), khi_shift, Lorene::Star::mp, nphi, nuf, nuq, omega, relativistic, set_der_0x0(), Lorene::Tensor::set_etat_zero(), Lorene::Tensor::set_triad(), ssjm1_dzeta, ssjm1_khi, ssjm1_nuf, ssjm1_nuq, ssjm1_tggg, ssjm1_wshift, Lorene::Scalar::std_spectral_base(), tggg, tkij, tnphi, unsurc2, uuu, and w_shift.

◆ Star_rot() [2/3]

Lorene::Star_rot::Star_rot ( const Star_rot et)

Copy constructor.

Definition at line 167 of file star_rot.C.

References set_der_0x0().

◆ Star_rot() [3/3]

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

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

Parameters
mp_iMapping on which the star is constructed
eos_iEquation of state of the stellar matter
fichinput file (must have been created by the function Star_rot::sauve )

Definition at line 200 of file star_rot.C.

References a_car, ak_car, b_car, bbb, dzeta, fait_nphi(), fait_shift(), Lorene::fread_be(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), khi_shift, Lorene::Star::mp, nuf, nuq, omega, relativistic, set_der_0x0(), Lorene::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()

Lorene::Star_rot::~Star_rot ( )
virtual

Destructor.

Definition at line 300 of file star_rot.C.

References del_deriv().

Member Function Documentation

◆ angu_mom()

◆ aplat()

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

Flatening r_pole/r_eq.

Definition at line 489 of file star_rot_global.C.

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

◆ area()

◆ del_deriv()

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

Deletes all the derived quantities.

Reimplemented from Lorene::Star.

Definition at line 310 of file star_rot.C.

References Lorene::Star::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_r_circ, p_r_circ_merid, p_r_isco, p_surf_grav, p_tsw, p_z_eqb, p_z_eqf, p_z_pole, and set_der_0x0().

◆ del_hydro_euler()

void Lorene::Star_rot::del_hydro_euler ( )
protectedvirtual

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

Reimplemented from Lorene::Star.

Definition at line 362 of file star_rot.C.

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

◆ display_poly()

◆ equation_of_state()

◆ equilibrium()

void Lorene::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.

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

Definition at line 77 of file star_rot_equil.C.

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

◆ equilibrium_spher()

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

Computes a spherical static configuration.

Parameters
ent_c[input] central value of the enthalpy
precis[input] threshold in the relative difference between the enthalpy fields of two consecutive steps to stop the iterative procedure (default value: 1.e-14)
ent_limit[input] : array of enthalpy values to be set at the boundaries between the domains; if set to 0x0 (default), the initial values will be kept.

Definition at line 101 of file star_equil_spher.C.

References a_car, Lorene::Map_et::adapt(), Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Param::add_tbl(), Lorene::Scalar::annule(), b_car, bbb, Lorene::diffrel(), Lorene::Scalar::dsdr(), Lorene::Map_af::dsdr(), dzeta, Lorene::Star::ener, Lorene::Star::ener_euler, Lorene::Star::ent, Lorene::Star::equation_of_state(), Lorene::exp(), Lorene::Star::gam_euler, Lorene::Star::gamma, Lorene::Map_af::get_alpha(), Lorene::Map_et::get_alpha(), Lorene::Map_af::get_beta(), Lorene::Map_et::get_beta(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map_af::homothetie(), Lorene::Scalar::integrale(), Lorene::Star::logn, Lorene::Star::mass_b(), Lorene::Star::mass_g(), Lorene::Star::mp, Lorene::Star::nn, Lorene::norme(), Lorene::Star::nzet, Lorene::Map_af::poisson(), Lorene::Star::press, Lorene::Star::s_euler, Lorene::Vector::set(), Lorene::Map_af::set_alpha(), Lorene::Map_af::set_beta(), Lorene::Scalar::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Tensor::set_etat_zero(), Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), Lorene::Star::u_euler, Lorene::Scalar::val_grid_point(), and Lorene::Map::val_r().

◆ espec_isco()

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

Energy of a particle on the ISCO.

Definition at line 292 of file star_rot_isco.C.

References p_espec_isco, and r_isco().

◆ extrinsic_curvature()

void Lorene::Star_rot::extrinsic_curvature ( )

Computes tkij and ak_car from shift , nnn and b_car .

Definition at line 113 of file star_rot_upmetr.C.

References Lorene::Scalar::dsdr(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Star::mp, nphi, Lorene::Tensor::set_etat_zero(), Lorene::Scalar::srdsdt(), and tkij.

◆ f_eq()

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

Orbital frequency at the equator.

Definition at line 310 of file star_rot_isco.C.

References p_f_eq, and r_isco().

◆ f_isco()

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

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

Definition at line 258 of file star_rot_isco.C.

References p_f_isco, and r_isco().

◆ fait_nphi()

void Lorene::Star_rot::fait_nphi ( )

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

Definition at line 779 of file star_rot.C.

References Lorene::Star::beta.

◆ fait_shift()

◆ get_a_car()

const Scalar& Lorene::Star_rot::get_a_car ( ) const
inline

Returns the square of the metric factor A.

Definition at line 333 of file star_rot.h.

References a_car.

◆ get_ak_car()

const Scalar& Lorene::Star_rot::get_ak_car ( ) const
inline

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

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

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

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

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

Definition at line 417 of file star_rot.h.

References ak_car.

◆ get_b_car()

const Scalar& Lorene::Star_rot::get_b_car ( ) const
inline

Returns the square of the metric factor B.

Definition at line 336 of file star_rot.h.

References b_car.

◆ get_bbb()

const Scalar& Lorene::Star_rot::get_bbb ( ) const
inline

Returns the metric factor B.

Definition at line 330 of file star_rot.h.

References bbb.

◆ get_beta()

const Vector& Lorene::Star::get_beta ( ) const
inlineinherited

Returns the shift vector $\beta^i$.

Definition at line 402 of file star.h.

References Lorene::Star::beta.

◆ get_dzeta()

const Scalar& Lorene::Star_rot::get_dzeta ( ) const
inline

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

Definition at line 360 of file star_rot.h.

References dzeta.

◆ get_ener()

const Scalar& Lorene::Star::get_ener ( ) const
inlineinherited

Returns the proper total energy density.

Definition at line 370 of file star.h.

References Lorene::Star::ener.

◆ get_ener_euler()

const Scalar& Lorene::Star::get_ener_euler ( ) const
inlineinherited

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

Definition at line 376 of file star.h.

References Lorene::Star::ener_euler.

◆ get_ent()

const Scalar& Lorene::Star::get_ent ( ) const
inlineinherited

Returns the enthalpy field.

Definition at line 364 of file star.h.

References Lorene::Star::ent.

◆ get_eos()

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

Returns the equation of state.

Definition at line 361 of file star.h.

References Lorene::Star::eos.

◆ get_gam_euler()

const Scalar& Lorene::Star::get_gam_euler ( ) const
inlineinherited

Returns the Lorentz factor between the fluid and Eulerian observers.

Definition at line 382 of file star.h.

References Lorene::Star::gam_euler.

◆ get_gamma()

const Metric& Lorene::Star::get_gamma ( ) const
inlineinherited

Returns the 3-metric $\gamma$.

Definition at line 409 of file star.h.

References Lorene::Star::gamma.

◆ get_khi_shift()

const Scalar& Lorene::Star_rot::get_khi_shift ( ) const
inline

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

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

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

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

Definition at line 391 of file star_rot.h.

References khi_shift.

◆ get_logn()

const Scalar& Lorene::Star::get_logn ( ) const
inlineinherited

Returns the logarithm of the lapse N.

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 396 of file star.h.

References Lorene::Star::logn.

◆ get_mp()

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

Returns the mapping.

Definition at line 355 of file star.h.

References Lorene::Star::mp.

◆ get_nbar()

const Scalar& Lorene::Star::get_nbar ( ) const
inlineinherited

Returns the proper baryon density.

Definition at line 367 of file star.h.

References Lorene::Star::nbar.

◆ get_nn()

const Scalar& Lorene::Star::get_nn ( ) const
inlineinherited

Returns the lapse function N.

Definition at line 399 of file star.h.

References Lorene::Star::nn.

◆ get_nphi()

const Scalar& Lorene::Star_rot::get_nphi ( ) const
inline

Returns the metric coefficient $N^\varphi$.

Definition at line 339 of file star_rot.h.

References nphi.

◆ get_nuf()

const Scalar& Lorene::Star_rot::get_nuf ( ) const
inline

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

Definition at line 352 of file star_rot.h.

References nuf.

◆ get_nuq()

const Scalar& Lorene::Star_rot::get_nuq ( ) const
inline

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

Definition at line 357 of file star_rot.h.

References nuq.

◆ get_nzet()

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

Returns the number of domains occupied by the star.

Definition at line 358 of file star.h.

References Lorene::Star::nzet.

◆ get_omega_c()

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

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

Definition at line 666 of file star_rot.C.

References omega.

◆ get_press()

const Scalar& Lorene::Star::get_press ( ) const
inlineinherited

Returns the fluid pressure.

Definition at line 373 of file star.h.

References Lorene::Star::press.

◆ get_s_euler()

const Scalar& Lorene::Star::get_s_euler ( ) const
inlineinherited

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

Definition at line 379 of file star.h.

References Lorene::Star::s_euler.

◆ get_stress_euler()

const Tensor& Lorene::Star::get_stress_euler ( ) const
inlineinherited

Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer.

Definition at line 390 of file star.h.

References Lorene::Star::stress_euler.

◆ get_tggg()

const Scalar& Lorene::Star_rot::get_tggg ( ) const
inline

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

Definition at line 363 of file star_rot.h.

References tggg.

◆ get_tkij()

const Sym_tensor& Lorene::Star_rot::get_tkij ( ) const
inline

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

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

Definition at line 398 of file star_rot.h.

References tkij.

◆ get_tnphi()

const Scalar& Lorene::Star_rot::get_tnphi ( ) const
inline

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

Definition at line 344 of file star_rot.h.

References tnphi.

◆ get_u_euler()

const Vector& Lorene::Star::get_u_euler ( ) const
inlineinherited

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

Definition at line 385 of file star.h.

References Lorene::Star::u_euler.

◆ get_uuu()

const Scalar& Lorene::Star_rot::get_uuu ( ) const
inline

Returns the norm of u_euler.

Definition at line 347 of file star_rot.h.

References uuu.

◆ get_w_shift()

const Vector& Lorene::Star_rot::get_w_shift ( ) const
inline

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

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

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

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

Definition at line 377 of file star_rot.h.

References w_shift.

◆ grv2()

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

Error on the virial identity GRV2.

This indicator is only valid for relativistic computations.

Definition at line 222 of file star_rot_global.C.

References a_car, ak_car, Lorene::Scalar::derive_cov(), Lorene::Star::ener_euler, Lorene::Map::flat_met_spher(), lambda_grv2(), Lorene::Star::logn, Lorene::Star::mp, Lorene::Star::nbar, p_grv2, Lorene::Star::press, relativistic, and uuu.

◆ grv3()

double Lorene::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.

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

Definition at line 253 of file star_rot_global.C.

References ak_car, bbb, Lorene::Scalar::derive_cov(), dzeta, Lorene::Map::flat_met_spher(), Lorene::Scalar::get_etat(), Lorene::log(), Lorene::Star::logn, Lorene::Star::mp, p_grv3, relativistic, Lorene::Scalar::srdsdt(), and Lorene::Scalar::std_spectral_base().

◆ hydro_euler()

void Lorene::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 Lorene::Star.

Definition at line 60 of file star_rot_hydro.C.

References Lorene::Tensor::annule_domain(), Lorene::Star::beta, Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Star::mp, Lorene::Star::nn, omega, Lorene::Vector::set(), Lorene::Tensor::set_etat_qcq(), Lorene::Tensor::set_triad(), Lorene::Vector::std_spectral_base(), Lorene::Star::u_euler, Lorene::Map::x, and Lorene::Map::y.

◆ is_relativistic()

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

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

Definition at line 322 of file star_rot.h.

References relativistic.

◆ l_surf()

const Itbl & Lorene::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 $(\theta', \phi')$.

The stellar surface is defined as the location where the enthalpy (member ent ) vanishes.

Reimplemented from Lorene::Star.

Definition at line 86 of file star_rot_global.C.

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

◆ lambda_grv2()

double Lorene::Star_rot::lambda_grv2 ( const Scalar sou_m,
const Scalar sou_q 
)
static

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

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

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

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

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

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

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

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

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

Definition at line 69 of file star_rot_lambda_grv2.C.

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

◆ lspec_isco()

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

Angular momentum of a particle on the ISCO.

Definition at line 275 of file star_rot_isco.C.

References p_lspec_isco, and r_isco().

◆ mass_b()

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

◆ mass_g()

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

◆ mean_radius()

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

Mean star radius from the area $ r_{\rm mean} = \sqrt{\mathcal{A}} / 4\pi$.

Definition at line 479 of file star_rot_global.C.

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

◆ mom_quad()

double Lorene::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 ${\bar I}_{zz}$ of the MTW (1973) reduced quadrupole moment ${\bar I}_{ij}$ by: $Q = -3/2 {\bar I}_{zz}$. Note that Q is the negative of the quadrupole moment defined by Laarakkers and Poisson, Astrophys. J. 512 , 282 (1999).

Definition at line 586 of file star_rot_global.C.

References a_car, ak_car, bbb, Lorene::Scalar::check_dzpuis(), Lorene::Scalar::derive_cov(), Lorene::Star::ener_euler, Lorene::Map::flat_met_spher(), Lorene::Scalar::get_etat(), Lorene::Scalar::inc_dzpuis(), Lorene::log(), Lorene::Star::logn, Lorene::Star::mp, Lorene::Scalar::mult_r(), Lorene::Star::nbar, p_mom_quad, relativistic, Lorene::Star::s_euler, and Lorene::Scalar::std_spectral_base().

◆ operator=()

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

◆ operator>>()

◆ partial_display()

◆ r_circ()

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

Circumferential equatorial radius.

Definition at line 360 of file star_rot_global.C.

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

◆ r_circ_merid()

◆ r_isco()

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

◆ ray_eq()

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

◆ ray_eq_3pis2()

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

◆ ray_eq_pi()

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

◆ ray_eq_pis2()

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

◆ ray_pole()

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

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

Definition at line 281 of file star_global.C.

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

◆ sauve()

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

◆ set_der_0x0()

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

Sets to 0x0 all the pointers on derived quantities.

Reimplemented from Lorene::Star.

Definition at line 337 of file star_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_r_circ, p_r_circ_merid, p_r_isco, p_surf_grav, p_tsw, p_z_eqb, p_z_eqf, p_z_pole, and Lorene::Star::set_der_0x0().

◆ set_enthalpy()

void Lorene::Star::set_enthalpy ( const Scalar ent_i)
inherited

Assignment of the enthalpy field.

Definition at line 382 of file star.C.

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

◆ set_mp()

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

Read/write of the mapping.

Definition at line 322 of file star.h.

References Lorene::Star::mp.

◆ surf_grav()

const Tbl & Lorene::Star_rot::surf_grav ( ) const
virtual

Surface gravity (table along the theta direction)

Definition at line 671 of file star_rot_global.C.

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

◆ tsw()

double Lorene::Star_rot::tsw ( ) const
virtual

◆ update_metric()

void Lorene::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 54 of file star_rot_upmetr.C.

References a_car, b_car, bbb, del_deriv(), Lorene::Scalar::div_rsint(), dzeta, Lorene::exp(), extrinsic_curvature(), Lorene::Star::gamma, Lorene::Map::get_bvect_spher(), Lorene::Star::logn, Lorene::Star::mp, Lorene::Star::nn, Lorene::Tensor::set(), Lorene::Scalar::std_spectral_base(), tggg, and unsurc2.

◆ xi_surf()

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

Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.

The stellar surface is defined as the location where the enthalpy (member ent) vanishes.

Definition at line 92 of file star_global.C.

References Lorene::Star::l_surf(), Lorene::Star::p_l_surf, and Lorene::Star::p_xi_surf.

◆ z_eqb()

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

Backward redshift factor at equator.

Definition at line 534 of file star_rot_global.C.

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

◆ z_eqf()

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

Forward redshift factor at equator.

Definition at line 506 of file star_rot_global.C.

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

◆ z_pole()

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

Redshift factor at North pole.

Definition at line 565 of file star_rot_global.C.

References Lorene::Star::nn, p_z_pole, Lorene::Star::ray_pole(), and Lorene::Scalar::val_point().

Member Data Documentation

◆ a_car

Scalar Lorene::Star_rot::a_car
protected

Square of the metric factor A.

Definition at line 116 of file star_rot.h.

◆ ak_car

Scalar Lorene::Star_rot::ak_car
protected

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

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

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

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

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

Definition at line 198 of file star_rot.h.

◆ b_car

Scalar Lorene::Star_rot::b_car
protected

Square of the metric factor B.

Definition at line 122 of file star_rot.h.

◆ bbb

Scalar Lorene::Star_rot::bbb
protected

Metric factor B.

Definition at line 119 of file star_rot.h.

◆ beta

Vector Lorene::Star::beta
protectedinherited

Shift vector.

Definition at line 228 of file star.h.

◆ dzeta

Scalar Lorene::Star_rot::dzeta
protected

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

Definition at line 146 of file star_rot.h.

◆ ener

Scalar Lorene::Star::ener
protectedinherited

Total energy density in the fluid frame.

Definition at line 193 of file star.h.

◆ ener_euler

Scalar Lorene::Star::ener_euler
protectedinherited

Total energy density in the Eulerian frame.

Definition at line 198 of file star.h.

◆ ent

Scalar Lorene::Star::ent
protectedinherited

Log-enthalpy.

Definition at line 190 of file star.h.

◆ eos

const Eos& Lorene::Star::eos
protectedinherited

Equation of state of the stellar matter.

Definition at line 185 of file star.h.

◆ gam_euler

Scalar Lorene::Star::gam_euler
protectedinherited

Lorentz factor between the fluid and Eulerian observers.

Definition at line 204 of file star.h.

◆ gamma

Metric Lorene::Star::gamma
protectedinherited

3-metric

Definition at line 235 of file star.h.

◆ khi_shift

Scalar Lorene::Star_rot::khi_shift
protected

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

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

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

Definition at line 172 of file star_rot.h.

◆ logn

Scalar Lorene::Star::logn
protectedinherited

Logarithm of the lapse N .

In the Newtonian case, this is the Newtonian gravitational potential (in units of $c^2$).

Definition at line 222 of file star.h.

◆ mp

Map& Lorene::Star::mp
protectedinherited

Mapping associated with the star.

Definition at line 180 of file star.h.

◆ nbar

Scalar Lorene::Star::nbar
protectedinherited

Baryon density in the fluid frame.

Definition at line 192 of file star.h.

◆ nn

Scalar Lorene::Star::nn
protectedinherited

Lapse function N .

Definition at line 225 of file star.h.

◆ nphi

Scalar Lorene::Star_rot::nphi
protected

Metric coefficient $N^\varphi$.

Definition at line 125 of file star_rot.h.

◆ nuf

Scalar Lorene::Star_rot::nuf
protected

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

Definition at line 138 of file star_rot.h.

◆ nuq

Scalar Lorene::Star_rot::nuq
protected

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

Definition at line 143 of file star_rot.h.

◆ nzet

int Lorene::Star::nzet
protectedinherited

Number of domains of *mp occupied by the star.

Definition at line 183 of file star.h.

◆ omega

double Lorene::Star_rot::omega
protected

Rotation angular velocity ([f_unit] )

Definition at line 113 of file star_rot.h.

◆ p_angu_mom

double* Lorene::Star_rot::p_angu_mom
mutableprotected

Angular momentum.

Definition at line 244 of file star_rot.h.

◆ p_aplat

double* Lorene::Star_rot::p_aplat
mutableprotected

Flatening r_pole/r_eq.

Definition at line 250 of file star_rot.h.

◆ p_area

double* Lorene::Star_rot::p_area
mutableprotected

Integrated surface area.

Definition at line 251 of file star_rot.h.

◆ p_espec_isco

double* Lorene::Star_rot::p_espec_isco
mutableprotected

Specific energy of a particle on the ISCO.

Definition at line 259 of file star_rot.h.

◆ p_f_eq

double* Lorene::Star_rot::p_f_eq
mutableprotected

Orbital frequency at the equator.

Definition at line 262 of file star_rot.h.

◆ p_f_isco

double* Lorene::Star_rot::p_f_isco
mutableprotected

Orbital frequency of the ISCO.

Definition at line 257 of file star_rot.h.

◆ p_grv2

double* Lorene::Star_rot::p_grv2
mutableprotected

Error on the virial identity GRV2.

Definition at line 246 of file star_rot.h.

◆ p_grv3

double* Lorene::Star_rot::p_grv3
mutableprotected

Error on the virial identity GRV3.

Definition at line 247 of file star_rot.h.

◆ p_l_surf

Itbl* Lorene::Star::p_l_surf
mutableprotectedinherited

Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surface at the collocation points in $(\theta', \phi')$.

Definition at line 260 of file star.h.

◆ p_lspec_isco

double* Lorene::Star_rot::p_lspec_isco
mutableprotected

Specific angular momentum of a particle on the ISCO.

Definition at line 261 of file star_rot.h.

◆ p_mass_b

double* Lorene::Star::p_mass_b
mutableprotectedinherited

Baryon mass.

Definition at line 268 of file star.h.

◆ p_mass_g

double* Lorene::Star::p_mass_g
mutableprotectedinherited

Gravitational mass.

Definition at line 269 of file star.h.

◆ p_mom_quad

double* Lorene::Star_rot::p_mom_quad
mutableprotected

Quadrupole moment.

Definition at line 255 of file star_rot.h.

◆ p_r_circ

double* Lorene::Star_rot::p_r_circ
mutableprotected

Circumferential radius (equator)

Definition at line 248 of file star_rot.h.

◆ p_r_circ_merid

double* Lorene::Star_rot::p_r_circ_merid
mutableprotected

Circumferential radius (meridian)

Definition at line 249 of file star_rot.h.

◆ p_r_isco

double* Lorene::Star_rot::p_r_isco
mutableprotected

Circumferential radius of the ISCO.

Definition at line 256 of file star_rot.h.

◆ p_ray_eq

double* Lorene::Star::p_ray_eq
mutableprotectedinherited

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

Definition at line 242 of file star.h.

◆ p_ray_eq_3pis2

double* Lorene::Star::p_ray_eq_3pis2
mutableprotectedinherited

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

Definition at line 251 of file star.h.

◆ p_ray_eq_pi

double* Lorene::Star::p_ray_eq_pi
mutableprotectedinherited

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

Definition at line 248 of file star.h.

◆ p_ray_eq_pis2

double* Lorene::Star::p_ray_eq_pis2
mutableprotectedinherited

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

Definition at line 245 of file star.h.

◆ p_ray_pole

double* Lorene::Star::p_ray_pole
mutableprotectedinherited

Coordinate radius at $\theta=0$.

Definition at line 254 of file star.h.

◆ p_surf_grav

Tbl* Lorene::Star_rot::p_surf_grav
mutableprotected

Surface gravity (along the theta direction)

Definition at line 263 of file star_rot.h.

◆ p_tsw

double* Lorene::Star_rot::p_tsw
mutableprotected

Ratio T/W.

Definition at line 245 of file star_rot.h.

◆ p_xi_surf

Tbl* Lorene::Star::p_xi_surf
mutableprotectedinherited

Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate $\xi$ on the surface at the collocation points in $(\theta', \phi')$.

Definition at line 266 of file star.h.

◆ p_z_eqb

double* Lorene::Star_rot::p_z_eqb
mutableprotected

Backward redshift factor at equator.

Definition at line 253 of file star_rot.h.

◆ p_z_eqf

double* Lorene::Star_rot::p_z_eqf
mutableprotected

Forward redshift factor at equator.

Definition at line 252 of file star_rot.h.

◆ p_z_pole

double* Lorene::Star_rot::p_z_pole
mutableprotected

Redshift factor at North pole.

Definition at line 254 of file star_rot.h.

◆ press

Scalar Lorene::Star::press
protectedinherited

Fluid pressure.

Definition at line 194 of file star.h.

◆ relativistic

bool Lorene::Star_rot::relativistic
protected

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

Definition at line 106 of file star_rot.h.

◆ s_euler

Scalar Lorene::Star::s_euler
protectedinherited

Trace of the stress scalar in the Eulerian frame.

Definition at line 201 of file star.h.

◆ ssjm1_dzeta

Scalar Lorene::Star_rot::ssjm1_dzeta
protected

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

Definition at line 215 of file star_rot.h.

◆ ssjm1_khi

Scalar Lorene::Star_rot::ssjm1_khi
protected

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

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

Definition at line 228 of file star_rot.h.

◆ ssjm1_nuf

Scalar Lorene::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 204 of file star_rot.h.

◆ ssjm1_nuq

Scalar Lorene::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 210 of file star_rot.h.

◆ ssjm1_tggg

Scalar Lorene::Star_rot::ssjm1_tggg
protected

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

Definition at line 220 of file star_rot.h.

◆ ssjm1_wshift

Vector Lorene::Star_rot::ssjm1_wshift
protected

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

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

Definition at line 237 of file star_rot.h.

◆ stress_euler

Sym_tensor Lorene::Star::stress_euler
protectedinherited

Spatial part of the stress-energy tensor with respect to the Eulerian observer.

Definition at line 212 of file star.h.

◆ tggg

Scalar Lorene::Star_rot::tggg
protected

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

Definition at line 149 of file star_rot.h.

◆ tkij

Sym_tensor Lorene::Star_rot::tkij
protected

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

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

Definition at line 179 of file star_rot.h.

◆ tnphi

Scalar Lorene::Star_rot::tnphi
protected

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

Definition at line 130 of file star_rot.h.

◆ u_euler

Vector Lorene::Star::u_euler
protectedinherited

Fluid 3-velocity with respect to the Eulerian observer.

Definition at line 207 of file star.h.

◆ unsurc2

double Lorene::Star_rot::unsurc2
protected

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

Definition at line 111 of file star_rot.h.

◆ uuu

Scalar Lorene::Star_rot::uuu
protected

Norm of u_euler.

Definition at line 133 of file star_rot.h.

◆ w_shift

Vector Lorene::Star_rot::w_shift
protected

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

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

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

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

Definition at line 162 of file star_rot.h.


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