LORENE
bin_bhns.C
1/*
2 * Methods of class Bin_bhns
3 *
4 * (see file bin_bhns.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2005-2007 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.C,v 1.5 2016/12/05 16:17:45 j_novak Exp $
32 * $Log: bin_bhns.C,v $
33 * Revision 1.5 2016/12/05 16:17:45 j_novak
34 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35 *
36 * Revision 1.4 2014/10/13 08:52:41 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.3 2014/10/06 15:13:00 j_novak
40 * Modified #include directives to use c++ syntax.
41 *
42 * Revision 1.2 2008/05/15 18:58:01 k_taniguchi
43 * Change of some parameters and introduction of new
44 * global quantities.
45 *
46 * Revision 1.1 2007/06/22 01:08:53 k_taniguchi
47 * *** empty log message ***
48 *
49 *
50 * $Header: /cvsroot/Lorene/C++/Source/Bin_bhns/bin_bhns.C,v 1.5 2016/12/05 16:17:45 j_novak Exp $
51 *
52 */
53
54// C++ headers
55//#include <>
56
57// C headers
58#include <cmath>
59
60// Lorene headers
61#include "bin_bhns.h"
62#include "eos.h"
63#include "utilitaires.h"
64#include "unites.h"
65
66 //---------------------//
67 // Constructor //
68 //---------------------//
69
70// Standard constructor
71// --------------------
72namespace Lorene {
73Bin_bhns::Bin_bhns(Map& mp_bh, Map& mp_ns, int nzet_i, const Eos& eos_i,
74 bool irrot_ns, bool kerrschild_i,
75 bool bc_lapconf_nd, bool bc_lapconf_fs, bool irrot_bh,
76 double massbh)
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)
80{
81
82 omega = 0. ;
83 separ = 0. ;
84 x_rot = 0. ;
85 y_rot = 0. ;
86
87 // Pointers of derived quantities are initialized to zero
88 set_der_0x0() ;
89
90}
91
92// Copy constructor
93// ----------------
95 : ref_triad(0., "Absolute frame Cartesian basis"),
96 hole(bibi.hole),
97 star(bibi.star),
98 omega(bibi.omega),
99 separ(bibi.separ),
100 x_rot(bibi.x_rot),
101 y_rot(bibi.y_rot)
102{
103
104 // Pointers of derived quantities are initialized to zero
105 set_der_0x0() ;
106
107}
108
109// Constructor from a file
110// -----------------------
111Bin_bhns::Bin_bhns(Map& mp_bh, Map& mp_ns, const Eos& eos_i, FILE* fich)
112 : ref_triad(0., "Absolute frame Cartesian basis"),
113 hole(mp_bh, fich),
114 star(mp_ns, eos_i, fich)
115{
116
117 // omega, separ, and x_rot are read from the file
118 fread_be(&omega, sizeof(double), 1, fich) ;
119 fread_be(&separ, sizeof(double), 1, fich) ;
120 fread_be(&x_rot, sizeof(double), 1, fich) ;
121 fread_be(&y_rot, sizeof(double), 1, fich) ;
122
123 // Pointers of derived quantities are initialized to zero
124 set_der_0x0() ;
125
126}
127
128
129 //--------------------//
130 // Destructor //
131 //--------------------//
132
134{
135
136 del_deriv() ;
137
138}
139
140
141 //------------------------------------------//
142 // Management of derived quantities //
143 //------------------------------------------//
144
146
147 if (p_mass_adm_bhns_surf != 0x0) delete p_mass_adm_bhns_surf ;
148 if (p_mass_adm_bhns_vol != 0x0) delete p_mass_adm_bhns_vol ;
149 if (p_mass_kom_bhns_surf != 0x0) delete p_mass_kom_bhns_surf ;
150 if (p_mass_kom_bhns_vol != 0x0) delete p_mass_kom_bhns_vol ;
151 if (p_line_mom_bhns != 0x0) delete p_line_mom_bhns ;
152 if (p_angu_mom_bhns != 0x0) delete p_angu_mom_bhns ;
153 if (p_virial_bhns_surf != 0x0) delete p_virial_bhns_surf ;
154 if (p_virial_bhns_vol != 0x0) delete p_virial_bhns_vol ;
155 if (p_xa_barycenter != 0x0) delete p_xa_barycenter ;
156 if (p_ya_barycenter != 0x0) delete p_ya_barycenter ;
157 if (p_omega_two_points != 0x0) delete p_omega_two_points ;
158 // if (p_ham_constr_bhns != 0x0) delete p_ham_constr_bhns ;
159 // if (p_mom_constr_bhns != 0x0) delete p_mom_constr_bhns ;
160
161 set_der_0x0() ;
162
163}
164
166
168 p_mass_adm_bhns_vol = 0x0 ;
170 p_mass_kom_bhns_vol = 0x0 ;
171 p_line_mom_bhns = 0x0 ;
172 p_angu_mom_bhns = 0x0 ;
173 p_virial_bhns_surf = 0x0 ;
174 p_virial_bhns_vol = 0x0 ;
175 p_xa_barycenter = 0x0 ;
176 p_ya_barycenter = 0x0 ;
177 p_omega_two_points = 0x0 ;
178 // p_ham_constr_bhns = 0x0 ;
179 // p_mom_constr_bhns = 0x0 ;
180
181}
182
183
184 //--------------------//
185 // Assignment //
186 //--------------------//
187
188// Assignment to anothe Bin_bhns
189// -----------------------------
190void Bin_bhns::operator=(const Bin_bhns& bibi) {
191
192 assert( bibi.ref_triad == ref_triad ) ;
193
194 hole = bibi.hole ;
195 star = bibi.star ;
196
197 omega = bibi.omega ;
198 separ = bibi.separ ;
199 x_rot = bibi.x_rot ;
200 y_rot = bibi.y_rot ;
201
202 del_deriv() ; // Deletes all derived quantities
203
204}
205
206
207 //-----------------//
208 // Outputs //
209 //-----------------//
210
211// Save in a file
212// --------------
213void Bin_bhns::sauve(FILE* fich) const {
214
215 hole.sauve(fich) ;
216 star.sauve(fich) ;
217
218 fwrite_be(&omega, sizeof(double), 1, fich) ;
219 fwrite_be(&separ, sizeof(double), 1, fich) ;
220 fwrite_be(&x_rot, sizeof(double), 1, fich) ;
221 fwrite_be(&y_rot, sizeof(double), 1, fich) ;
222
223}
224
225// Printing
226// --------
227ostream& operator<<(ostream& ost, const Bin_bhns& bibi) {
228
229 bibi >> ost ;
230 return ost ;
231
232}
233
234ostream& Bin_bhns::operator>>(ostream& ost) const {
235
236 using namespace Unites ;
237
238 ost << endl ;
239 ost << "BH-NS binary system" << endl ;
240 ost << "===================" << endl ;
241
242 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 ;
249
250 ost << endl << "Black hole : " ;
251 ost << endl << "========== " << endl ;
252
253 int nt = (hole.get_mp()).get_mg()->get_nt(1) ;
254
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 ;
271
272 ost << endl << "Neutron star : " ;
273 ost << endl << "============ " << endl ;
274
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. : "
280 << star.mass_b_bhns(hole.is_kerrschild(),hole.get_mass_bh(),separ)/msol
281 << " [Mo]" << endl ;
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 ;
304
305
306 return ost ;
307
308}
309
310// Display in polytropic unites
311// ----------------------------
312void Bin_bhns::display_poly(ostream& ost) const {
313
314 using namespace Unites ;
315
316 const Eos* p_eos = &( star.get_eos() ) ;
317 const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos ) ;
318
319 if (p_eos_poly != 0x0) {
320
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) ) ;
324
325 // Polytropic unit of length in terms of r_unit :
326 double r_poly = kap_ns2 / sqrt(ggrav) ;
327
328 // Polytropic unit of time in terms of t_unit :
329 double t_poly = r_poly ;
330
331 // Polytropic unit of mass in terms of m_unit :
332 double m_poly = r_poly / ggrav ;
333
334 // Polytropic unit of angular momentum in terms of j_unit :
335 // double j_poly = r_poly * r_poly / ggrav ;
336
337 double r0 = 0.5 * ( star.ray_eq() + star.ray_eq_pi() ) ;
338 // = 1 in Baumgarte et al.
339 // double d_ns = separ + star.ray_eq() - r0 ;
340 // Orbital separation of Baumgarte et al.
341
342 double y_ns = star.get_mp().get_ori_y() ;
343 double d_coord = sqrt(separ*separ + y_ns*y_ns) ;
344
345 ost.precision(16) ;
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
359 << endl ;
360 ost << " M_bar(NS) : "
361 << star.mass_b_bhns(hole.is_kerrschild(),hole.get_mass_bh(),separ)
362 / m_poly << endl ;
363 ost << " M_g(NS)iso : " << star.mass_g_bhns() / m_poly << endl ;
364 ost << " R_0(NS) : " << r0 / r_poly << endl ;
365
366 }
367
368}
369}
void display_poly(ostream &) const
Display in polytropic units.
Definition bin_bhns.C:312
void del_deriv() const
Deletes all the derived quantities.
Definition bin_bhns.C:145
double * p_mass_adm_bhns_surf
Total ADM mass of the system calculated by the surface integral at infinity.
Definition bin_bhns.h:97
friend ostream & operator<<(ostream &, const Bin_bhns &)
Display.
Definition bin_bhns.C:227
virtual void sauve(FILE *) const
Save in a file.
Definition bin_bhns.C:213
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.
Definition bin_bhns.C:73
double y_rot
Absolute Y coordinate of the rotation axis.
Definition bin_bhns.h:89
const Base_vect_cart ref_triad
Cartesian triad of the absolute reference frame.
Definition bin_bhns.h:69
double * p_mass_kom_bhns_surf
Total Komar mass of the system calculated by the surface integral at infinity.
Definition bin_bhns.h:107
Hole_bhns hole
Black hole.
Definition bin_bhns.h:72
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition bin_bhns.C:165
double omega
Angular velocity with respect to an asymptotically inertial observer.
Definition bin_bhns.h:80
double * p_mass_adm_bhns_vol
Total ADM mass of the system calculated by the volume integral and the surface integral at the appare...
Definition bin_bhns.h:102
double * p_virial_bhns_vol
Virial theorem error calculated by the ADM mass and the Komar mass of the volume integral.
Definition bin_bhns.h:128
double * p_omega_two_points
Orbital angular velocity derived from another method.
Definition bin_bhns.h:137
double separ
Absolute orbital separation between two centers of BH and NS.
Definition bin_bhns.h:83
Tbl * p_angu_mom_bhns
Total angular momentum of the system.
Definition bin_bhns.h:118
double * p_ya_barycenter
Absolute coordinate Y of the barycenter of the baryon density.
Definition bin_bhns.h:134
double x_rot
Absolute X coordinate of the rotation axis.
Definition bin_bhns.h:86
virtual ~Bin_bhns()
Destructor.
Definition bin_bhns.C:133
void operator=(const Bin_bhns &)
Assignment to another Bin_bhns.
Definition bin_bhns.C:190
double * p_virial_bhns_surf
Virial theorem error calculated by the ADM mass and the Komar mass of the surface integral at infinit...
Definition bin_bhns.h:123
double * p_mass_kom_bhns_vol
Total Komar mass of the system calculated by the volume integral and the surface integral at the appa...
Definition bin_bhns.h:112
Star_bhns star
Neutron star.
Definition bin_bhns.h:75
double * p_xa_barycenter
Absolute coordinate X of the barycenter of the baryon density.
Definition bin_bhns.h:131
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<).
Definition bin_bhns.C:234
Tbl * p_line_mom_bhns
Total linear momentum of the system.
Definition bin_bhns.h:115
Polytropic equation of state (relativistic case).
Definition eos.h:809
double get_gam() const
Returns the adiabatic index (cf. Eq. (3)).
Definition eos_poly.C:271
double get_kap() const
Returns the pressure coefficient (cf.
Definition eos_poly.C:275
Equation of state base class.
Definition eos.h:206
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
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
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
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142
Standard units of space, time and mass.