LORENE
Lorene::Star_bhns Class Reference

Class for stars in black hole-neutron star binaries. More...

#include <star_bhns.h>

Inheritance diagram for Lorene::Star_bhns:
Lorene::Star

Public Member Functions

 Star_bhns (Map &mp_i, int nzet_i, const Eos &eos_i, bool irrot_i)
 Standard constructor. More...
 
 Star_bhns (const Star_bhns &)
 Copy constructor. More...
 
 Star_bhns (Map &mp_i, const Eos &eos_i, FILE *fich)
 Constructor from a file (see sauve(FILE*) ) More...
 
virtual ~Star_bhns ()
 Destructor. More...
 
void operator= (const Star_bhns &)
 Assignment to another Star_bhns. More...
 
Scalarset_pot_centri ()
 Read/write the centrifugal potential. More...
 
Scalarset_lapconf_auto ()
 Read/write of the lapconf function generated by the neutron star. More...
 
Scalarset_lapconf_comp ()
 Read/write of the lapconf function generated by the companion black hole. More...
 
Vectorset_shift_auto ()
 Read/write of the shift vector generated by the neutron star. More...
 
Vectorset_shift_comp ()
 Read/write of the shift vector generated by the companion black hole. More...
 
Scalarset_confo_auto ()
 Read/write of the conformal factor generated by the neutron star. More...
 
Scalarset_confo_comp ()
 Read/write of the conformal factor generated by the companion black hole. More...
 
bool is_irrotational () const
 Returns true for an irrotational motion, false for a corotating one. More...
 
const Scalarget_psi0 () const
 Returns the non-translational part of the velocity potential. More...
 
const Vectorget_d_psi () const
 Returns the covariant derivative of the velocity potential (Spherical components with respect to the mapping of the star) More...
 
const Vectorget_wit_w () const
 Returns the spatial projection of the fluid 3-velocity with respect to the co-orbiting observer. More...
 
const Scalarget_loggam () const
 Returns the logarithm of the Lorentz factor between the fluid and the co-orbiting observer. More...
 
const Vectorget_bsn () const
 Returns the shift vector, divided by N , of the rotating coordinates, $\beta^i/N$. More...
 
const Scalarget_gam () const
 Returns the Lorentz factor gam. More...
 
const Scalarget_gam0 () const
 Returns the Lorentz factor gam0. More...
 
const Scalarget_pot_centri () const
 Returns the centrifugal potential. More...
 
const Scalarget_lapconf_auto () const
 Returns the part of the lapconf function generated by the star. More...
 
const Scalarget_lapconf_comp () const
 Returns the part of the lapconf function generated by the companion black hole. More...
 
const Scalarget_lapconf_tot () const
 Returns the total lapconf function. More...
 
const Scalarget_lapse_auto () const
 
const Scalarget_lapse_tot () const
 Returns the total lapse function. More...
 
const Vectorget_d_lapconf_auto () const
 Returns the derivative of the lapse function generated by the star. More...
 
const Vectorget_d_lapconf_comp () const
 Returns the derivative of the lapse function generated by the companion black hole. More...
 
const Vectorget_shift_auto () const
 Returns the part of the shift vector generated by the star. More...
 
const Vectorget_shift_comp () const
 Returns the part of the shift vector generated by the companion black hole. More...
 
const Vectorget_shift_tot () const
 Returns the part total shift vector. More...
 
const Tensorget_d_shift_auto () const
 Returns the derivative of the shift vector generated by the star. More...
 
const Tensorget_d_shift_comp () const
 Returns the derivative of the shift vector generated by the companion black hole. More...
 
const Scalarget_confo_auto () const
 Returns the part of the conformal factor generated by the star. More...
 
const Scalarget_confo_comp () const
 Returns the part of the conformal factor generated by the companion black hole. More...
 
const Scalarget_confo_tot () const
 Returns the total conformal factor. More...
 
const Vectorget_d_confo_auto () const
 Returns the derivative of the conformal factor generated by the star. More...
 
const Vectorget_d_confo_comp () const
 Returns the derivative of the conformal factor generated by the companion black hole. More...
 
const Scalarget_psi4 () const
 Returns the fourth power of the total conformal factor. More...
 
const Sym_tensorget_taij_auto () const
 Returns the part of the extrinsic curvature tensor $\tilde A^{ij}$ generated by the neutron star part. More...
 
const Scalarget_taij_quad_auto () const
 Returns the part of the scalar $\eta_{ik} \eta_{jl} A^{ij} A^{kl}$ generated by $A_{ij}^{\rm auto}$. More...
 
virtual void sauve (FILE *) const
 Save in a file. More...
 
virtual double mass_b () const
 Baryon mass. More...
 
virtual double mass_b_bhns (bool kerrschild, const double &mass_bh, const double &sepa) const
 
virtual double mass_g () const
 Gravitational mass. More...
 
virtual double mass_g_bhns () const
 
void hydro_euler_bhns (bool kerrschild, const double &mass_bh, const double &sepa)
 Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame, as well as wit_w and loggam . More...
 
void update_metric_bhns (const Hole_bhns &hole, const Star_bhns &star_prev, double relax)
 Computes metric coefficients from known potentials with relaxation when the companion is a black hole. More...
 
void update_met_der_comp_bhns (const Hole_bhns &hole)
 Computes derivative of metric quantities from the companion black hole. More...
 
