LORENE
|
Class for stars in black hole-neutron star binaries. More...
#include <star_bhns.h>
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... | |
Scalar & | set_pot_centri () |
Read/write the centrifugal potential. More... | |
Scalar & | set_lapconf_auto () |
Read/write of the lapconf function generated by the neutron star. More... | |
Scalar & | set_lapconf_comp () |
Read/write of the lapconf function generated by the companion black hole. More... | |
Vector & | set_shift_auto () |
Read/write of the shift vector generated by the neutron star. More... | |
Vector & | set_shift_comp () |
Read/write of the shift vector generated by the companion black hole. More... | |
Scalar & | set_confo_auto () |
Read/write of the conformal factor generated by the neutron star. More... | |
Scalar & | set_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 Scalar & | get_psi0 () const |
Returns the non-translational part of the velocity potential. More... | |
const Vector & | get_d_psi () const |
Returns the covariant derivative of the velocity potential (Spherical components with respect to the mapping of the star) More... | |
const Vector & | get_wit_w () const |
Returns the spatial projection of the fluid 3-velocity with respect to the co-orbiting observer. More... | |
const Scalar & | get_loggam () const |
Returns the logarithm of the Lorentz factor between the fluid and the co-orbiting observer. More... | |
const Vector & | get_bsn () const |
Returns the shift vector, divided by N , of the rotating coordinates, . More... | |
const Scalar & | get_gam () const |
Returns the Lorentz factor gam. More... | |
const Scalar & | get_gam0 () const |
Returns the Lorentz factor gam0. More... | |
const Scalar & | get_pot_centri () const |
Returns the centrifugal potential. More... | |
const Scalar & | get_lapconf_auto () const |
Returns the part of the lapconf function generated by the star. More... | |
const Scalar & | get_lapconf_comp () const |
Returns the part of the lapconf function generated by the companion black hole. More... | |
const Scalar & | get_lapconf_tot () const |
Returns the total lapconf function. More... | |
const Scalar & | get_lapse_auto () const |
const Scalar & | get_lapse_tot () const |
Returns the total lapse function. More... | |
const Vector & | get_d_lapconf_auto () const |
Returns the derivative of the lapse function generated by the star. More... | |
const Vector & | get_d_lapconf_comp () const |
Returns the derivative of the lapse function generated by the companion black hole. More... | |
const Vector & | get_shift_auto () const |
Returns the part of the shift vector generated by the star. More... | |
const Vector & | get_shift_comp () const |
Returns the part of the shift vector generated by the companion black hole. More... | |
const Vector & | get_shift_tot () const |
Returns the part total shift vector. More... | |
const Tensor & | get_d_shift_auto () const |
Returns the derivative of the shift vector generated by the star. More... | |
const Tensor & | get_d_shift_comp () const |
Returns the derivative of the shift vector generated by the companion black hole. More... | |
const Scalar & | get_confo_auto () const |
Returns the part of the conformal factor generated by the star. More... | |
const Scalar & | get_confo_comp () const |
Returns the part of the conformal factor generated by the companion black hole. More... | |
const Scalar & | get_confo_tot () const |
Returns the total conformal factor. More... | |
const Vector & | get_d_confo_auto () const |
Returns the derivative of the conformal factor generated by the star. More... | |
const Vector & | get_d_confo_comp () const |
Returns the derivative of the conformal factor generated by the companion black hole. More... | |
const Scalar & | get_psi4 () const |
Returns the fourth power of the total conformal factor. More... | |
const Sym_tensor & | get_taij_auto () const |
Returns the part of the extrinsic curvature tensor generated by the neutron star part. More... | |
const Scalar & | get_taij_quad_auto () const |
Returns the part of the scalar generated by . 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 . 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 by solving the continuity equation. More... | |
double | chi_rp (double radius, double phi) |
Sensitive indicator of the mass-shedding to the direction of , , . More... | |
double | radius_p (double phi) |
Radius of the star to the direction of and . 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... | |
Map & | set_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 Map & | get_mp () const |
Returns the mapping. More... | |
int | get_nzet () const |
Returns the number of domains occupied by the star. More... | |
const Eos & | get_eos () const |
Returns the equation of state. More... | |
const Scalar & | get_ent () const |
Returns the enthalpy field. More... | |
const Scalar & | get_nbar () const |
Returns the proper baryon density. More... | |
const Scalar & | get_ener () const |
Returns the proper total energy density. More... | |
const Scalar & | get_press () const |
Returns the fluid pressure. More... | |
const Scalar & | get_ener_euler () const |
Returns the total energy density with respect to the Eulerian observer. More... | |
const Scalar & | get_s_euler () const |
Returns the trace of the stress tensor in the Eulerian frame. More... | |
const Scalar & | get_gam_euler () const |
Returns the Lorentz factor between the fluid and Eulerian observers. More... | |
const Vector & | get_u_euler () const |
Returns the fluid 3-velocity with respect to the Eulerian observer. More... | |
const Tensor & | get_stress_euler () const |
Returns the spatial part of the stress-energy tensor with respect to the Eulerian observer. More... | |
const Scalar & | get_logn () const |
Returns the logarithm of the lapse N. More... | |
const Scalar & | get_nn () const |
Returns the lapse function N. More... | |
const Vector & | get_beta () const |
Returns the shift vector . More... | |
const Scalar & | get_lnq () const |
const Metric & | get_gamma () const |
Returns the 3-metric . More... | |
double | ray_eq () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_eq_pis2 () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_eq_pi () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_eq_3pis2 () const |
Coordinate radius at , [r_unit]. More... | |
double | ray_pole () const |
Coordinate radius at [r_unit]. More... | |
virtual const Itbl & | l_surf () const |
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on the surface at the collocation points in . More... | |
const Tbl & | xi_surf () const |
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinate on the surface at the collocation points in . More... | |
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 of the non-translational part of fluid 4-velocity (in the irrotational case) More... | |
Vector | d_psi |
Gradient of (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, . 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 . 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 . 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 . 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 generated by shift_auto , lapse_auto , and confo_auto . More... | |
Scalar | taij_quad_auto |
Part of the scalar generated by . 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 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 by means of Map_et::poisson . More... | |
double * | p_mass_b_bhns |
Baryon mass. More... | |
double * | p_mass_g_bhns |
Gravitational mass. More... | |
Map & | mp |
Mapping associated with the star. More... | |
int | nzet |
Number of domains of *mp occupied by the star. More... | |
const Eos & | eos |
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 , . More... | |
double * | p_ray_eq_pis2 |
Coordinate radius at , . More... | |
double * | p_ray_eq_pi |
Coordinate radius at , . More... | |
double * | p_ray_eq_3pis2 |
Coordinate radius at , . More... | |
double * | p_ray_pole |
Coordinate radius at . More... | |
Itbl * | p_l_surf |
Description of the stellar surface: 2-D Itbl containing the values of the domain index l on the surface at the collocation points in . More... | |
Tbl * | p_xi_surf |
Description of the stellar surface: 2-D Tbl containing the values of the radial coordinate on the surface at the collocation points in . More... | |
double * | p_mass_b |
Baryon mass. More... | |
double * | p_mass_g |
Gravitational mass. More... | |
Friends | |
class | Bin_bhns |
Standard constructor.
mp_i | Mapping on which the star will be defined |
nzet_i | Number of domains occupied by the star |
eos_i | Equation of state of the stellar matter |
irrot_i | should 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.
Lorene::Star_bhns::Star_bhns | ( | const Star_bhns & | star | ) |
Constructor from a file (see sauve(FILE*)
)
mp_i | Mapping on which the star will be defined |
eos_i | Equation of state of the stellar matter |
fich | input file (must have been created by the function sauve ) |
Definition at line 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.
|
virtual |
double Lorene::Star_bhns::chi_rp | ( | double | radius, |
double | phi | ||
) |
Sensitive indicator of the mass-shedding to the direction of , , .
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().
|
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().
|
protectedvirtualinherited |
Sets to ETATNONDEF
(undefined state) the hydrodynamical quantities relative to the Eulerian observer.
Reimplemented in Lorene::Star_bin_xcts, Lorene::Star_bin, Lorene::Star_rot, Lorene::Star_rot_Dirac, and Lorene::Star_rot_CFC.
Definition at line 333 of file star.C.
References Lorene::Star::del_deriv(), Lorene::Star::ener_euler, Lorene::Star::gam_euler, Lorene::Star::s_euler, Lorene::Scalar::set_etat_nondef(), Lorene::Tensor::set_etat_nondef(), Lorene::Star::stress_euler, and Lorene::Star::u_euler.
|
inherited |
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition at line 465 of file star.C.
References Lorene::Param::add_int(), Lorene::Scalar::allocate_all(), Lorene::Star::del_deriv(), Lorene::Star::ener, Lorene::Eos::ener_ent(), Lorene::Star::ent, Lorene::Star::eos, Lorene::Mg3d::get_grille3d(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Star::mp, Lorene::Star::nbar, Lorene::Eos::nbar_ent(), Lorene::Star::nzet, Lorene::Star::press, Lorene::Eos::press_ent(), Lorene::Mtbl::set(), Lorene::Scalar::set_domain(), Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::std_spectral_base(), Lorene::Mtbl::t, and Lorene::Grille3d::x.
void Lorene::Star_bhns::equil_spher_bhns | ( | double | ent_c, |
double | precis | ||
) |
Computes a spherical configuration.
ent_c | [input] Central enthalpy |
precis | [input] precision |
Definition at line 67 of file star_bhns_spher.C.
References Lorene::Scalar::annule(), confo_auto, confo_tot, Lorene::diffrel(), Lorene::Scalar::dsdr(), Lorene::Star::ener, Lorene::Star::ener_euler, Lorene::Star::ent, Lorene::Star::equation_of_state(), Lorene::exp(), Lorene::Star::gam_euler, Lorene::Scalar::get_dzpuis(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map_af::homothetie(), Lorene::Scalar::inc_dzpuis(), Lorene::Scalar::integrale(), lapconf_auto, lapconf_tot, lapse_auto, lapse_tot, Lorene::log(), Lorene::Star::mp, Lorene::norme(), Lorene::Star::nzet, Lorene::Map_af::poisson(), Lorene::pow(), Lorene::Star::press, psi4, Lorene::Star::s_euler, Lorene::Vector::set(), Lorene::Cmp::set_etat_qcq(), shift_auto, Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), Lorene::Star::u_euler, Lorene::Scalar::val_grid_point(), and Lorene::Map::val_r().
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.
ent_c | [input] Central enthalpy |
mass_bh | [input] BH mass in the background metric |
sepa | [input] Separation between NS and BH |
kerrschild | should 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.
|
virtualinherited |
Computes a spherical static configuration.
ent_c | [input] central value of the enthalpy |
precis | [input] threshold in the relative difference between the enthalpy fields of two consecutive steps to stop the iterative procedure (default value: 1.e-14) |
ent_limit | [input] : array of enthalpy values to be set at the boundaries between the domains; if set to 0x0 (default), the initial values will be kept. |
Definition at line 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().
void Lorene::Star_bhns::extr_curv_bhns | ( | ) |
Computes taij_auto
, taij_quad_auto
from shift_auto
, lapse_auto
, confo_auto
.
Definition at line 63 of file star_bhns_extr_curv.C.
References Lorene::Metric_flat::con(), confo_auto, Lorene::Metric_flat::cov(), d_shift_auto, flat, Lorene::Map::get_bvect_cart(), lapconf_auto, Lorene::Star::mp, Lorene::pow(), Lorene::Tensor::set_etat_qcq(), Lorene::Tensor::std_spectral_base(), Lorene::Scalar::std_spectral_base(), taij_auto, and taij_quad_auto.
void Lorene::Star_bhns::fait_d_psi_bhns | ( | ) |
Computes the gradient of the total velocity potential .
Definition at line 538 of file star_bhns.C.
References bsn, d_psi, Lorene::Scalar::deriv(), Lorene::Star::ent, Lorene::exp(), Lorene::Star::gam_euler, Lorene::Map::get_bvect_cart(), irrotational, Lorene::Star::mp, psi0, psi4, Lorene::Vector::set(), Lorene::Tensor::set_etat_nondef(), Lorene::Tensor::set_etat_qcq(), and Lorene::Scalar::std_spectral_base().
|
inlineinherited |
|
inline |
Returns the shift vector, divided by N , of the rotating coordinates, .
(Spherical components with respect to the mapping of the star)
Definition at line 327 of file star_bhns.h.
References bsn.
|
inline |
Returns the part of the conformal factor generated by the star.
Definition at line 383 of file star_bhns.h.
References confo_auto.
|
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.
|
inline |
Returns the total conformal factor.
Definition at line 391 of file star_bhns.h.
References confo_tot.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
|
inlineinherited |
Returns the proper total energy density.
Definition at line 370 of file star.h.
References Lorene::Star::ener.
|
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.
|
inlineinherited |
|
inlineinherited |
|
inline |
|
inline |
|
inlineinherited |
Returns the Lorentz factor between the fluid and Eulerian observers.
Definition at line 382 of file star.h.
References Lorene::Star::gam_euler.
|
inlineinherited |
|
inline |
Returns the part of the lapconf function generated by the star.
Definition at line 339 of file star_bhns.h.
References lapconf_auto.
|
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.
|
inline |
Returns the total lapconf function.
Definition at line 347 of file star_bhns.h.
References lapconf_tot.
|
inline |
|
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.
|
inlineinherited |
Returns the logarithm of the lapse N.
In the Newtonian case, this is the Newtonian gravitational potential (in units of ).
Definition at line 396 of file star.h.
References Lorene::Star::logn.
|
inlineinherited |
|
inlineinherited |
Returns the proper baryon density.
Definition at line 367 of file star.h.
References Lorene::Star::nbar.
|
inlineinherited |
|
inlineinherited |
Returns the number of domains occupied by the star.
Definition at line 358 of file star.h.
References Lorene::Star::nzet.
|
inline |
Returns the centrifugal potential.
Definition at line 336 of file star_bhns.h.
References pot_centri.
|
inlineinherited |
|
inline |
Returns the non-translational part of the velocity potential.
Definition at line 305 of file star_bhns.h.
References psi0.
|
inline |
Returns the fourth power of the total conformal factor.
Definition at line 404 of file star_bhns.h.
References psi4.
|
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.
|
inline |
Returns the part of the shift vector generated by the star.
Definition at line 364 of file star_bhns.h.
References shift_auto.
|
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.
|
inline |
Returns the part total shift vector.
Definition at line 372 of file star_bhns.h.
References shift_tot.
|
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.
|
inline |
Returns the part of the extrinsic curvature tensor generated by the neutron star part.
Definition at line 409 of file star_bhns.h.
References taij_auto.
|
inline |
Returns the part of the scalar generated by .
Definition at line 415 of file star_bhns.h.
References taij_quad_auto.
|
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.
|
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.
|
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.
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
.
kerrschild | should be true for a Kerr-Schild background, false for a Conformally flat one |
mass_bh | BH mass in the background metric |
sepa | Separation 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.
|
inline |
Returns true
for an irrotational motion, false
for a corotating one.
Definition at line 302 of file star_bhns.h.
References irrotational.
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.
kerrschild | should be true for a Kerr-Schild background, false for a Conformally flat one |
mass_bh | BH mass in the background metric |
sepa | Separation between NS and BH |
omega | angular velocity with respect to an asymptotically inertial observer |
x_rot | absolute X coordinate of the rotation axis |
y_rot | absolute 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.
|
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 .
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.
|
virtual |
Baryon mass.
Implements Lorene::Star.
Definition at line 66 of file star_bhns_global.C.
References confo_tot, Lorene::Star::gam_euler, Lorene::Scalar::integrale(), Lorene::Star::nbar, Lorene::Star::p_mass_b, Lorene::pow(), and Lorene::Scalar::std_spectral_base().
|
virtual |
void Lorene::Star_bhns::operator= | ( | const Star_bhns & | star | ) |
Assignment to another Star_bhns.
Definition at line 371 of file star_bhns.C.
References bsn, 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, del_deriv(), flat, gam, gam0, irrotational, lapconf_auto, lapconf_comp, lapconf_tot, lapse_auto, lapse_tot, loggam, mp_aff, Lorene::Star::operator=(), pot_centri, psi0, psi4, shift_auto, shift_comp, shift_tot, ssjm1_confo, ssjm1_khi, ssjm1_lapconf, ssjm1_wshift, taij_auto, taij_quad_auto, and wit_w.
|
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().
double Lorene::Star_bhns::phi_local_min | ( | double | phi_ini | ) |
Azimuthal angle when the indicator of the mass-shedding takes its local minimum.
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().
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().
double Lorene::Star_bhns::radius_p | ( | double | phi | ) |
Radius of the star to the direction of and .
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().
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 111 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq.
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 236 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq_3pis2.
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 189 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq_pi.
|
inherited |
Coordinate radius at , [r_unit].
Definition at line 141 of file star_global.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), Lorene::Star::mp, and Lorene::Star::p_ray_eq_pis2.
|
inherited |
Coordinate radius at [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.
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
.
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.
|
virtual |
Save in a file.
Reimplemented from Lorene::Star.
Definition at line 473 of file star_bhns.C.
References confo_auto, Lorene::Star::gam_euler, irrotational, lapconf_auto, psi0, Lorene::Star::sauve(), Lorene::Tensor::sauve(), Lorene::Scalar::sauve(), shift_auto, ssjm1_confo, ssjm1_khi, ssjm1_lapconf, and ssjm1_wshift.
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().
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().
|
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().
|
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().
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.
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.
|
inlineinherited |
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.
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.
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.
void Lorene::Star_bhns::update_met_der_comp_bhns | ( | const Hole_bhns & | hole | ) |
Computes derivative of metric quantities from the companion black hole.
hole | companion black hole |
Definition at line 63 of file star_bhns_upmetr_der.C.
References Lorene::Hole_bhns::get_d_lapconf_auto().
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
,
hole | companion black hole |
star_prev | previous value of the star |
relax | relaxation parameter |
Definition at line 63 of file star_bhns_upmetr.C.
References Lorene::Hole_bhns::get_lapconf_auto(), and Lorene::Black_hole::get_mp().
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 by solving the continuity equation.
mass_bh | [input] BH mass in the background metric |
sepa | [input] Separation between NS and BH |
kerrschild | should 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 between two successive steps is lower than precis . |
relax | [input] Relaxation factor. |
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().
|
inherited |
Description of the stellar surface: returns a 2-D Tbl
containing the values of the radial coordinate on the surface at the collocation points in .
The stellar surface is defined as the location where the enthalpy (member ent
) vanishes.
Definition at line 92 of file star_global.C.
References Lorene::Star::l_surf(), Lorene::Star::p_l_surf, and Lorene::Star::p_xi_surf.
|
protected |
3-vector shift, divided by N , of the rotating coordinates, .
(Spherical components with respect to the mapping of the star)
Definition at line 99 of file star_bhns.h.
|
protected |
Conformal factor generated by the star.
Definition at line 157 of file star_bhns.h.
|
protected |
Conformal factor generated by the companion black hole.
Definition at line 160 of file star_bhns.h.
|
protected |
Total conformal factor.
Definition at line 163 of file star_bhns.h.
|
protected |
Derivative of the conformal factor generated by the star .
Definition at line 168 of file star_bhns.h.
|
protected |
Derivative of the conformal factor generated by the companion black hole.
Definition at line 173 of file star_bhns.h.
|
protected |
Derivative of the lapconf function generated by the star .
Definition at line 130 of file star_bhns.h.
|
protected |
Derivative of the lapconf function generated by the companion black hole.
Definition at line 135 of file star_bhns.h.
|
protected |
Gradient of (in the irrotational case) (Spherical components with respect to the mapping of the star)
Definition at line 82 of file star_bhns.h.
|
protected |
Derivative of the shift vector generated by the star .
Definition at line 149 of file star_bhns.h.
|
protected |
Derivative of the shift vector generated by the companion black hole.
Definition at line 154 of file star_bhns.h.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
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.
|
protected |
Lorentz factor between the fluid and the co-orbiting observer.
Definition at line 102 of file star_bhns.h.
|
protected |
Lorentz factor between the co-orbiting observer and the Eulerian one.
Definition at line 107 of file star_bhns.h.
|
protectedinherited |
|
protected |
true
for an irrotational star, false
for a corotating one
Definition at line 72 of file star_bhns.h.
|
protected |
Lapconf function generated by the star.
Definition at line 113 of file star_bhns.h.
|
protected |
Lapconf function generated by the companion black hole.
Definition at line 116 of file star_bhns.h.
|
protected |
Total lapconf function.
Definition at line 119 of file star_bhns.h.
|
protected |
Lapse function generated by the "star".
Definition at line 122 of file star_bhns.h.
|
protected |
Total lapse function.
Definition at line 125 of file star_bhns.h.
|
protected |
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition at line 93 of file star_bhns.h.
|
protectedinherited |
|
protectedinherited |
|
protected |
Affine mapping for solving poisson's equations of metric quantities.
Definition at line 67 of file star_bhns.h.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotected |
Baryon mass.
Definition at line 225 of file star_bhns.h.
|
mutableprotectedinherited |
|
mutableprotected |
Gravitational mass.
Definition at line 226 of file star_bhns.h.
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
protected |
Centrifugal potential.
Definition at line 110 of file star_bhns.h.
|
protectedinherited |
|
protected |
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case)
Definition at line 77 of file star_bhns.h.
|
protected |
Fourth power of the total conformal factor.
Definition at line 176 of file star_bhns.h.
|
protectedinherited |
|
protected |
Shift vector generated by the star.
Definition at line 138 of file star_bhns.h.
|
protected |
Shift vector generated by the companion black hole.
Definition at line 141 of file star_bhns.h.
|
protected |
Total shift vector.
Definition at line 144 of file star_bhns.h.
|
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.
|
protected |
Effective source at the previous step for the resolution of the Poisson equation for the scalar by means of Map_et::poisson
.
is an intermediate quantity for the resolution of the elliptic equation for the shift vector
Definition at line 210 of file star_bhns.h.
|
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.
|
protected |
Effective source at the previous step for the resolution of the vector Poisson equation for by means of Map_et::poisson
.
is an intermediate quantity for the resolution of the elliptic equation for the shift vector (Components with respect to the Cartesian triad associated with the mapping mp
)
Definition at line 220 of file star_bhns.h.
|
protectedinherited |
|
protected |
Part of the extrinsic curvature tensor generated by shift_auto
, lapse_auto
, and confo_auto
.
Definition at line 182 of file star_bhns.h.
|
protected |
Part of the scalar generated by .
Definition at line 187 of file star_bhns.h.
|
protectedinherited |
|
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.