28char bin_bhns_extr_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Bin_bhns_extr/bin_bhns_extr.C,v 1.10 2014/10/13 08:52:41 j_novak Exp $" ;
73#include "bin_bhns_extr.h"
75#include "utilitaires.h"
86 bool relat,
bool kerrs,
bool multi)
87 :
ref_triad(0.,
"Absolute frame Cartesian basis"),
103 :
ref_triad(0.,
"Absolute frame Cartesian basis"),
118 :
ref_triad(0.,
"Absolute frame Cartesian basis"),
215 ost <<
"Binary BH-NS system" << endl ;
216 ost <<
"===================" << endl ;
218 if (
star.in_kerrschild()) {
219 ost <<
"Kerr-Schild background metric" << endl ;
220 ost <<
"-----------------------------" << endl ;
223 ost <<
"Conformally flat background metric" << endl ;
224 ost <<
"----------------------------------" << endl ;
226 if (
star.with_multipole()) {
227 ost <<
"Multipole falloff boundary condition" << endl ;
228 ost <<
"------------------------------------" << endl ;
231 ost <<
"1/r falloff boundary condition" << endl ;
232 ost <<
"------------------------------" << endl ;
235 <<
"Orbital angular velocity : "
236 <<
omega * f_unit <<
" rad/s" << endl ;
237 ost <<
"Coordinate separation between BH and NS : "
238 <<
separ / km <<
" km" << endl ;
240 ost << endl <<
"Neutron star : " ;
241 ost << endl <<
"============ " << endl ;
243 if (
star.is_relativistic()) {
244 ost <<
"Relativistic star" << endl ;
245 ost <<
"-----------------" << endl ;
248 ost <<
"WARNING : BH-NS binary should be relativistic !!!" << endl ;
249 ost <<
"-------------------------------------------------" << endl ;
251 ost <<
"Number of domains occupied by the star : " <<
star.get_nzet()
253 ost <<
"Equation of state : " << endl ;
254 ost <<
star.get_eos() << endl ;
257 <<
"Enthalpy at the coordinate origin : "
258 <<
star.get_ent()()(0,0,0,0) <<
" c^2" << endl ;
259 ost <<
"Proper baryon density at the coordinate origin : "
260 <<
star.get_nbar()()(0,0,0,0) <<
" x 0.1 fm^-3" << endl ;
261 ost <<
"Proper energy density at the coordinate origin : "
262 <<
star.get_ener()()(0,0,0,0) <<
" rho_nuc c^2" << endl ;
263 ost <<
"Pressure at the coordinate origin : "
264 <<
star.get_press()()(0,0,0,0) <<
" rho_nuc c^2" << endl ;
266 <<
"Lapse N at the coordinate origin : "
267 <<
star.get_nnn()()(0,0,0,0) << endl ;
268 ost <<
"Conformal factor A^2 at the coordinate origin : "
269 <<
star.get_a_car()()(0,0,0,0) << endl ;
272 <<
"Equatorial radius (to BH) a_to : "
273 <<
star.ray_eq_pi()/km <<
" km" << endl ;
274 ost <<
"Equatorial radius (opp. to BH) a_opp : "
275 <<
star.ray_eq()/km <<
" km" << endl ;
278 <<
"Baryon mass in isolation : " <<
star.mass_b() / msol
280 <<
"Gravitational mass in isolation : " <<
star.mass_g() / msol
282 <<
"Baryon mass in a binary system : " <<
mass_b_extr() / msol
286 ost <<
"Star in a binary system" << endl ;
287 ost <<
"-----------------------" << endl ;
289 if (
star.is_irrotational()) {
290 ost <<
"irrotational configuration" << endl ;
293 ost <<
"corotating configuration" << endl ;
296 ost <<
"Absolute abscidia of the stellar center: "
297 <<
star.get_mp().get_ori_x()/km <<
" km" << endl ;
298 ost <<
"Orientation with respect to the absolute frame : "
299 <<
star.get_mp().get_rot_phi() <<
" rad" << endl ;
301 double r0 = 0.5 * (
star.ray_eq() +
star.ray_eq_pi() ) ;
303 double d_ns =
separ +
star.ray_eq() - r0 ;
305 double d_tilde = d_ns / r0 ;
307 ost << endl <<
"Comparison with those by Baumgarte et al. :" << endl ;
308 ost <<
" Radius r0 : " << r0 / km <<
" km " << endl ;
309 ost <<
" Separation d : " << d_ns / km <<
" km " << endl ;
310 ost <<
" Normalized sep. (d/r0) : " << d_tilde << endl ;
312 ost << endl <<
"Black hole : " ;
313 ost << endl <<
"========== " << endl ;
314 ost <<
"Gravitational mass of BH : "
315 <<
mass_bh / msol <<
" M_sol" << endl ;
327 const Eos* p_eos = &(
star.get_eos() ) ;
330 if (p_eos_poly != 0x0) {
332 double kappa = p_eos_poly->
get_kap() ;
333 double gamma = p_eos_poly->
get_gam() ;
334 double kap_ns2 =
pow( kappa, 0.5 /(gamma-1) ) ;
337 double r_poly = kap_ns2 /
sqrt(ggrav) ;
340 double t_poly = r_poly ;
343 double m_poly = r_poly / ggrav ;
348 double r0 = 0.5 * (
star.ray_eq() +
star.ray_eq_pi() ) ;
350 double d_ns =
separ +
star.ray_eq() - r0 ;
354 ost << endl <<
"Quantities in polytropic units : " ;
355 ost << endl <<
"==============================" << endl ;
356 ost <<
" ( r_poly = " << r_poly / km <<
" km )" << endl ;
357 ost <<
" d_e_max : " <<
separ / r_poly << endl ;
360 ost <<
" d_bh/M_bh : " << d_ns/r_poly / (
mass_bh/m_poly) << endl ;
361 ost <<
" Omega : " <<
omega * t_poly << endl ;
362 ost <<
" Omega M_bh : " <<
omega * t_poly *
mass_bh / m_poly
364 ost <<
" M_bar(NS) : " <<
mass_b_extr() / m_poly << endl ;
365 ost <<
" M_bar(NS_0) : " <<
star.mass_b() / m_poly << endl ;
367 << 0.5 * (
star.ray_eq() +
star.ray_eq_pi()) / r_poly << endl ;
368 ost <<
" M_grv(BH) : " <<
mass_bh / m_poly << endl ;
const Base_vect_cart ref_triad
Cartesian triad of the absolute reference frame.
double * p_ya_barycenter_extr
Absolute coordinate Y of the barycenter of the baryon density in the Kerr-Schild background metric.
double * p_mass_b_extr
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat one.
double ya_barycenter_extr() const
in the Kerr-Schild background metric
double mass_b_extr() const
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat.
Et_bin_bhns_extr star
Neutron star.
void operator=(const Bin_bhns_extr &)
Assignment to another Bin_bhns_extr.
Bin_bhns_extr(Map &mp, int nzet, const Eos &eos, bool irrot, bool relat, bool kerrs, bool multi)
Standard constructor.
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
double separ
Absolute orbital separation between two centers of BH and NS.
friend ostream & operator<<(ostream &, const Bin_bhns_extr &)
Display.
double mass_bh
Gravitational mass of BH.
double xa_barycenter_extr() const
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or...
void display_poly(ostream &) const
Display in polytropic units.
void sauve(FILE *) const
Save in a file.
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<).
~Bin_bhns_extr()
Destructor.
double omega
Angular velocity with respect to an asymptotically inertial observer.
double * p_xa_barycenter_extr
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or...
void del_deriv() const
Deletes all the derived quantities.
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.