void kinema_bhns (bool kerrschild, const double &mass_bh, const double &sepa, double omega, double x_rot, double y_rot)
 Computes the quantities bsn and pot_centri . More...
 
void fait_d_psi_bhns ()
 Computes the gradient of the total velocity potential $\psi$. More...
 
void extr_curv_bhns ()
 Computes taij_auto , taij_quad_auto from shift_auto , lapse_auto , confo_auto . More...
 
void equilibrium_bhns (double ent_c, const double &mass_bh, const double &sepa, bool kerrschild, int mer, int mermax_ns, int mermax_potvit, int mermax_poisson, int filter_r, int filter_r_s, int filter_p_s, double relax_poisson, double relax_potvit, double thres_adapt, double resize_ns, const Tbl &fact_resize, Tbl &diff)
 Computes an equilibrium configuration. More...
 
double velo_pot_bhns (const double &mass_bh, const double &sepa, bool kerrschild, int mermax, double precis, double relax)
 Computes the non-translational part of the velocity scalar potential $\psi0$ by solving the continuity equation. More...
 
double chi_rp (double radius, double phi)
 Sensitive indicator of the mass-shedding to the direction of $r$, $\theta=\pi/2$, $\phi$. More...
 
double radius_p (double phi)
 Radius of the star to the direction of $\theta=\pi/2$ and $\phi$. More...
 
double phi_min ()
 Azimuthal angle when the indicator of the mass-shedding takes its minimum chi_min. More...
 
double phi_local_min (double phi_ini)
 Azimuthal angle when the indicator of the mass-shedding takes its local minimum. More...
 
void relax_bhns (const Star_bhns &star_prev, double relax_ent, double relax_met, int mer, int fmer_met)
 Performs a relaxation on ent , lapse_auto , shift_auto , confo_auto . More...
 
void equil_spher_bhns (double ent_c, double precis)
 Computes a spherical 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 hydro_euler ()
 Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame (nbar, ener and press). 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...
 
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...
 
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...
 

Protected Member Functions

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

Protected Attributes

Map_af mp_aff
 Affine mapping for solving poisson's equations of metric quantities. More...
 
bool irrotational
 true for an irrotational star, false for a corotating one More...
 
Scalar psi0
 Scalar potential $\Psi_0$ of the non-translational part of fluid 4-velocity (in the irrotational case) More...
 
Vector d_psi
 Gradient of $\Psi$ (in the irrotational case) (Spherical components with respect to the mapping of the star) More...
 
Vector wit_w
 Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer. More...
 
Scalar loggam
 Logarithm of the Lorentz factor between the fluid and the co-orbiting observer. More...
 
Vector bsn
 3-vector shift, divided by N , of the rotating coordinates, $\beta^i/N$. More...
 
Scalar gam
 Lorentz factor between the fluid and the co-orbiting observer. More...
 
Scalar gam0
 Lorentz factor between the co-orbiting observer and the Eulerian one. More...
 
Scalar pot_centri
 Centrifugal potential. More...
 
Scalar lapconf_auto
 Lapconf function generated by the star. More...
 
Scalar lapconf_comp
 Lapconf function generated by the companion black hole. More...
 
Scalar lapconf_tot
 Total lapconf function. More...
 
Scalar lapse_auto
 Lapse function generated by the "star". More...
 
Scalar lapse_tot
 Total lapse function. More...
 
Vector d_lapconf_auto
 Derivative of the lapconf function generated by the star $ \partial_j \alpha $. More...
 
Vector d_lapconf_comp
 Derivative of the lapconf function generated by the companion black hole. More...
 
Vector shift_auto
 Shift vector generated by the star. More...
 
Vector shift_comp
 Shift vector generated by the companion black hole. More...
 
Vector shift_tot
 Total shift vector. More...
 
Tensor d_shift_auto
 Derivative of the shift vector generated by the star $ \eta^{ik} \partial_k \beta^j $. More...
 
Tensor d_shift_comp
 Derivative of the shift vector generated by the companion black hole. More...
 
Scalar confo_auto
 Conformal factor generated by the star. More...
 
Scalar confo_comp
 Conformal factor generated by the companion black hole. More...
 
Scalar confo_tot
 Total conformal factor. More...
 
Vector d_confo_auto
 Derivative of the conformal factor generated by the star $ \partial_j \psi $. More...
 
Vector d_confo_comp
 Derivative of the conformal factor generated by the companion black hole. More...
 
Scalar psi4
 Fourth power of the total conformal factor. More...
 
Sym_tensor taij_auto
 Part of the extrinsic curvature tensor $ A^{ij}$ generated by shift_auto , lapse_auto , and confo_auto . More...
 
Scalar taij_quad_auto
 Part of the scalar $\eta_{ik} \eta_{jl} A^{ij} A^{kl}$ generated by $A_{ij}^{\rm auto}$. More...
 
Metric_flat flat
 Flat metric defined on the mapping (Spherical components with respect to the mapping of the star ). More...
 
Scalar ssjm1_lapconf
 Effective source at the previous step for the resolution of the Poisson equation for lapconf_auto . More...
 
Scalar ssjm1_confo
 Effective source at the previous step for the resolution of the Poisson equation for confo_auto . 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$ by means of Map_et::poisson . More...
 
double * p_mass_b_bhns
 Baryon mass. More...
 
double * p_mass_g_bhns
 Gravitational mass. 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 Bin_bhns
 

Detailed Description

