LORENE
bin_bhns_extr.C
1 /*
2  * Methods of class Bin_bhns_extr
3  *
4  * (see file bin_bhns_extr.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2004-2005 Keisuke Taniguchi
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License version 2
15  * as published by the Free Software Foundation.
16  *
17  * LORENE is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with LORENE; if not, write to the Free Software
24  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25  *
26  */
27 
28 
29 
30 /*
31  * $Id: bin_bhns_extr.C,v 1.11 2016/12/05 16:17:46 j_novak Exp $
32  * $Log: bin_bhns_extr.C,v $
33  * Revision 1.11 2016/12/05 16:17:46 j_novak
34  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35  *
36  * Revision 1.10 2014/10/13 08:52:41 j_novak
37  * Lorene classes and functions now belong to the namespace Lorene.
38  *
39  * Revision 1.9 2014/10/06 15:13:00 j_novak
40  * Modified #include directives to use c++ syntax.
41  *
42  * Revision 1.8 2005/02/28 23:06:16 k_taniguchi
43  * Modification of some arguments to include the case of the conformally
44  * flat background metric.
45  *
46  * Revision 1.7 2004/12/13 21:08:59 k_taniguchi
47  * Addition of some outputs in display_poly.
48  *
49  * Revision 1.6 2004/12/03 20:17:52 k_taniguchi
50  * Correction of "display_poly">
51  *
52  * Revision 1.5 2004/12/02 22:43:35 k_taniguchi
53  * Modification of "display_poly".
54  *
55  * Revision 1.4 2004/12/02 17:36:42 k_taniguchi
56  * Modification of "display_poly()".
57  *
58  * Revision 1.3 2004/12/01 16:37:11 k_taniguchi
59  * Modification of the output again.
60  *
61  * Revision 1.2 2004/12/01 16:27:09 k_taniguchi
62  * Modification of the output.
63  *
64  * Revision 1.1 2004/11/30 20:45:49 k_taniguchi
65  * *** empty log message ***
66  *
67  *
68  * $Header: /cvsroot/Lorene/C++/Source/Bin_bhns_extr/bin_bhns_extr.C,v 1.11 2016/12/05 16:17:46 j_novak Exp $
69  *
70  */
71 
72 // C headers
73 #include <cmath>
74 
75 // Lorene headers
76 #include "bin_bhns_extr.h"
77 #include "eos.h"
78 #include "utilitaires.h"
79 #include "unites.h"
80 
81  //--------------------------------//
82  // Constructors //
83  //--------------------------------//
84 
85 // Standard constructor
86 // --------------------
87 namespace Lorene {
88 Bin_bhns_extr::Bin_bhns_extr(Map& mp, int nzet, const Eos& eos, bool irrot,
89  bool relat, bool kerrs, bool multi)
90  : ref_triad(0., "Absolute frame Cartesian basis"),
91  star(mp, nzet, relat, eos, irrot, ref_triad, kerrs, multi)
92 {
93 
94  omega = 0. ;
95  separ = 0. ;
96  mass_bh = 0. ;
97 
98  // Pointers of derived quantities initialized to zero :
99  set_der_0x0() ;
100 
101 }
102 
103 // Copy constructor
104 // ----------------
106  : ref_triad(0., "Absolute frame Cartesian basis"),
107  star(bibi.star),
108  omega(bibi.omega),
109  separ(bibi.separ),
110  mass_bh(bibi.mass_bh)
111 {
112 
113  // Pointers of derived quantities initialized to zero :
114  set_der_0x0() ;
115 
116 }
117 
118 // Constructor from a file
119 // -----------------------
120 Bin_bhns_extr::Bin_bhns_extr(Map& mp, const Eos& eos, FILE* fich)
121  : ref_triad(0., "Absolute frame Cartesian basis"),
122  star(mp, eos, ref_triad, fich)
123 {
124 
125  // omega and x_axe are read in the file :
126  fread_be(&omega, sizeof(double), 1, fich) ;
127  fread_be(&separ, sizeof(double), 1, fich) ;
128  fread_be(&mass_bh, sizeof(double), 1, fich) ;
129 
130  // Pointers of derived quantities initialized to zero :
131  set_der_0x0() ;
132 
133 }
134 
135  //------------------------------//
136  // Destructor //
137  //------------------------------//
138 
140 {
141  del_deriv() ;
142 }
143 
144  //----------------------------------------------------//
145  // Management of derived quantities //
146  //----------------------------------------------------//
147 
149 
150  if (p_xa_barycenter_extr != 0x0) delete p_xa_barycenter_extr ;
151  if (p_ya_barycenter_extr != 0x0) delete p_ya_barycenter_extr ;
152  if (p_mass_b_extr != 0x0) delete p_mass_b_extr ;
153 
154  set_der_0x0() ;
155 
156 }
157 
159 
160  p_xa_barycenter_extr = 0x0 ;
161  p_ya_barycenter_extr = 0x0 ;
162  p_mass_b_extr = 0x0 ;
163 
164 }
165 
166  //------------------------------//
167  // Assignment //
168  //------------------------------//
169 
170 // Assignment to another Bin_bhns_extr
171 // -----------------------------------
173 
174  assert( bibi.ref_triad == ref_triad ) ;
175 
176  star = bibi.star ;
177 
178  omega = bibi.omega ;
179  separ = bibi.separ ;
180  mass_bh = bibi.mass_bh ;
181 
182  // ref_triad remains unchanged
183 
184  del_deriv() ; // Deletes all derived quantities
185 
186 }
187 
188  //---------------------------//
189  // Outputs //
190  //---------------------------//
191 
192 // Save in a file
193 // --------------
194 void Bin_bhns_extr::sauve(FILE* fich) const {
195 
196  star.sauve(fich) ;
197 
198  fwrite_be(&omega, sizeof(double), 1, fich) ;
199  fwrite_be(&separ, sizeof(double), 1, fich) ;
200  fwrite_be(&mass_bh, sizeof(double), 1, fich) ;
201 
202 }
203 
204 // Printing
205 // --------
206 ostream& operator<<(ostream& ost, const Bin_bhns_extr& bibi) {
207 
208  bibi >> ost ;
209  return ost ;
210 
211 }
212 
213 ostream& Bin_bhns_extr::operator>>(ostream& ost) const {
214 
215  using namespace Unites ;
216 
217  ost << endl ;
218  ost << "Binary BH-NS system" << endl ;
219  ost << "===================" << endl ;
220  ost << endl ;
221  if (star.in_kerrschild()) {
222  ost << "Kerr-Schild background metric" << endl ;
223  ost << "-----------------------------" << endl ;
224  }
225  else {
226  ost << "Conformally flat background metric" << endl ;
227  ost << "----------------------------------" << endl ;
228  }
229  if (star.with_multipole()) {
230  ost << "Multipole falloff boundary condition" << endl ;
231  ost << "------------------------------------" << endl ;
232  }
233  else {
234  ost << "1/r falloff boundary condition" << endl ;
235  ost << "------------------------------" << endl ;
236  }
237  ost << endl
238  << "Orbital angular velocity : "
239  << omega * f_unit << " rad/s" << endl ;
240  ost << "Coordinate separation between BH and NS : "
241  << separ / km << " km" << endl ;
242 
243  ost << endl << "Neutron star : " ;
244  ost << endl << "============ " << endl ;
245  ost << endl ;
246  if (star.is_relativistic()) {
247  ost << "Relativistic star" << endl ;
248  ost << "-----------------" << endl ;
249  }
250  else {
251  ost << "WARNING : BH-NS binary should be relativistic !!!" << endl ;
252  ost << "-------------------------------------------------" << endl ;
253  }
254  ost << "Number of domains occupied by the star : " << star.get_nzet()
255  << endl ;
256  ost << "Equation of state : " << endl ;
257  ost << star.get_eos() << endl ;
258 
259  ost << endl
260  << "Enthalpy at the coordinate origin : "
261  << star.get_ent()()(0,0,0,0) << " c^2" << endl ;
262  ost << "Proper baryon density at the coordinate origin : "
263  << star.get_nbar()()(0,0,0,0) << " x 0.1 fm^-3" << endl ;
264  ost << "Proper energy density at the coordinate origin : "
265  << star.get_ener()()(0,0,0,0) << " rho_nuc c^2" << endl ;
266  ost << "Pressure at the coordinate origin : "
267  << star.get_press()()(0,0,0,0) << " rho_nuc c^2" << endl ;
268  ost << endl
269  << "Lapse N at the coordinate origin : "
270  << star.get_nnn()()(0,0,0,0) << endl ;
271  ost << "Conformal factor A^2 at the coordinate origin : "
272  << star.get_a_car()()(0,0,0,0) << endl ;
273 
274  ost << endl
275  << "Equatorial radius (to BH) a_to : "
276  << star.ray_eq_pi()/km << " km" << endl ;
277  ost << "Equatorial radius (opp. to BH) a_opp : "
278  << star.ray_eq()/km << " km" << endl ;
279 
280  ost << endl
281  << "Baryon mass in isolation : " << star.mass_b() / msol
282  << " Mo" << endl
283  << "Gravitational mass in isolation : " << star.mass_g() / msol
284  << " Mo" << endl
285  << "Baryon mass in a binary system : " << mass_b_extr() / msol
286  << " Mo" << endl ;
287 
288  ost << endl ;
289  ost << "Star in a binary system" << endl ;
290  ost << "-----------------------" << endl ;
291 
292  if (star.is_irrotational()) {
293  ost << "irrotational configuration" << endl ;
294  }
295  else {
296  ost << "corotating configuration" << endl ;
297  }
298 
299  ost << "Absolute abscidia of the stellar center: "
300  << star.get_mp().get_ori_x()/km << " km" << endl ;
301  ost << "Orientation with respect to the absolute frame : "
302  << star.get_mp().get_rot_phi() << " rad" << endl ;
303 
304  double r0 = 0.5 * ( star.ray_eq() + star.ray_eq_pi() ) ;
305  // = 1 in Baumgarte et al.
306  double d_ns = separ + star.ray_eq() - r0 ;
307  // Orbital separation of Baumgarte et al.
308  double d_tilde = d_ns / r0 ; // Normalized separation of Baumgarte et al.
309 
310  ost << endl << "Comparison with those by Baumgarte et al. :" << endl ;
311  ost << " Radius r0 : " << r0 / km << " km " << endl ;
312  ost << " Separation d : " << d_ns / km << " km " << endl ;
313  ost << " Normalized sep. (d/r0) : " << d_tilde << endl ;
314 
315  ost << endl << "Black hole : " ;
316  ost << endl << "========== " << endl ;
317  ost << "Gravitational mass of BH : "
318  << mass_bh / msol << " M_sol" << endl ;
319 
320  return ost ;
321 
322 }
323 
324 // Display in polytropic units
325 // ---------------------------
326 void Bin_bhns_extr::display_poly(ostream& ost) const {
327 
328  using namespace Unites ;
329 
330  const Eos* p_eos = &( star.get_eos() ) ;
331  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos ) ;
332 
333  if (p_eos_poly != 0x0) {
334 
335  double kappa = p_eos_poly->get_kap() ;
336  double gamma = p_eos_poly->get_gam() ;
337  double kap_ns2 = pow( kappa, 0.5 /(gamma-1) ) ;
338 
339  // Polytropic unit of length in terms of r_unit :
340  double r_poly = kap_ns2 / sqrt(ggrav) ;
341 
342  // Polytropic unit of time in terms of t_unit :
343  double t_poly = r_poly ;
344 
345  // Polytropic unit of mass in terms of m_unit :
346  double m_poly = r_poly / ggrav ;
347 
348  // Polytropic unit of angular momentum in terms of j_unit :
349  // double j_poly = r_poly * r_poly / ggrav ;
350 
351  double r0 = 0.5 * ( star.ray_eq() + star.ray_eq_pi() ) ;
352  // = 1 in Baumgarte et al.
353  double d_ns = separ + star.ray_eq() - r0 ;
354  // Orbital separation of Baumgarte et al.
355 
356  ost.precision(16) ;
357  ost << endl << "Quantities in polytropic units : " ;
358  ost << endl << "==============================" << endl ;
359  ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
360  ost << " d_e_max : " << separ / r_poly << endl ;
361  ost << " d_G_x : " << xa_barycenter_extr() / r_poly << endl
362  << " d_G_y : " << ya_barycenter_extr() / r_poly << endl ;
363  ost << " d_bh/M_bh : " << d_ns/r_poly / (mass_bh/m_poly) << endl ;
364  ost << " Omega : " << omega * t_poly << endl ;
365  ost << " Omega M_bh : " << omega * t_poly * mass_bh / m_poly
366  << endl ;
367  ost << " M_bar(NS) : " << mass_b_extr() / m_poly << endl ;
368  ost << " M_bar(NS_0) : " << star.mass_b() / m_poly << endl ;
369  ost << " R_0(NS) : "
370  << 0.5 * (star.ray_eq() + star.ray_eq_pi()) / r_poly << endl ;
371  ost << " M_grv(BH) : " << mass_bh / m_poly << endl ;
372 
373  }
374 
375 }
376 }
bool with_multipole() const
Returns true for the multipole falloff boundary condition, false for the one.
double * p_xa_barycenter_extr
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or...
Definition: bin_bhns_extr.h:86
const Map & get_mp() const
Returns the mapping.
Definition: etoile.h:662
const Tenseur & get_a_car() const
Returns the total conformal factor .
Definition: etoile.h:736
void operator=(const Bin_bhns_extr &)
Assignment to another Bin_bhns_extr.
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
double xa_barycenter_extr() const
Absolute coordinate X of the barycenter of the baryon density in the Kerr-Schild background metric or...
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
Equation of state base class.
Definition: eos.h:209
double ray_eq() const
Coordinate radius at , [r_unit].
Base class for coordinate mappings.
Definition: map.h:682
double get_ori_x() const
Returns the x coordinate of the origin.
Definition: map.h:780
int get_nzet() const
Returns the number of domains occupied by the star.
Definition: etoile.h:665
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<)
bool is_irrotational() const
Returns true for an irrotational motion, false for a corotating one.
Definition: etoile.h:1095
virtual double mass_g() const
Gravitational mass.
double get_rot_phi() const
Returns the angle between the x –axis and X –axis.
Definition: map.h:787
Et_bin_bhns_extr star
Neutron star.
Definition: bin_bhns_extr.h:66
virtual void sauve(FILE *) const
Save in a file.
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
Definition: eos_poly.C:271
const Tenseur & get_press() const
Returns the fluid pressure.
Definition: etoile.h:685
bool in_kerrschild() const
Returns true for the Kerr-Schild background metric, false for the conformally flat one...
void del_deriv() const
Deletes all the derived quantities.
void display_poly(ostream &) const
Display in polytropic units.
double ya_barycenter_extr() const
in the Kerr-Schild background metric
void sauve(FILE *) const
Save in a file.
double get_kap() const
Returns the pressure coefficient (cf.
Definition: eos_poly.C:275
const Base_vect_cart ref_triad
Cartesian triad of the absolute reference frame.
Definition: bin_bhns_extr.h:63
Polytropic equation of state (relativistic case).
Definition: eos.h:812
const Tenseur & get_nnn() const
Returns the total lapse function N.
Definition: etoile.h:730
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.
Definition: fwrite_be.C:73
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition: fread_be.C:72
double ray_eq_pi() const
Coordinate radius at , [r_unit].
const Tenseur & get_ener() const
Returns the proper total energy density.
Definition: etoile.h:682
Bin_bhns_extr(Map &mp, int nzet, const Eos &eos, bool irrot, bool relat, bool kerrs, bool multi)
Standard constructor.
Definition: bin_bhns_extr.C:88
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
~Bin_bhns_extr()
Destructor.
double separ
Absolute orbital separation between two centers of BH and NS.
Definition: bin_bhns_extr.h:74
double omega
Angular velocity with respect to an asymptotically inertial observer.
Definition: bin_bhns_extr.h:71
double * p_mass_b_extr
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat one...
Definition: bin_bhns_extr.h:96
virtual double mass_b() const
Baryon mass.
const Tenseur & get_nbar() const
Returns the proper baryon density.
Definition: etoile.h:679
bool is_relativistic() const
Returns true for a relativistic star, false for a Newtonian one.
Definition: etoile.h:670
double mass_bh
Gravitational mass of BH.
Definition: bin_bhns_extr.h:77
const Eos & get_eos() const
Returns the equation of state.
Definition: etoile.h:673
const Tenseur & get_ent() const
Returns the enthalpy field.
Definition: etoile.h:676
double * p_ya_barycenter_extr
Absolute coordinate Y of the barycenter of the baryon density in the Kerr-Schild background metric...
Definition: bin_bhns_extr.h:91
double mass_b_extr() const
Baryon mass of the neutron star in the Kerr-Schild background metric or in the conformally flat...
Class for computing a Black hole - Neutron star binary system with an extreme mass ratio...
Definition: bin_bhns_extr.h:57