26 #ifndef __STAR_BHNS_H_ 27 #define __STAR_BHNS_H_ 420 virtual void sauve(FILE *)
const ;
424 virtual ostream&
operator>>(ostream& )
const ;
430 virtual double mass_b()
const ;
432 virtual double mass_b_bhns(
bool kerrschild,
const double& mass_bh,
433 const double& sepa)
const ;
436 virtual double mass_g()
const ;
438 virtual double mass_g_bhns()
const ;
462 const double& sepa) ;
508 void kinema_bhns(
bool kerrschild,
const double& mass_bh,
509 const double& sepa,
double omega,
510 double x_rot,
double y_rot) ;
549 const double& sepa,
bool kerrschild,
550 int mer,
int mermax_ns,
int mermax_potvit,
551 int mermax_poisson,
int filter_r,
int filter_r_s,
552 int filter_p_s,
double relax_poisson,
553 double relax_potvit,
double thres_adapt,
555 const Tbl& fact_resize,
Tbl& diff) ;
574 double velo_pot_bhns(
const double& mass_bh,
const double& sepa,
576 int mermax,
double precis,
double relax) ;
584 double chi_rp(
double radius,
double phi) ;
617 double relax_met,
int mer,
int fmer_met) ;
Scalar psi4
Fourth power of the total conformal factor.
double phi_local_min(double phi_ini)
Azimuthal angle when the indicator of the mass-shedding takes its local minimum.
void fait_d_psi_bhns()
Computes the gradient of the total velocity potential .
Vector d_lapconf_auto
Derivative of the lapconf function generated by the star .
Vector shift_tot
Total shift vector.
Scalar lapconf_tot
Total lapconf function.
Class for stars in black hole-neutron star binaries.
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 fram...
double * p_mass_g_bhns
Gravitational mass.
const Scalar & get_psi0() const
Returns the non-translational part of the velocity potential.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
const Scalar & get_psi4() const
Returns the fourth power of the total conformal factor.
const Vector & get_d_lapconf_comp() const
Returns the derivative of the lapse function generated by the companion black hole.
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
virtual void del_deriv() const
Deletes all the derived quantities.
const Scalar & get_lapconf_comp() const
Returns the part of the lapconf function generated by the companion black hole.
Vector bsn
3-vector shift, divided by N , of the rotating coordinates, .
Equation of state base class.
Flat metric for tensor calculation.
Vector d_lapconf_comp
Derivative of the lapconf function generated by the companion black hole.
Tensor field of valence 0 (or component of a tensorial field).
void update_met_der_comp_bhns(const Hole_bhns &hole)
Computes derivative of metric quantities from the companion black hole.
Vector wit_w
Spatial projection of the fluid 3-velocity with respect to the co-orbiting observer.
Base class for coordinate mappings.
const Vector & get_bsn() const
Returns the shift vector, divided by N , of the rotating coordinates, .
Scalar loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
const Scalar & get_lapconf_auto() const
Returns the part of the lapconf function generated by the star.
double * p_mass_b_bhns
Baryon mass.
const Vector & get_shift_comp() const
Returns the part of the shift vector generated by the companion black hole.
const Sym_tensor & get_taij_auto() const
Returns the part of the extrinsic curvature tensor generated by the neutron star part...
Scalar & set_pot_centri()
Read/write the centrifugal potential.
const Scalar & get_lapconf_tot() const
Returns the total lapconf function.
Vector shift_auto
Shift vector generated by the star.
Tensor d_shift_auto
Derivative of the shift vector generated by the star .
void extr_curv_bhns()
Computes taij_auto , taij_quad_auto from shift_auto , lapse_auto , confo_auto .
const Vector & get_d_confo_auto() const
Returns the derivative of the conformal factor generated by the star.
Tensor field of valence 1.
Scalar lapse_tot
Total lapse function.
const Tensor & get_d_shift_comp() const
Returns the derivative of the shift vector generated by the companion black hole. ...
Scalar gam
Lorentz factor between the fluid and the co-orbiting observer.
Map_af mp_aff
Affine mapping for solving poisson's equations of metric quantities.
Vector d_psi
Gradient of (in the irrotational case) (Spherical components with respect to the mapping of the star...
const Vector & get_wit_w() const
Returns the spatial projection of the fluid 3-velocity with respect to the co-orbiting observer...
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star )...
void operator=(const Star_bhns &)
Assignment to another Star_bhns.
Scalar confo_tot
Total conformal factor.
const Vector & get_shift_auto() const
Returns the part of the shift vector generated by the star.
const Vector & get_d_psi() const
Returns the covariant derivative of the velocity potential (Spherical components with respect to the ...
virtual ~Star_bhns()
Destructor.
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 .
const Scalar & get_confo_tot() const
Returns the total conformal factor.
Vector d_confo_comp
Derivative of the conformal factor generated by the companion black hole.
double radius_p(double phi)
Radius of the star to the direction of and .
Vector & set_shift_comp()
Read/write of the shift vector generated by the companion black hole.
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 .
Star_bhns(Map &mp_i, int nzet_i, const Eos &eos_i, bool irrot_i)
Standard constructor.
void equil_spher_bhns(double ent_c, double precis)
Computes a spherical configuration.
Scalar & set_confo_auto()
Read/write of the conformal factor generated by the neutron star.
Scalar lapconf_comp
Lapconf function generated by the companion black hole.
Vector & set_shift_auto()
Read/write of the shift vector generated by the neutron star.
Scalar taij_quad_auto
Part of the scalar generated by .
const Scalar & get_confo_comp() const
Returns the part of the conformal factor generated by the companion black hole.
const Scalar & get_confo_auto() const
Returns the part of the conformal factor generated by the star.
Vector shift_comp
Shift vector generated by the companion black hole.
Scalar gam0
Lorentz factor between the co-orbiting observer and the Eulerian one.
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for by means...
Class for black holes in black hole-neutron star binaries.
Scalar pot_centri
Centrifugal potential.
const Vector & get_d_lapconf_auto() const
Returns the derivative of the lapse function generated by the star.
double phi_min()
Azimuthal angle when the indicator of the mass-shedding takes its minimum chi_min.
const Scalar & get_gam0() const
Returns the Lorentz factor gam0.
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...
double chi_rp(double radius, double phi)
Sensitive indicator of the mass-shedding to the direction of , , .
Scalar lapconf_auto
Lapconf function generated by the star.
bool irrotational
true for an irrotational star, false for a corotating one
Scalar & set_lapconf_comp()
Read/write of the lapconf function generated by the companion black hole.
Vector d_confo_auto
Derivative of the conformal factor generated by the star .
Scalar confo_auto
Conformal factor generated by the star.
const Scalar & get_lapse_tot() const
Returns the total lapse function.
Scalar confo_comp
Conformal factor generated by the companion black hole.
Sym_tensor taij_auto
Part of the extrinsic curvature tensor generated by shift_auto , lapse_auto , and confo_auto ...
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.
const Vector & get_d_confo_comp() const
Returns the derivative of the conformal factor generated by the companion black hole.
Scalar ssjm1_lapconf
Effective source at the previous step for the resolution of the Poisson equation for lapconf_auto ...
Class for computing a black hole - neutron star binary system with comparable mass () ...
Scalar & set_confo_comp()
Read/write of the conformal factor generated by the companion black hole.
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Scalar & set_lapconf_auto()
Read/write of the lapconf function generated by the neutron star.
bool is_irrotational() const
Returns true for an irrotational motion, false for a corotating one.
const Scalar & get_loggam() const
Returns the logarithm of the Lorentz factor between the fluid and the co-orbiting observer...
Scalar ssjm1_confo
Effective source at the previous step for the resolution of the Poisson equation for confo_auto ...
virtual void sauve(FILE *) const
Save in a file.
const Vector & get_shift_tot() const
Returns the part total shift vector.
Scalar lapse_auto
Lapse function generated by the "star".
virtual double mass_g() const
Gravitational mass.
Scalar psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case) ...
const Scalar & get_taij_quad_auto() const
Returns the part of the scalar generated by .
const Tensor & get_d_shift_auto() const
Returns the derivative of the shift vector generated by the star.
Class intended to describe valence-2 symmetric tensors.
Tensor d_shift_comp
Derivative of the shift vector generated by the companion black hole.
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 equat...
const Scalar & get_gam() const
Returns the Lorentz factor gam.
const Scalar & get_pot_centri() const
Returns the centrifugal potential.
virtual double mass_b() const
Baryon mass.