170 bool irrot_ns,
bool kerrschild,
171 bool bc_lapse_nd,
bool bc_lapse_fs,
bool irrot_bh,
252 virtual void sauve(FILE *)
const ;
271 double mass_adm_bhns_vol()
const ;
276 double mass_kom_bhns_vol()
const ;
341 void orbit_omega(
double fact_omeg_min,
double fact_omeg_max) ;
359 void rotation_axis_y(
double thres_rot,
double rot_exp_y,
double fact) ;
369 void shift_analytic(
double reduce_shift_bh,
double reduce_shift_ns) ;
373ostream& operator<<(ostream& ,
const Bin_bhns& ) ;
Cartesian vectorial bases (triads).
Class for computing a black hole - neutron star binary system with comparable mass ()
void display_poly(ostream &) const
Display in polytropic units.
void del_deriv() const
Deletes all the derived quantities.
double * p_mass_adm_bhns_surf
Total ADM mass of the system calculated by the surface integral at infinity.
virtual void sauve(FILE *) const
Save in a file.
double y_rot
Absolute Y coordinate of the rotation axis.
const Base_vect_cart ref_triad
Cartesian triad of the absolute reference frame.
double get_separ() const
Returns the coordinate separation of the binary system [{\tt r_unit}].
double * p_mass_kom_bhns_surf
Total Komar mass of the system calculated by the surface integral at infinity.
Hole_bhns hole
Black hole.
double mass_kom_bhns_surf() const
Total Komar mass.
double get_omega() const
Returns the orbital angular velocity [{\tt f_unit}].
double ya_barycenter() const
Absolute coordinate Y of the barycenter of the baryon density.
void shift_analytic(double reduce_shift_bh, double reduce_shift_ns)
Sets some analytical template for the initial shift vector.
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
double omega
Angular velocity with respect to an asymptotically inertial observer.
double & set_y_rot()
Sets the absolute coordinate Y of the rotation axis [{\tt r_unit}].
double * p_mass_adm_bhns_vol
Total ADM mass of the system calculated by the volume integral and the surface integral at the appare...
double * p_virial_bhns_vol
Virial theorem error calculated by the ADM mass and the Komar mass of the volume integral.
const Tbl & line_mom_bhns() const
Total linear momentum.
friend ostream & operator<<(ostream &, const Bin_bhns &)
Display.
double * p_omega_two_points
Orbital angular velocity derived from another method.
double xa_barycenter() const
Absolute coordinate X of the barycenter of the baryon density.
double separ
Absolute orbital separation between two centers of BH and NS.
void rotation_axis_y(double thres_rot, double rot_exp_y, double fact)
Computes the position of the rotation axis Y.
Tbl * p_angu_mom_bhns
Total angular momentum of the system.
double get_y_rot() const
Returns the absolute coordinate Y of the rotation axis [{\tt r_unit}].
const Hole_bhns & get_bh() const
Returns a reference to the black hole.
double virial_bhns_surf() const
Estimates the relative error on the virial theorem $|1 - M_{\rm Komar} / M_{\rm ADM}|$.
void orbit_omega(double fact_omeg_min, double fact_omeg_max)
Estimates the relative error on the Hamiltonian constraint.
const Tbl & angu_mom_bhns() const
Total angular momentum.
double virial_bhns_vol() const
Estimates the relative error on the virial theorem $|1 - M_{\rm Komar} / M_{\rm ADM}|$.
double * p_ya_barycenter
Absolute coordinate Y of the barycenter of the baryon density.
double omega_two_points() const
Orbital angular velocity derived from another method.
double x_rot
Absolute X coordinate of the rotation axis.
virtual ~Bin_bhns()
Destructor.
Hole_bhns & set_bh()
Read/write of the black hole.
void operator=(const Bin_bhns &)
Assignment to another Bin_bhns.
const Star_bhns & get_ns() const
Returns a reference to the neutron star.
double & set_x_rot()
Sets the absolute coordinate X of the rotation axis [{\tt r_unit}].
double * p_virial_bhns_surf
Virial theorem error calculated by the ADM mass and the Komar mass of the surface integral at infinit...
double * p_mass_kom_bhns_vol
Total Komar mass of the system calculated by the volume integral and the surface integral at the appa...
Star_bhns star
Neutron star.
double * p_xa_barycenter
Absolute coordinate X of the barycenter of the baryon density.
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<)
void rotation_axis_x(double rot_exp_x)
Computes the position of the rotation axis X.
double & set_separ()
Sets the orbital separation [{\tt r_unit}].
Tbl * p_line_mom_bhns
Total linear momentum of the system.
double mass_adm_bhns_surf() const
Total ADM mass.
double & set_omega()
Sets the orbital angular velocity [{\tt f_unit}].
double get_x_rot() const
Returns the absolute coordinate X of the rotation axis [{\tt r_unit}].
Star_bhns & set_ns()
Read/write of the neutron star.
Equation of state base class.
Class for black holes in black hole-neutron star binaries.
Base class for coordinate mappings.
Class for stars in black hole-neutron star binaries.