Class for stars in black hole-neutron star binaries.

()

Definition at line 59 of file star_bhns.h.

Constructor & Destructor Documentation

◆ Star_bhns() [1/3]

Lorene::Star_bhns::Star_bhns ( Map mp_i,
int  nzet_i,
const Eos eos_i,
bool  irrot_i 
)

Standard constructor.

Parameters
mp_iMapping on which the star will be defined
nzet_iNumber of domains occupied by the star
eos_iEquation of state of the stellar matter
irrot_ishould be true for an irrotational star, false for a corotating one

Definition at line 76 of file star_bhns.C.

References bsn, Lorene::Vector::change_triad(), confo_auto, confo_comp, confo_tot, d_confo_auto, d_confo_comp, d_lapconf_auto, d_lapconf_comp, d_psi, d_shift_auto, d_shift_comp, gam, gam0, Lorene::Map::get_bvect_cart(), lapconf_auto, lapconf_comp, lapconf_tot, lapse_auto, lapse_tot, loggam, pot_centri, psi0, psi4, set_der_0x0(), Lorene::Tensor::set_etat_zero(), shift_auto, shift_comp, shift_tot, ssjm1_confo, ssjm1_khi, ssjm1_lapconf, ssjm1_wshift, Lorene::Scalar::std_spectral_base(), taij_auto, taij_quad_auto, Lorene::Star::u_euler, and wit_w.

◆ Star_bhns() [2/3]

Lorene::Star_bhns::Star_bhns ( const Star_bhns star)

Copy constructor.

Definition at line 181 of file star_bhns.C.

References set_der_0x0().

◆ Star_bhns() [3/3]

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

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

Parameters
mp_iMapping on which the star will be defined
eos_iEquation of state of the stellar matter
fichinput file (must have been created by the function sauve )

Definition at line 224 of file star_bhns.C.

References bsn, Lorene::Vector::change_triad(), confo_comp, confo_tot, d_confo_auto, d_confo_comp, d_lapconf_auto, d_lapconf_comp, d_psi, d_shift_auto, d_shift_comp, gam, gam0, Lorene::Star::gam_euler, Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), irrotational, lapconf_comp, lapconf_tot, lapse_auto, lapse_tot, loggam, Lorene::Star::mp, pot_centri, psi0, psi4, set_der_0x0(), Lorene::Tensor::set_etat_zero(), shift_comp, shift_tot, Lorene::Scalar::std_spectral_base(), taij_auto, taij_quad_auto, Lorene::Star::u_euler, and wit_w.

◆ ~Star_bhns()

Lorene::Star_bhns::~Star_bhns ( )
virtual

Destructor.

Definition at line 332 of file star_bhns.C.

References del_deriv().

Member Function Documentation

◆ chi_rp()

double Lorene::Star_bhns::chi_rp ( double  radius,
double  phi 
)

Sensitive indicator of the mass-shedding to the direction of $r$, $\theta=\pi/2$, $\phi$.

Parameters
radius[input] Radial coordinate
phi[input] Azimuthal angle

Definition at line 64 of file star_bhns_chi.C.

References Lorene::Scalar::dsdr(), Lorene::Star::ent, Lorene::Star::ray_pole(), and Lorene::Scalar::val_point().

◆ del_deriv()

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

Deletes all the derived quantities.

Reimplemented from Lorene::Star.

Definition at line 344 of file star_bhns.C.

References Lorene::Star::del_deriv(), p_mass_b_bhns, p_mass_g_bhns, and set_der_0x0().

◆ del_hydro_euler()

void Lorene::Star::del_hydro_euler ( )
protectedvirtualinherited

◆ equation_of_state()

◆ equil_spher_bhns()

◆ equilibrium_bhns()

void Lorene::Star_bhns::equilibrium_bhns ( double  ent_c,
const double &  mass_bh,
const double &  sepa,
bool  kerrschild,
int  mer,
int  mermax_ns,
int  mermax_potvit,
int  mermax_poisson,
int  filter_r,
int  filter_r_s,
int  filter_p_s,
double  relax_poisson,
double  relax_potvit,
double  thres_adapt,
double  resize_ns,
const Tbl fact_resize,
Tbl diff 
)

Computes an equilibrium configuration.

Parameters
ent_c[input] Central enthalpy
mass_bh[input] BH mass in the background metric
sepa[input] Separation between NS and BH
kerrschildshould be true for a Kerr-Schild background, false for a Conformally flat one
mer[input] Number of iteration
mermax_ns[input] Maximum number of iteration steps
mermax_potvit[input] Maximum number of steps in Map_radial::poisson_compact
mermax_poisson[input] Maximum number of steps in poisson scalar
filter_r[input] No. points to be deleted for (r): lap,conf
filter_r_s[input] No. points to be deleted for (r): shift
filter_p_s[input] No. points to be deleted for (phi): shift
relax_poisson[input] Relaxation factor in poisson scalar
relax_potvit[input] Relaxation factor in Map_radial::poisson_compact
thres_adapt[input] Threshold on dH/dr for the adaptation of the mapping
resize_ns[input] Resize factor for the first shell
fact_resize[input] 1-D Tbl for the input of some factors : \ fact(0) : A resizing factor for the first shell
diff[output] 1-D Tbl for the storage of some error indicators :

Definition at line 73 of file star_bhns_equilibrium.C.

