LORENE
binary_xcts.C
1/*
2 * Methods of class Binary_xcts
3 * (see file binary_xcts.h for documentation)
4 */
5
6/*
7 * Copyright (c) 2010 Michal Bejger
8 *
9 * This file is part of LORENE.
10 *
11 * LORENE is free software; you can redistribute it and/or modify
12 * it under the terms of the GNU General Public License version 2
13 * as published by the Free Software Foundation.
14 *
15 * LORENE is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with LORENE; if not, write to the Free Software
22 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23 *
24 */
25
26
27
28/*
29 * $Id: binary_xcts.C,v 1.7 2016/12/05 16:17:47 j_novak Exp $
30 * $Log: binary_xcts.C,v $
31 * Revision 1.7 2016/12/05 16:17:47 j_novak
32 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
33 *
34 * Revision 1.6 2014/10/13 08:52:45 j_novak
35 * Lorene classes and functions now belong to the namespace Lorene.
36 *
37 * Revision 1.5 2014/10/06 15:12:59 j_novak
38 * Modified #include directives to use c++ syntax.
39 *
40 * Revision 1.4 2010/12/20 09:56:02 m_bejger
41 * Pointer to the linear momentum added
42 *
43 * Revision 1.3 2010/12/09 10:38:46 m_bejger
44 * virial_vol() added, fait_decouple() removed
45 *
46 * Revision 1.2 2010/10/28 12:00:07 m_bejger
47 * Mass-shedding indicators added to the output in Binary_xcts::write_global
48 *
49 * Revision 1.1 2010/05/04 07:35:54 m_bejger
50 * Initial version
51 *
52 * $Header: /cvsroot/Lorene/C++/Source/Binary_xcts/binary_xcts.C,v 1.7 2016/12/05 16:17:47 j_novak Exp $
53 *
54 */
55
56// Headers C
57#include <cmath>
58
59// Headers Lorene
60#include "binary_xcts.h"
61#include "eos.h"
62#include "utilitaires.h"
63#include "graphique.h"
64#include "param.h"
65#include "unites.h"
66
67 //--------------//
68 // Constructors //
69 //--------------//
70
71// Standard constructor
72// --------------------
73
74namespace Lorene {
76 int nzet1,
77 const Eos& eos1,
78 int irrot1,
79 Map& mp2,
80 int nzet2,
81 const Eos& eos2,
82 int irrot2)
83 : star1(mp1, nzet1, eos1, irrot1),
84 star2(mp2, nzet2, eos2, irrot2)
85{
86
87 et[0] = &star1 ;
88 et[1] = &star2 ;
89
90 omega = 0 ;
91 x_axe = 0 ;
92
93 // Pointers of derived quantities initialized to zero :
94 set_der_0x0() ;
95}
96
97// Copy constructor
98// ----------------
100 : star1(bibi.star1),
101 star2(bibi.star2),
102 omega(bibi.omega),
103 x_axe(bibi.x_axe)
104{
105 et[0] = &star1 ;
106 et[1] = &star2 ;
107
108 // Pointers of derived quantities initialized to zero :
109 set_der_0x0() ;
110}
111
112// Constructor from a file
113// -----------------------
115 const Eos& eos1,
116 Map& mp2,
117 const Eos& eos2,
118 FILE* fich)
119 : star1(mp1, eos1, fich),
120 star2(mp2, eos2, fich)
121{
122
123 et[0] = &star1 ;
124 et[1] = &star2 ;
125
126 // omega and x_axe are read in the file:
127 fread_be(&omega, sizeof(double), 1, fich) ;
128 fread_be(&x_axe, 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 //----------------------------------//
146 // Management of derived quantities //
147 //----------------------------------//
148
150
151 if (p_mass_adm != 0x0) delete p_mass_adm ;
152 if (p_mass_kom != 0x0) delete p_mass_kom ;
153 if (p_angu_mom != 0x0) delete p_angu_mom ;
154 if (p_lin_mom != 0x0) delete p_lin_mom ;
155 if (p_total_ener != 0x0) delete p_total_ener ;
156 if (p_virial != 0x0) delete p_virial ;
157 if (p_virial_vol != 0x0) delete p_virial_vol ;
158
159 set_der_0x0() ;
160}
161
162
163
164
166
167 p_mass_adm = 0x0 ;
168 p_mass_kom = 0x0 ;
169 p_angu_mom = 0x0 ;
170 p_lin_mom = 0x0 ;
171 p_total_ener = 0x0 ;
172 p_virial = 0x0 ;
173 p_virial_vol = 0x0 ;
174
175}
176
177
178 //--------------//
179 // Assignment //
180 //--------------//
181
182// Assignment to another Binary_xcts
183// --------------------------------
184
186
187 star1 = bibi.star1 ;
188 star2 = bibi.star2 ;
189
190 omega = bibi.omega ;
191 x_axe = bibi.x_axe ;
192
193 del_deriv() ; // Deletes all derived quantities
194
195}
196
197 //--------------//
198 // Outputs //
199 //--------------//
200
201// Save in a file
202// --------------
203void Binary_xcts::sauve(FILE* fich) const {
204
205 star1.sauve(fich) ;
206 star2.sauve(fich) ;
207
208 fwrite_be(&omega, sizeof(double), 1, fich) ;
209 fwrite_be(&x_axe, sizeof(double), 1, fich) ;
210
211}
212
213// Printing
214// --------
215ostream& operator<<(ostream& ost, const Binary_xcts& bibi) {
216
217 bibi >> ost ;
218 return ost ;
219}
220
221
222ostream& Binary_xcts::operator>>(ostream& ost) const {
223
224 using namespace Unites ;
225
226 ost << endl ;
227 ost << "Binary neutron stars" << endl ;
228 ost << "=============" << endl ;
229 ost << endl <<
230 "Orbital angular velocity : " << omega * f_unit << " rad/s" << endl ;
231 ost << endl <<
232 "Coordinate separation between the two stellar centers : "
233 << separation() / km << " km" << endl ;
234 ost <<
235 "Absolute coordinate X of the rotation axis : " << x_axe / km
236 << " km" << endl ;
237 ost << endl << "Star 1 : " << endl ;
238 ost << "====== " << endl ;
239 ost << star1 << endl ;
240 ost << "Star 2 : " << endl ;
241 ost << "====== " << endl ;
242 ost << star2 << endl ;
243 return ost ;
244}
245
246// Display in polytropic units
247// ---------------------------
248
249void Binary_xcts::display_poly(ostream& ost) const {
250
251 using namespace Unites ;
252
253 const Eos* p_eos1 = &( star1.get_eos() ) ;
254 const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;
255
256 if (p_eos_poly != 0x0) {
257
258 assert( star1.get_eos() == star2.get_eos() ) ;
259
260 double kappa = p_eos_poly->get_kap() ;
261 double gamma = p_eos_poly->get_gam() ; ;
262 double kap_ns2 = pow( kappa, 0.5 /(gamma-1) ) ;
263
264 // Polytropic unit of length in terms of r_unit :
265 double r_poly = kap_ns2 / sqrt(ggrav) ;
266
267 // Polytropic unit of time in terms of t_unit :
268 double t_poly = r_poly ;
269
270 // Polytropic unit of mass in terms of m_unit :
271 double m_poly = r_poly / ggrav ;
272
273 // Polytropic unit of angular momentum in terms of j_unit :
274 double j_poly = r_poly * r_poly / ggrav ;
275
276 ost.precision(10) ;
277 ost << endl << "Quantities in polytropic units : " << endl ;
278 ost << "==============================" << endl ;
279 ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
280 ost << " d_e_max : " << separation() / r_poly << endl ;
281 ost << " d_G : "
282 << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly
283 << endl ;
284 ost << " Omega : " << omega * t_poly << endl ;
285 ost << " J : " << angu_mom()(2) / j_poly << endl ;
286 ost << " M_ADM : " << mass_adm() / m_poly << endl ;
287 ost << " M_Komar : " << mass_kom() / m_poly << endl ;
288 ost << " M_bar(star 1) : " << star1.mass_b() / m_poly << endl ;
289 ost << " M_bar(star 2) : " << star2.mass_b() / m_poly << endl ;
290 ost << " R_0(star 1) : " <<
291 0.5 * ( star1.ray_eq() + star1.ray_eq_pi() ) / r_poly << endl ;
292 ost << " R_0(star 2) : " <<
293 0.5 * ( star2.ray_eq() + star2.ray_eq_pi() ) / r_poly << endl ;
294
295 }
296
297}
298
299void Binary_xcts::write_global(ostream& ost) const {
300
301 using namespace Unites ;
302
303 const Map& mp1 = star1.get_mp() ;
304 const Mg3d* mg1 = mp1.get_mg() ;
305 int nz1 = mg1->get_nzone() ;
306
307 // Mass-shedding indicators
308 double dent1_eq = (star1.ent).dsdr().val_point(star1.ray_eq(),M_PI/2.,0.) ;
309 double dent1_pole = (star1.ent).dsdr().val_point(star1.ray_pole(),0.,0.) ;
310 double chi1 = fabs( dent1_eq / dent1_pole ) ;
311
312 double dent2_eq = (star2.ent).dsdr().val_point(star2.ray_eq(),M_PI/2.,0.) ;
313 double dent2_pole = (star2.ent).dsdr().val_point(star2.ray_pole(),0.,0.) ;
314 double chi2 = fabs( dent2_eq / dent2_pole ) ;
315
316 ost.precision(5) ;
317 ost << "# Grid 1 : " << nz1 << "x"
318 << mg1->get_nr(0) << "x" << mg1->get_nt(0) << "x" << mg1->get_np(0)
319 << " R_out(l) [km] : " ;
320 for (int l=0; l<nz1; l++) {
321 ost << " " << mp1.val_r(l, 1., M_PI/2, 0) / km ;
322 }
323 ost << endl ;
324
325 ost << "# VE(M) VE(M) [vol]" << endl ;
326
327
328 ost.setf(ios::scientific) ;
329 ost.width(14) ;
330 ost << virial() << " " << virial_vol() << endl ;
331
332 ost << "# d [km] "
333 << " d_G [km] "
334 << " d/(a1 +a1') "
335 << " f [Hz] "
336 << " M_ADM [M_sol] "
337 << " M_ADM_vol [M_sol] "
338 << " M_Komar [M_sol] "
339 << " M_Komar_vol [M_sol] "
340 << " J [G M_sol^2/c] " << endl ;
341
342 ost.precision(14) ;
343 ost.width(20) ;
344 ost << separation() / km ; ost.width(22) ;
345 ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / km ; ost.width(22) ;
346 ost << separation() / (star1.ray_eq() + star2.ray_eq()) ; ost.width(22) ;
347 ost << omega / (2*M_PI)* f_unit ; ost.width(22) ;
348 ost << mass_adm() / msol ; ost.width(22) ;
349 ost << mass_adm_vol() / msol ; ost.width(22) ;
350 ost << mass_kom() / msol ; ost.width(22) ;
351 ost << mass_kom_vol() / msol ; ost.width(22) ;
352 ost << angu_mom()(2)/ ( qpig / (4* M_PI) * msol*msol) << endl ;
353
354 ost << "# H_c(1)[c^2] "
355 << " e_c(1)[rho_nuc] "
356 << " M_B(1) [M_sol] "
357 << " r_eq(1) [km] "
358 << " a2/a1(1) "
359 << " a3/a1(1) "
360 << " chi1 " << endl ;
361
362 ost.width(20) ;
363 ost << star1.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
364 ost << star1.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
365 ost << star1.mass_b() / msol ; ost.width(22) ;
366 ost << star1.ray_eq() / km ; ost.width(22) ;
367 ost << star1.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
368 ost << star1.ray_pole() / star1.ray_eq() ; ost.width(22) ;
369 ost << chi1 << endl ;
370
371 ost << "# H_c(2)[c^2] "
372 << " e_c(2)[rho_nuc] "
373 << " M_B(2) [M_sol] "
374 << " r_eq(2) [km] "
375 << " a2/a1(2) "
376 << " a3/a1(2) "
377 << " chi2 " << endl ;
378
379
380 ost.width(20) ;
381 ost << star2.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
382 ost << star2.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
383 ost << star2.mass_b() / msol ; ost.width(22) ;
384 ost << star2.ray_eq() / km ; ost.width(22) ;
385 ost << star2.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
386 ost << star2.ray_pole() / star1.ray_eq() ; ost.width(22) ;
387 ost << chi2 << endl ;
388
389 // Quantities in polytropic units if the EOS is a polytropic one
390 // -------------------------------------------------------------
391 const Eos* p_eos1 = &( star1.get_eos() ) ;
392 const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;
393
394 if ((p_eos_poly != 0x0) && ( star1.get_eos() == star2.get_eos() )) {
395
396 double kappa = p_eos_poly->get_kap() ;
397 double gamma = p_eos_poly->get_gam() ; ;
398 double kap_ns2 = pow( kappa, 0.5 /(gamma-1.) ) ;
399
400 // Polytropic unit of length in terms of r_unit :
401 double r_poly = kap_ns2 / sqrt(ggrav) ;
402
403 // Polytropic unit of time in terms of t_unit :
404 double t_poly = r_poly ;
405
406 // Polytropic unit of mass in terms of m_unit :
407 double m_poly = r_poly / ggrav ;
408
409 // Polytropic unit of angular momentum in terms of j_unit :
410 double j_poly = r_poly * r_poly / ggrav ;
411
412 ost << "# d [poly] "
413 << " d_G [poly] "
414 << " Omega [poly] "
415 << " M_ADM [poly] "
416 << " J [poly] "
417 << " M_B(1) [poly] "
418 << " M_B(2) [poly] " << endl ;
419
420 ost.width(20) ;
421 ost << separation() / r_poly ; ost.width(22) ;
422 ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly ; ost.width(22) ;
423 ost << omega * t_poly ; ost.width(22) ;
424 ost << mass_adm() / m_poly ; ost.width(22) ;
425 ost << angu_mom()(2) / j_poly ; ost.width(22) ;
426 ost << star1.mass_b() / m_poly ; ost.width(22) ;
427 ost << star2.mass_b() / m_poly << endl ;
428
429 }
430
431}
432
433 //-------------------------------//
434 // Miscellaneous //
435 //-------------------------------//
436
438
439 double dx = star1.mp.get_ori_x() - star2.mp.get_ori_x() ;
440 double dy = star1.mp.get_ori_y() - star2.mp.get_ori_y() ;
441 double dz = star1.mp.get_ori_z() - star2.mp.get_ori_z() ;
442
443 return sqrt( dx*dx + dy*dy + dz*dz ) ;
444
445}
446}
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
double virial_vol() const
Estimates the relative error on the virial theorem (volume version).
double * p_mass_kom
Total Komar mass of the system.
Definition binary_xcts.h:94
double x_axe
Absolute X coordinate of the rotation axis.
Definition binary_xcts.h:84
void write_global(ostream &) const
Write global quantities in a formatted file.
Star_bin_xcts star2
Second star of the system.
Definition binary_xcts.h:69
double mass_adm_vol() const
Total ADM mass (computed by a volume integral).
double * p_virial_vol
Virial theorem error (volume version).
Tbl * p_lin_mom
Total linear momentum of the system.
void operator=(const Binary_xcts &)
Assignment to another Binary_xcts.
double mass_kom_vol() const
Total Komar mass (computed by a volume integral).
Star_bin_xcts star1
First star of the system.
Definition binary_xcts.h:66
void sauve(FILE *) const
Save in a file.
double mass_adm() const
Total ADM mass.
double separation() const
Returns the coordinate separation of the two stellar centers [r_unit].
void del_deriv() const
Deletes all the derived quantities.
Binary_xcts(Map &mp1, int nzet1, const Eos &eos1, int irrot1, Map &mp2, int nzet2, const Eos &eos2, int irrot2)
Standard constructor.
Definition binary_xcts.C:75
double * p_mass_adm
Total ADM mass of the system.
Definition binary_xcts.h:91
double * p_total_ener
Total energy of the system.
friend ostream & operator<<(ostream &, const Binary_xcts &)
Display.
void display_poly(ostream &) const
Display in polytropic units.
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<).
double mass_kom() const
Total Komar mass.
~Binary_xcts()
Destructor.
double virial() const
Estimates the relative error on the virial theorem.
Tbl * p_angu_mom
Total angular momentum of the system.
Definition binary_xcts.h:97
Star_bin_xcts * et[2]
Array of the two stars (to perform loops on the stars): et[0] contains the address of star1 and et[1]...
Definition binary_xcts.h:75
double * p_virial
Virial theorem error.
const Tbl & angu_mom() const
Total angular momentum.
double omega
Angular velocity with respect to an asymptotically inertial observer.
Definition binary_xcts.h:80
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
Multi-domain grid.
Definition grilles.h:279
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:479
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:474
int get_nzone() const
Returns the number of domains.
Definition grilles.h:465
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:469
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
virtual void dsdr(const Cmp &ci, Cmp &resu) const =0
Computes of a Cmp .
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142
Standard units of space, time and mass.