63#include "utilitaires.h"
74 bool irrot_ns,
bool kerrschild_i,
75 bool bc_lapconf_nd,
bool bc_lapconf_fs,
bool irrot_bh,
77 :
ref_triad(0.,
"Absolute frame Cartesian basis"),
78 hole(mp_bh,kerrschild_i,bc_lapconf_nd,bc_lapconf_fs,irrot_bh, massbh),
79 star(mp_ns, nzet_i, eos_i, irrot_ns)
95 :
ref_triad(0.,
"Absolute frame Cartesian basis"),
112 :
ref_triad(0.,
"Absolute frame Cartesian basis"),
114 star(mp_ns, eos_i, fich)
239 ost <<
"BH-NS binary system" << endl ;
240 ost <<
"===================" << endl ;
243 <<
"Orbital angular velocity : "
244 <<
omega * f_unit <<
" [rad/s]" << endl ;
245 ost <<
"Coordinate separation between BH and NS : "
246 <<
separ / km <<
" [km]" << endl ;
247 ost <<
"Position of the rotation axis : "
248 <<
x_rot / km <<
" [km]" << endl ;
250 ost << endl <<
"Black hole : " ;
251 ost << endl <<
"========== " << endl ;
253 int nt = (
hole.get_mp()).get_mg()->get_nt(1) ;
255 ost <<
"Irreducible mass of BH : "
256 <<
hole.mass_irr_bhns() / msol <<
" [Mo]" << endl ;
257 ost <<
"Mass in the background : "
258 <<
hole.get_mass_bh() / msol <<
" [Mo]" << endl ;
259 ost <<
"Radius of the apparent horison : "
260 <<
hole.rad_ah() / km <<
" [km]" << endl ;
261 ost <<
"Lapse function on the AH : "
262 << (
hole.get_lapse_tot()).val_grid_point(1,0,nt-1,0) << endl ;
263 ost <<
"Conformal factor on the AH : "
264 << (
hole.get_confo_tot()).val_grid_point(1,0,nt-1,0) << endl ;
265 ost <<
"shift(1) on the AH : "
266 << (
hole.get_shift_tot()(1)).val_grid_point(1,0,nt-1,0) << endl ;
267 ost <<
"shift(2) on the AH : "
268 << (
hole.get_shift_tot()(2)).val_grid_point(1,0,nt-1,0) << endl ;
269 ost <<
"shift(3) on the AH : "
270 << (
hole.get_shift_tot()(3)).val_grid_point(1,0,nt-1,0) << endl ;
272 ost << endl <<
"Neutron star : " ;
273 ost << endl <<
"============ " << endl ;
275 ost <<
"Baryon mass of NS in isolation : "
276 <<
star.mass_b()/msol <<
" [Mo]" << endl ;
277 ost <<
"Gravitational mass of NS : "
278 <<
star.mass_g_bhns()/msol <<
" [Mo]" << endl ;
279 ost <<
"Baryon mass of NS in BH bg. : "
282 ost <<
"Coordinate radius R_eq_tow : "
283 <<
star.ray_eq_pi() / km <<
" [km]" << endl ;
284 ost <<
"Coordinate radius R_eq_opp : "
285 <<
star.ray_eq() / km <<
" [km]" << endl ;
286 ost <<
"Coordinate radius R_eq_orb : "
287 <<
star.ray_eq_pis2() / km <<
" [km]" << endl ;
288 ost <<
"Coordinate radius R_eq_orb_opp : "
289 <<
star.ray_eq_3pis2() / km <<
" [km]" << endl ;
290 ost <<
"Coordinate radius R_pole : "
291 <<
star.ray_pole() / km <<
" [km]" << endl ;
292 ost <<
"Central enthalpy H_ent : "
293 << (
star.get_ent()).val_grid_point(0,0,0,0) << endl ;
294 ost <<
"Lapse function at the center of NS : "
295 << (
star.get_lapse_tot()).val_grid_point(0,0,0,0) << endl ;
296 ost <<
"Conformal factor at the center of NS : "
297 << (
star.get_confo_tot()).val_grid_point(0,0,0,0) << endl ;
298 ost <<
"shift(1) at the center of NS : "
299 << (
star.get_shift_tot()(1)).val_grid_point(0,0,0,0) << endl ;
300 ost <<
"shift(2) at the center of NS : "
301 << (
star.get_shift_tot()(2)).val_grid_point(0,0,0,0) << endl ;
302 ost <<
"shift(3) at the center of NS : "
303 << (
star.get_shift_tot()(3)).val_grid_point(0,0,0,0) << endl ;
316 const Eos* p_eos = &(
star.get_eos() ) ;
319 if (p_eos_poly != 0x0) {
321 double kappa = p_eos_poly->
get_kap() ;
322 double gamma = p_eos_poly->
get_gam() ;
323 double kap_ns2 =
pow( kappa, 0.5 /(gamma-1) ) ;
326 double r_poly = kap_ns2 /
sqrt(ggrav) ;
329 double t_poly = r_poly ;
332 double m_poly = r_poly / ggrav ;
337 double r0 = 0.5 * (
star.ray_eq() +
star.ray_eq_pi() ) ;
342 double y_ns =
star.get_mp().get_ori_y() ;
346 ost << endl <<
"Quantities in polytropic units : " ;
347 ost << endl <<
"==============================" << endl ;
348 ost <<
" ( r_poly = " << r_poly / km <<
" km )" << endl ;
349 ost <<
" d_separ : " <<
separ / r_poly << endl ;
350 ost <<
" d_coord : " << d_coord / r_poly << endl ;
351 ost <<
" Omega : " <<
omega * t_poly << endl ;
352 ost <<
" Omega M_irr : "
353 <<
omega * t_poly *
hole.mass_irr_bhns() / m_poly << endl ;
354 ost <<
" Omega_spin : " <<
hole.get_omega_spin() * t_poly << endl ;
355 ost <<
" M_irr : " <<
hole.mass_irr_bhns() / m_poly << endl ;
356 ost <<
" M_bh : " <<
hole.get_mass_bh() / m_poly << endl ;
357 ost <<
" R_ah : " <<
hole.rad_ah() / r_poly << endl ;
358 ost <<
" M_bar(NS)iso : " <<
star.mass_b() / m_poly
360 ost <<
" M_bar(NS) : "
363 ost <<
" M_g(NS)iso : " <<
star.mass_g_bhns() / m_poly << endl ;
364 ost <<
" R_0(NS) : " << r0 / r_poly << endl ;
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.
friend ostream & operator<<(ostream &, const Bin_bhns &)
Display.
virtual void sauve(FILE *) const
Save in a file.
Bin_bhns(Map &mp_bh, Map &mp_ns, int nzet, const Eos &eos, bool irrot_ns, bool kerrschild, bool bc_lapse_nd, bool bc_lapse_fs, bool irrot_bh, double mass_bh)
Relative error on the Hamiltonian constraint.
double y_rot
Absolute Y coordinate of the rotation axis.
const Base_vect_cart ref_triad
Cartesian triad of the absolute reference frame.
double * p_mass_kom_bhns_surf
Total Komar mass of the system calculated by the surface integral at infinity.
Hole_bhns hole
Black hole.
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 * 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.
double * p_omega_two_points
Orbital angular velocity derived from another method.
double separ
Absolute orbital separation between two centers of BH and NS.
Tbl * p_angu_mom_bhns
Total angular momentum of the system.
double * p_ya_barycenter
Absolute coordinate Y of the barycenter of the baryon density.
double x_rot
Absolute X coordinate of the rotation axis.
virtual ~Bin_bhns()
Destructor.
void operator=(const Bin_bhns &)
Assignment to another Bin_bhns.
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 <<).
Tbl * p_line_mom_bhns
Total linear momentum of the system.
Polytropic equation of state (relativistic case).
double get_gam() const
Returns the adiabatic index (cf. Eq. (3)).
double get_kap() const
Returns the pressure coefficient (cf.
Equation of state base class.
Cmp sqrt(const Cmp &)
Square root.
Cmp pow(const Cmp &, int)
Power .
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Standard units of space, time and mass.