References Lorene::Map::adapt(), Lorene::Param::add_cmp_mod(), Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Param::add_tbl(), Lorene::Param::add_tenseur_mod(), Lorene::Scalar::annule(), chi_rp(), confo_auto, confo_tot, Lorene::Scalar::dsdr(), Lorene::Scalar::dsdx(), Lorene::Scalar::dsdy(), Lorene::Star::ener_euler, Lorene::Star::ent, Lorene::Star::equation_of_state(), Lorene::exp(), gam, gam0, Lorene::Map::get_bvect_cart(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map_et::homothetie(), hydro_euler_bhns(), Lorene::Scalar::inc_dzpuis(), irrotational, lapconf_auto, lapconf_comp, Lorene::log(), loggam, Lorene::Star::mp, Lorene::Star::nzet, phi_min(), pot_centri, Lorene::pow(), radius_p(), Lorene::Star::ray_eq_pi(), Lorene::Star::ray_pole(), Lorene::Map::reevaluate(), Lorene::Map::resize(), Lorene::Star::s_euler, Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_etat_qcq(), Lorene::Tensor::set_etat_qcq(), Lorene::Tenseur::set_etat_qcq(), Lorene::Cmp::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::smooth(), Lorene::sqrt(), ssjm1_confo, ssjm1_khi, ssjm1_lapconf, ssjm1_wshift, Lorene::Scalar::std_spectral_base(), taij_quad_auto, Lorene::Scalar::val_grid_point(), Lorene::Scalar::val_point(), Lorene::Map::val_r(), velo_pot_bhns(), Lorene::Map::x, and Lorene::Map::y.

◆ 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 Lorene::Star_rot::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(), Lorene::Star_rot::b_car, Lorene::Star_rot::bbb, Lorene::diffrel(), Lorene::Scalar::dsdr(), Lorene::Map_af::dsdr(), Lorene::Star_rot::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().

◆ extr_curv_bhns()

void Lorene::Star_bhns::extr_curv_bhns ( )

◆ fait_d_psi_bhns()

◆ 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_bsn()

const Vector& Lorene::Star_bhns::get_bsn ( ) const
inline

Returns the shift vector, divided by N , of the rotating coordinates, $\beta^i/N$.

(Spherical components with respect to the mapping of the star)

Definition at line 327 of file star_bhns.h.

References bsn.

◆ get_confo_auto()

const Scalar& Lorene::Star_bhns::get_confo_auto ( ) const
inline

Returns the part of the conformal factor generated by the star.

Definition at line 383 of file star_bhns.h.

References confo_auto.

◆ get_confo_comp()

const Scalar& Lorene::Star_bhns::get_confo_comp ( ) const
inline

Returns the part of the conformal factor generated by the companion black hole.

Definition at line 388 of file star_bhns.h.

References confo_comp.

◆ get_confo_tot()

const Scalar& Lorene::Star_bhns::get_confo_tot ( ) const
inline

Returns the total conformal factor.

Definition at line 391 of file star_bhns.h.

References confo_tot.

◆ get_d_confo_auto()

const Vector& Lorene::Star_bhns::get_d_confo_auto ( ) const
inline

Returns the derivative of the conformal factor generated by the star.

Definition at line 396 of file star_bhns.h.

References d_confo_auto.

◆ get_d_confo_comp()

const Vector& Lorene::Star_bhns::get_d_confo_comp ( ) const
inline

Returns the derivative of the conformal factor generated by the companion black hole.

Definition at line 401 of file star_bhns.h.

References d_confo_comp.

◆ get_d_lapconf_auto()

const Vector& Lorene::Star_bhns::get_d_lapconf_auto ( ) const
inline

Returns the derivative of the lapse function generated by the star.

Definition at line 356 of file star_bhns.h.

References d_lapconf_auto.

◆ get_d_lapconf_comp()

const Vector& Lorene::Star_bhns::get_d_lapconf_comp ( ) const
inline

Returns the derivative of the lapse function generated by the companion black hole.

Definition at line 361 of file star_bhns.h.

References d_lapconf_comp.

◆ get_d_psi()

const Vector& Lorene::Star_bhns::get_d_psi ( ) const
inline

Returns the covariant derivative of the velocity potential (Spherical components with respect to the mapping of the star)

Definition at line 310 of file star_bhns.h.

References d_psi.

◆ get_d_shift_auto()

const Tensor& Lorene::Star_bhns::get_d_shift_auto ( ) const
inline

Returns the derivative of the shift vector generated by the star.

Definition at line 375 of file star_bhns.h.

References d_shift_auto.

◆ get_d_shift_comp()

const Tensor& Lorene::Star_bhns::get_d_shift_comp ( ) const
inline

Returns the derivative of the shift vector generated by the companion black hole.

Definition at line 380 of file star_bhns.h.

References d_shift_comp.

◆ 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()

const Scalar& Lorene::Star_bhns::get_gam ( ) const
inline

Returns the Lorentz factor gam.

Definition at line 330 of file star_bhns.h.

References gam.

◆ get_gam0()

const Scalar& Lorene::Star_bhns::get_gam0 ( ) const
inline

Returns the Lorentz factor gam0.

Definition at line 333 of file star_bhns.h.

References gam0.

◆ 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_lapconf_auto()

const Scalar& Lorene::Star_bhns::get_lapconf_auto ( ) const
inline

Returns the part of the lapconf function generated by the star.

Definition at line 339 of file star_bhns.h.

References lapconf_auto.

◆ get_lapconf_comp()

const Scalar& Lorene::Star_bhns::get_lapconf_comp ( ) const
inline

Returns the part of the lapconf function generated by the companion black hole.

Definition at line 344 of file star_bhns.h.

References lapconf_comp.

◆ get_lapconf_tot()

const Scalar& Lorene::Star_bhns::get_lapconf_tot ( ) const
inline

Returns the total lapconf function.

Definition at line 347 of file star_bhns.h.

References lapconf_tot.

◆ get_lapse_tot()

const Scalar& Lorene::Star_bhns::get_lapse_tot ( ) const
inline

Returns the total lapse function.

Definition at line 353 of file star_bhns.h.

References lapse_tot.

◆ get_loggam()

const Scalar& Lorene::Star_bhns::get_loggam ( ) const
inline

Returns the logarithm of the Lorentz factor between the fluid and the co-orbiting observer.

Definition at line 321 of file star_bhns.h.

References loggam.

◆ 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_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_pot_centri()

const Scalar& Lorene::Star_bhns::get_pot_centri ( ) const
inline

Returns the centrifugal potential.

Definition at line 336 of file star_bhns.h.

References pot_centri.

◆ 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_psi0()

const Scalar& Lorene::Star_bhns::get_psi0 ( ) const
inline

Returns the non-translational part of the velocity potential.

Definition at line 305 of file star_bhns.h.

References psi0.

◆ get_psi4()

const Scalar& Lorene::Star_bhns::get_psi4 ( ) const
inline

Returns the fourth power of the total conformal factor.

Definition at line 404 of file star_bhns.h.

References psi4.

◆ 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_shift_auto()

const Vector& Lorene::Star_bhns::get_shift_auto ( ) const
inline

Returns the part of the shift vector generated by the star.

Definition at line 364 of file star_bhns.h.

References shift_auto.

◆ get_shift_comp()

const Vector& Lorene::Star_bhns::get_shift_comp ( ) const
inline

Returns the part of the shift vector generated by the companion black hole.

Definition at line 369 of file star_bhns.h.

References shift_comp.

◆ get_shift_tot()

const Vector& Lorene::Star_bhns::get_shift_tot ( ) const
inline

Returns the part total shift vector.

Definition at line 372 of file star_bhns.h.

References shift_tot.

◆ 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_taij_auto()

const Sym_tensor& Lorene::Star_bhns::get_taij_auto ( ) const
inline

Returns the part of the extrinsic curvature tensor $\tilde A^{ij}$ generated by the neutron star part.

Definition at line 409 of file star_bhns.h.

References taij_auto.

◆ get_taij_quad_auto()

const Scalar& Lorene::Star_bhns::get_taij_quad_auto ( ) const
inline

Returns the part of the scalar $\eta_{ik} \eta_{jl} A^{ij} A^{kl}$ generated by $A_{ij}^{\rm auto}$.

Definition at line 415 of file star_bhns.h.

References taij_quad_auto.

◆ 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_wit_w()

const Vector& Lorene::Star_bhns::get_wit_w ( ) const
inline

Returns the spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.

(Spherical components with respect to the mapping of the star)

Definition at line 316 of file star_bhns.h.

References wit_w.

◆ hydro_euler()

void Lorene::Star::hydro_euler ( )
virtualinherited

Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame (nbar, ener and press).

Reimplemented in Lorene::Star_bin_xcts, Lorene::Star_bin, Lorene::Star_rot, Lorene::Star_rot_Dirac, Lorene::Star_rot_CFC, and Lorene::Star_rot_Dirac_diff.

Definition at line 580 of file star.C.

◆ hydro_euler_bhns()

void Lorene::Star_bhns::hydro_euler_bhns ( bool  kerrschild,
const double &  mass_bh,
const double &  sepa 
)

Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid frame, as well as wit_w and loggam .

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

Parameters
kerrschildshould be true for a Kerr-Schild background, false for a Conformally flat one
mass_bhBH mass in the background metric
sepaSeparation between NS and BH

Definition at line 62 of file star_bhns_hydro.C.

References Lorene::Tensor::annule_domain(), bsn, d_psi, del_deriv(), Lorene::Star::ener, Lorene::Star::ener_euler, Lorene::Star::ent, Lorene::exp(), gam, gam0, Lorene::Star::gam_euler, Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_y(), Lorene::Tensor::get_triad(), irrotational, Lorene::log(), loggam, Lorene::Star::mp, Lorene::Star::press, psi4, Lorene::Star::s_euler, Lorene::Vector::set(), Lorene::Tensor::set(), Lorene::Scalar::set_dzpuis(), Lorene::Tensor::set_etat_qcq(), Lorene::Tensor::set_etat_zero(), Lorene::sqrt(), Lorene::Vector::std_spectral_base(), Lorene::Scalar::std_spectral_base(), Lorene::Star::u_euler, wit_w, Lorene::Map::x, Lorene::Map::y, and Lorene::Map::z.

◆ is_irrotational()

bool Lorene::Star_bhns::is_irrotational ( ) const
inline

Returns true for an irrotational motion, false for a corotating one.

Definition at line 302 of file star_bhns.h.

References irrotational.

◆ kinema_bhns()

void Lorene::Star_bhns::kinema_bhns ( bool  kerrschild,
const double &  mass_bh,
const double &  sepa,
double  omega,
double  x_rot,
double  y_rot 
)

Computes the quantities bsn and pot_centri .

The calculation is performed starting from the quantities lapse_tot , shift_tot , which are supposed to be up to date.

Parameters
kerrschildshould be true for a Kerr-Schild background, false for a Conformally flat one
mass_bhBH mass in the background metric
sepaSeparation between NS and BH
omegaangular velocity with respect to an asymptotically inertial observer
x_rotabsolute X coordinate of the rotation axis
y_rotabsolute Y coordinate of the rotation axis

Definition at line 64 of file star_bhns_kinema.C.

References Lorene::Tensor::annule_domain(), bsn, confo_tot, del_deriv(), gam0, Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_y(), Lorene::Map::get_rot_phi(), lapconf_tot, Lorene::log(), Lorene::Star::mp, pot_centri, psi4, Lorene::Vector::set(), Lorene::Tensor::set(), Lorene::Tensor::set_etat_qcq(), shift_tot, Lorene::sqrt(), Lorene::Vector::std_spectral_base(), Lorene::Scalar::std_spectral_base(), Lorene::Map::x, Lorene::Map::xa, Lorene::Map::y, Lorene::Map::ya, and Lorene::Map::z.

◆ l_surf()

const Itbl & Lorene::Star::l_surf ( ) const
virtualinherited

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 in Lorene::Star_rot.

Definition at line 66 of file star_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.

◆ mass_b()

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

◆ mass_g()

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

Gravitational mass.

Implements Lorene::Star.

Definition at line 141 of file star_bhns_global.C.

◆ operator=()

◆ operator>>()

ostream & Lorene::Star_bhns::operator>> ( ostream &  ost) const
protectedvirtual

Operator >> (virtual function called by the operator <<).

Reimplemented from Lorene::Star.

Definition at line 498 of file star_bhns.C.

References confo_tot, Lorene::Star::ent, lapse_tot, Lorene::Star::ray_eq(), Lorene::Star::ray_eq_pi(), Lorene::Star::ray_eq_pis2(), Lorene::Star::ray_pole(), shift_tot, and Lorene::Scalar::val_grid_point().

◆ phi_local_min()

double Lorene::Star_bhns::phi_local_min ( double  phi_ini)

Azimuthal angle when the indicator of the mass-shedding takes its local minimum.

Parameters
phi_ini[input] Initial azumuthal angle to search minimum

Definition at line 148 of file star_bhns_chi.C.

References chi_rp(), and radius_p().

◆ phi_min()

double Lorene::Star_bhns::phi_min ( )

Azimuthal angle when the indicator of the mass-shedding takes its minimum chi_min.

Definition at line 86 of file star_bhns_chi.C.

References chi_rp(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Star::mp, phi_local_min(), radius_p(), Lorene::Tbl::set(), and Lorene::Tbl::set_etat_qcq().

◆ radius_p()

double Lorene::Star_bhns::radius_p ( double  phi)

Radius of the star to the direction of $\theta=\pi/2$ and $\phi$.

Parameters
phi[input] Azimuthal angle

Definition at line 77 of file star_bhns_chi.C.

References Lorene::Star::mp, Lorene::Star::nzet, and Lorene::Map::val_r().

◆ 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.

◆ relax_bhns()

void Lorene::Star_bhns::relax_bhns ( const Star_bhns star_prev,
double  relax_ent,
double  relax_met,
int  mer,
int  fmer_met 
)

Performs a relaxation on ent , lapse_auto , shift_auto , confo_auto .

Parameters
star_prev[input] star at the previous step
relax_ent[input] Relaxation factor for ent
relax_met[input] Relaxation factor for lapse_auto , shift_auto , confo_auto , only if (merfmer_met=0) .
mer[input] Step number
fmer_met[input] Step interval between metric updates

Definition at line 597 of file star_bhns.C.

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

◆ sauve()

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

◆ set_confo_auto()

Scalar & Lorene::Star_bhns::set_confo_auto ( )

Read/write of the conformal factor generated by the neutron star.

Definition at line 452 of file star_bhns.C.

References confo_auto, and del_deriv().

◆ set_confo_comp()

Scalar & Lorene::Star_bhns::set_confo_comp ( )

Read/write of the conformal factor generated by the companion black hole.

Definition at line 459 of file star_bhns.C.

References confo_comp, and del_deriv().

◆ set_der_0x0()

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

Sets to 0x0 all the pointers on derived quantities.

Reimplemented from Lorene::Star.

Definition at line 355 of file star_bhns.C.

References p_mass_b_bhns, p_mass_g_bhns, 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_lapconf_auto()

Scalar & Lorene::Star_bhns::set_lapconf_auto ( )

Read/write of the lapconf function generated by the neutron star.

Definition at line 424 of file star_bhns.C.

References del_deriv(), and lapconf_auto.

◆ set_lapconf_comp()

Scalar & Lorene::Star_bhns::set_lapconf_comp ( )

Read/write of the lapconf function generated by the companion black hole.

Definition at line 431 of file star_bhns.C.

References del_deriv(), and lapconf_comp.

◆ 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.

◆ set_pot_centri()

Scalar & Lorene::Star_bhns::set_pot_centri ( )

Read/write the centrifugal potential.

Definition at line 417 of file star_bhns.C.

References del_deriv(), and pot_centri.

◆ set_shift_auto()

Vector & Lorene::Star_bhns::set_shift_auto ( )

Read/write of the shift vector generated by the neutron star.

Definition at line 438 of file star_bhns.C.

References del_deriv(), and shift_auto.

◆ set_shift_comp()

Vector & Lorene::Star_bhns::set_shift_comp ( )

Read/write of the shift vector generated by the companion black hole.

Definition at line 445 of file star_bhns.C.

References del_deriv(), and shift_comp.

◆ update_met_der_comp_bhns()

void Lorene::Star_bhns::update_met_der_comp_bhns ( const Hole_bhns hole)

Computes derivative of metric quantities from the companion black hole.

Parameters
holecompanion black hole

Definition at line 63 of file star_bhns_upmetr_der.C.

References Lorene::Hole_bhns::get_d_lapconf_auto().

◆ update_metric_bhns()

void Lorene::Star_bhns::update_metric_bhns ( const Hole_bhns hole,
const Star_bhns star_prev,
double  relax 
)

Computes metric coefficients from known potentials with relaxation when the companion is a black hole.

The calculation is performed starting from the quantities lapse_auto , shift_auto , confo_auto , comp.lapse_auto , comp.confo_auto which are supposed to be up to date. From these, the following fields are updated: lapse_comp , lapse_tot, confo_comp, confo_tot , psi4 ,

Parameters
holecompanion black hole
star_prevprevious value of the star
relaxrelaxation parameter

Definition at line 63 of file star_bhns_upmetr.C.

References Lorene::Hole_bhns::get_lapconf_auto(), and Lorene::Black_hole::get_mp().

◆ velo_pot_bhns()

double Lorene::Star_bhns::velo_pot_bhns ( const double &  mass_bh,
const double &  sepa,
bool  kerrschild,
int  mermax,
double  precis,
double  relax 
)

Computes the non-translational part of the velocity scalar potential $\psi0$ by solving the continuity equation.

Parameters
mass_bh[input] BH mass in the background metric
sepa[input] Separation between NS and BH
kerrschildshould be true for a Kerr-Schild background, false for a Conformally flat one
mermax[input] Maximum number of steps in the iteration
precis[input] Required precision: the iteration will be stopped when the relative difference on $\psi0$ between two successive steps is lower than precis .
relax[input] Relaxation factor.
Returns
Relative error of the resolution obtained by comparing the operator acting on the solution with the source.

Definition at line 73 of file star_bhns_vel_pot.C.

References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Scalar::annule(), bsn, Lorene::Vector::change_triad(), confo_tot, Lorene::contract(), d_confo_auto, d_confo_comp, d_lapconf_auto, d_lapconf_comp, Lorene::Eos::der_nbar_ent(), Lorene::Scalar::deriv(), Lorene::Scalar::derive_con(), Lorene::Scalar::derive_cov(), Lorene::Star::ent, Lorene::Star::eos, Lorene::exp(), flat, Lorene::Map::flat_met_spher(), Lorene::Star::gam_euler, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::get_triad(), lapconf_tot, Lorene::log(), Lorene::Star::mp, Lorene::Star::nzet, psi0, psi4, Lorene::Scalar::set_domain(), Lorene::Tensor::set_etat_qcq(), Lorene::Tensor::set_triad(), and Lorene::Scalar::std_spectral_base().

◆ 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.

Member Data Documentation

◆ beta

Vector Lorene::Star::beta
protectedinherited

Shift vector.

Definition at line 228 of file star.h.

◆ bsn

Vector Lorene::Star_bhns::bsn
protected

3-vector shift, divided by N , of the rotating coordinates, $\beta^i/N$.

(Spherical components with respect to the mapping of the star)

Definition at line 99 of file star_bhns.h.

◆ confo_auto

Scalar Lorene::Star_bhns::confo_auto
protected

Conformal factor generated by the star.

Definition at line 157 of file star_bhns.h.

◆ confo_comp

Scalar Lorene::Star_bhns::confo_comp
protected

Conformal factor generated by the companion black hole.

Definition at line 160 of file star_bhns.h.

◆ confo_tot

Scalar Lorene::Star_bhns::confo_tot
protected

Total conformal factor.

Definition at line 163 of file star_bhns.h.

◆ d_confo_auto

Vector Lorene::Star_bhns::d_confo_auto
protected

Derivative of the conformal factor generated by the star $ \partial_j \psi $.

Definition at line 168 of file star_bhns.h.

◆ d_confo_comp

Vector Lorene::Star_bhns::d_confo_comp
protected

Derivative of the conformal factor generated by the companion black hole.

Definition at line 173 of file star_bhns.h.

◆ d_lapconf_auto

Vector Lorene::Star_bhns::d_lapconf_auto
protected

Derivative of the lapconf function generated by the star $ \partial_j \alpha $.

Definition at line 130 of file star_bhns.h.

◆ d_lapconf_comp

Vector Lorene::Star_bhns::d_lapconf_comp
protected

Derivative of the lapconf function generated by the companion black hole.

Definition at line 135 of file star_bhns.h.

◆ d_psi

Vector Lorene::Star_bhns::d_psi
protected

Gradient of $\Psi$ (in the irrotational case) (Spherical components with respect to the mapping of the star)

Definition at line 82 of file star_bhns.h.

◆ d_shift_auto

Tensor Lorene::Star_bhns::d_shift_auto
protected

Derivative of the shift vector generated by the star $ \eta^{ik} \partial_k \beta^j $.

Definition at line 149 of file star_bhns.h.

◆ d_shift_comp

Tensor Lorene::Star_bhns::d_shift_comp
protected

Derivative of the shift vector generated by the companion black hole.

Definition at line 154 of file star_bhns.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.

◆ flat

Metric_flat Lorene::Star_bhns::flat
protected

Flat metric defined on the mapping (Spherical components with respect to the mapping of the star ).

Definition at line 192 of file star_bhns.h.

◆ gam

Scalar Lorene::Star_bhns::gam
protected

Lorentz factor between the fluid and the co-orbiting observer.

Definition at line 102 of file star_bhns.h.

◆ gam0

Scalar Lorene::Star_bhns::gam0
protected

Lorentz factor between the co-orbiting observer and the Eulerian one.

Definition at line 107 of file star_bhns.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.

◆ irrotational

bool Lorene::Star_bhns::irrotational
protected

true for an irrotational star, false for a corotating one

Definition at line 72 of file star_bhns.h.

◆ lapconf_auto

Scalar Lorene::Star_bhns::lapconf_auto
protected

Lapconf function generated by the star.

Definition at line 113 of file star_bhns.h.

◆ lapconf_comp

Scalar Lorene::Star_bhns::lapconf_comp
protected

Lapconf function generated by the companion black hole.

Definition at line 116 of file star_bhns.h.

◆ lapconf_tot

Scalar Lorene::Star_bhns::lapconf_tot
protected

Total lapconf function.

Definition at line 119 of file star_bhns.h.

◆ lapse_auto

Scalar Lorene::Star_bhns::lapse_auto
protected

Lapse function generated by the "star".

Definition at line 122 of file star_bhns.h.

◆ lapse_tot

Scalar Lorene::Star_bhns::lapse_tot
protected

Total lapse function.

Definition at line 125 of file star_bhns.h.

◆ loggam

Scalar Lorene::Star_bhns::loggam
protected

Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.

Definition at line 93 of file star_bhns.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.

◆ mp_aff

Map_af Lorene::Star_bhns::mp_aff
protected

Affine mapping for solving poisson's equations of metric quantities.

Definition at line 67 of file star_bhns.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.

◆ nzet

int Lorene::Star::nzet
protectedinherited

Number of domains of *mp occupied by the star.

Definition at line 183 of file star.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_mass_b

double* Lorene::Star::p_mass_b
mutableprotectedinherited

Baryon mass.

Definition at line 268 of file star.h.

◆ p_mass_b_bhns

double* Lorene::Star_bhns::p_mass_b_bhns
mutableprotected

Baryon mass.

Definition at line 225 of file star_bhns.h.

◆ p_mass_g

double* Lorene::Star::p_mass_g
mutableprotectedinherited

Gravitational mass.

Definition at line 269 of file star.h.

◆ p_mass_g_bhns

double* Lorene::Star_bhns::p_mass_g_bhns
mutableprotected

Gravitational mass.

Definition at line 226 of file star_bhns.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_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.

◆ pot_centri

Scalar Lorene::Star_bhns::pot_centri
protected

Centrifugal potential.

Definition at line 110 of file star_bhns.h.

◆ press

Scalar Lorene::Star::press
protectedinherited

Fluid pressure.

Definition at line 194 of file star.h.

◆ psi0

Scalar Lorene::Star_bhns::psi0
protected

Scalar potential $\Psi_0$ of the non-translational part of fluid 4-velocity (in the irrotational case)

Definition at line 77 of file star_bhns.h.

◆ psi4

Scalar Lorene::Star_bhns::psi4
protected

Fourth power of the total conformal factor.

Definition at line 176 of file star_bhns.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.

◆ shift_auto

Vector Lorene::Star_bhns::shift_auto
protected

Shift vector generated by the star.

Definition at line 138 of file star_bhns.h.

◆ shift_comp

Vector Lorene::Star_bhns::shift_comp
protected

Shift vector generated by the companion black hole.

Definition at line 141 of file star_bhns.h.

◆ shift_tot

Vector Lorene::Star_bhns::shift_tot
protected

Total shift vector.

Definition at line 144 of file star_bhns.h.

◆ ssjm1_confo

Scalar Lorene::Star_bhns::ssjm1_confo
protected

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

Definition at line 202 of file star_bhns.h.

◆ ssjm1_khi

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

◆ ssjm1_lapconf

Scalar Lorene::Star_bhns::ssjm1_lapconf
protected

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

Definition at line 197 of file star_bhns.h.

◆ ssjm1_wshift

Vector Lorene::Star_bhns::ssjm1_wshift
protected

Effective source at the previous step for the resolution of the vector Poisson equation for $W^i$ by means of Map_et::poisson .

$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 220 of file star_bhns.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.

◆ taij_auto

Sym_tensor Lorene::Star_bhns::taij_auto
protected

Part of the extrinsic curvature tensor $ A^{ij}$ generated by shift_auto , lapse_auto , and confo_auto .

Definition at line 182 of file star_bhns.h.

◆ taij_quad_auto

Scalar Lorene::Star_bhns::taij_quad_auto
protected

Part of the scalar $\eta_{ik} \eta_{jl} A^{ij} A^{kl}$ generated by $A_{ij}^{\rm auto}$.

Definition at line 187 of file star_bhns.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.

◆ wit_w

Vector Lorene::Star_bhns::wit_w
protected

Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.

(Spherical components with respect to the mapping of the star)

Definition at line 88 of file star_bhns.h.


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