File: binary.C

package info (click to toggle)
lorene 0.0.0~cvs20161116%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 26,472 kB
  • sloc: cpp: 212,946; fortran: 21,645; makefile: 1,750; sh: 4
file content (613 lines) | stat: -rw-r--r-- 18,115 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
/*
 * Methods of class Binary
 *
 */

/*
 *   Copyright (c) 2004 Francois Limousin
 *
 *   This file is part of LORENE.
 *
 *   LORENE is free software; you can redistribute it and/or modify
 *   it under the terms of the GNU General Public License as published by
 *   the Free Software Foundation; either version 2 of the License, or
 *   (at your option) any later version.
 *
 *   LORENE is distributed in the hope that it will be useful,
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *   GNU General Public License for more details.
 *
 *   You should have received a copy of the GNU General Public License
 *   along with LORENE; if not, write to the Free Software
 *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
 *
 */
char Binary_C[] = "$Header: /cvsroot/Lorene/C++/Source/Binary/binary.C,v 1.16 2014/10/13 08:52:44 j_novak Exp $" ;

/*
 * $Id: binary.C,v 1.16 2014/10/13 08:52:44 j_novak Exp $
 * $Log: binary.C,v $
 * Revision 1.16  2014/10/13 08:52:44  j_novak
 * Lorene classes and functions now belong to the namespace Lorene.
 *
 * Revision 1.15  2014/10/06 15:12:59  j_novak
 * Modified #include directives to use c++ syntax.
 *
 * Revision 1.14  2005/09/15 14:39:14  e_gourgoulhon
 * Added printing of angular momentum in display_poly.
 *
 * Revision 1.13  2005/09/13 19:38:31  f_limousin
 * Reintroduction of the resolution of the equations in cartesian coordinates.
 *
 * Revision 1.12  2005/02/24 17:31:27  f_limousin
 * Update of the function decouple().
 *
 * Revision 1.11  2005/02/18 13:14:06  j_novak
 * Changing of malloc/free to new/delete + suppression of some unused variables
 * (trying to avoid compilation warnings).
 *
 * Revision 1.10  2004/07/21 11:47:10  f_limousin
 * Add / Delete comments.
 *
 * Revision 1.9  2004/05/25 14:25:12  f_limousin
 * Add the virial theorem for conformally flat configurations.
 *
 * Revision 1.8  2004/03/25 10:29:01  j_novak
 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
 *
 * Revision 1.7  2004/03/23 10:00:47  f_limousin
 * Minor changes.
 *
 * Revision 1.6  2004/02/27 09:59:33  f_limousin
 * Modification in the routine decouple().
 *
 * Revision 1.5  2004/02/21 17:05:12  e_gourgoulhon
 * Method Scalar::point renamed Scalar::val_grid_point.
 * Method Scalar::set_point renamed Scalar::set_grid_point.
 *
 * Revision 1.4  2004/01/20 15:21:12  f_limousin
 * First version
 *
 *
 * $Header: /cvsroot/Lorene/C++/Source/Binary/binary.C,v 1.16 2014/10/13 08:52:44 j_novak Exp $
 *
 */

// Headers C
#include <cmath>

// Headers Lorene
#include "binary.h"
#include "eos.h"
#include "utilitaires.h"
#include "graphique.h"
#include "param.h"
#include "unites.h"	    

			    //--------------//
			    // Constructors //
			    //--------------//

// Standard constructor
// --------------------

namespace Lorene {
Binary::Binary(Map& mp1, int nzet1, const Eos& eos1, int irrot1, 
	       Map& mp2, int nzet2, const Eos& eos2, int irrot2, 
	       int conf_flat) 
                 : star1(mp1, nzet1, eos1, irrot1, conf_flat), 
		   star2(mp2, nzet2, eos2, irrot2, conf_flat)
{

    et[0] = &star1 ; 
    et[1] = &star2 ; 
    
    omega = 0 ; 
    x_axe = 0 ; 

    // Pointers of derived quantities initialized to zero :
    set_der_0x0() ;
}

// Copy constructor
// ----------------
Binary::Binary(const Binary& bibi) 
		: star1(bibi.star1), 
		  star2(bibi.star2),
		  omega(bibi.omega), 
		  x_axe(bibi.x_axe) 
{
    et[0] = &star1 ; 
    et[1] = &star2 ; 

    // Pointers of derived quantities initialized to zero :
    set_der_0x0() ;    
}

// Constructor from a file
// -----------------------
Binary::Binary(Map& mp1, const Eos& eos1, Map& mp2, const Eos& eos2, 
	       FILE* fich)
		: star1(mp1, eos1, fich), 
		  star2(mp2, eos2, fich) 
{
    et[0] = &star1 ; 
    et[1] = &star2 ; 

    // omega and x_axe are read in the file:
    fread_be(&omega, sizeof(double), 1, fich) ;		
    fread_be(&x_axe, sizeof(double), 1, fich) ;		

    // Pointers of derived quantities initialized to zero : 
    set_der_0x0() ;    
    
}			

			    //------------//
			    // Destructor //
			    //------------//

Binary::~Binary(){

    del_deriv() ; 

}

			//----------------------------------//
			// Management of derived quantities //
			//----------------------------------//

void Binary::del_deriv() const {

    if (p_mass_adm != 0x0) delete p_mass_adm ; 
    if (p_mass_kom != 0x0) delete p_mass_kom ; 
    if (p_angu_mom != 0x0) delete p_angu_mom ; 
    if (p_total_ener != 0x0) delete p_total_ener ; 
    if (p_virial != 0x0) delete p_virial ; 
    if (p_ham_constr != 0x0) delete p_ham_constr ; 
    if (p_mom_constr != 0x0) delete p_mom_constr ; 

    set_der_0x0() ; 
}			    




void Binary::set_der_0x0() const {

    p_mass_adm = 0x0 ; 
    p_mass_kom = 0x0 ; 
    p_angu_mom = 0x0 ; 
    p_total_ener = 0x0 ; 
    p_virial = 0x0 ; 
    p_ham_constr = 0x0 ; 
    p_mom_constr = 0x0 ; 

}			    


			    //--------------//
			    //  Assignment  //
			    //--------------//

// Assignment to another Binary
// --------------------------------

void Binary::operator=(const Binary& bibi) {

    star1 = bibi.star1 ; 
    star2 = bibi.star2 ; 
    
    omega = bibi.omega ; 
    x_axe = bibi.x_axe ; 
    
    del_deriv() ;  // Deletes all derived quantities
    
}

			    //--------------//
			    //	  Outputs   //
			    //--------------//

// Save in a file
// --------------
void Binary::sauve(FILE* fich) const {
    
    star1.sauve(fich) ; 
    star2.sauve(fich) ; 
    
    fwrite_be(&omega, sizeof(double), 1, fich) ;		
    fwrite_be(&x_axe, sizeof(double), 1, fich) ;		
    
}

// Printing
// --------
ostream& operator<<(ostream& ost, const Binary& bibi)  {
    bibi >> ost ;
    return ost ;
}
    

ostream& Binary::operator>>(ostream& ost) const {

  using namespace Unites ;

    ost << endl ; 
    ost << "Binary neutron stars" << endl ; 
    ost << "=============" << endl ; 
    ost << endl << 
	"Orbital angular velocity : " << omega * f_unit << " rad/s" << endl ; 
    ost << endl << 
	"Coordinate separation between the two stellar centers : " 
	<< separation() / km  << " km" << endl ; 
    ost << 
	"Absolute coordinate X of the rotation axis : " << x_axe / km 
	    << " km" << endl ; 
    ost << endl << "Star 1 : " << endl ; 
    ost << "======   " << endl ; 
    ost << star1 << endl ; 
    ost << "Star 2 : " << endl ; 
    ost << "======   " << endl ; 
    ost << star2 << endl ; 
    return ost ;
}

// Display in polytropic units
// ---------------------------

void Binary::display_poly(ostream& ost) const {

  using namespace Unites ;

    const Eos* p_eos1 = &( star1.get_eos() ) ; 
    const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ; 	  

    if (p_eos_poly != 0x0) {

	assert( star1.get_eos() == star2.get_eos() ) ; 

	double kappa = p_eos_poly->get_kap() ; 
	double gamma = p_eos_poly->get_gam() ;  ; 
	double kap_ns2 = pow( kappa,  0.5 /(gamma-1) ) ; 
    
	// Polytropic unit of length in terms of r_unit : 
	double r_poly = kap_ns2 / sqrt(ggrav) ; 
    
	// Polytropic unit of time in terms of t_unit :
	double t_poly = r_poly ; 

	// Polytropic unit of mass in terms of m_unit :
	double m_poly = r_poly / ggrav ; 
    
	// Polytropic unit of angular momentum in terms of j_unit :
	double j_poly = r_poly * r_poly / ggrav ; 
    
	ost.precision(10) ; 
	ost << endl << "Quantities in polytropic units : " << endl ; 
	ost	 << "==============================" << endl ; 
	ost << " ( r_poly = " << r_poly / km << " km )" << endl ; 
	ost << "  d_e_max	: " << separation() / r_poly << endl ; 
	ost << "  d_G		: " 
	     << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly 
	     << endl ; 
	ost << "  Omega	  : " << omega * t_poly << endl ; 
	ost << "  J	  : " << angu_mom()(2) / j_poly << endl ; 
	ost << "  M_ADM   : " << mass_adm() / m_poly << endl ;      
	ost << "  M_Komar : " << mass_kom() / m_poly << endl ; 
	ost << "  M_bar(star 1) : " << star1.mass_b() / m_poly << endl ; 
	ost << "  M_bar(star 2) : " << star2.mass_b() / m_poly << endl ; 
	ost << "  R_0(star 1)	: " << 
	0.5 * ( star1.ray_eq() + star1.ray_eq_pi() ) / r_poly << endl ;  
	ost << "  R_0(star 2)	: " << 
	0.5 * ( star2.ray_eq() + star2.ray_eq_pi() ) / r_poly << endl ;  
    
    }
    

} 


void Binary::fait_decouple () {

    int nz_un = star1.get_mp().get_mg()->get_nzone() ;
    int nz_deux = star2.get_mp().get_mg()->get_nzone() ;
    
    // We determine R_limite :
    double distance = star2.get_mp().get_ori_x() - star1.get_mp().get_ori_x() ;
    cout << "distance = " << distance << endl ;
    double lim_un = distance/2. ;
    double lim_deux = distance/2. ;
    double int_un = distance/6. ;
    double int_deux = distance/6. ;
    
    // The functions used.
    Scalar fonction_f_un (star1.get_mp()) ;
    fonction_f_un = 0.5*pow(
	cos((star1.get_mp().r-int_un)*M_PI/2./(lim_un-int_un)), 2.)+0.5 ;
    fonction_f_un.std_spectral_base();
    
    Scalar fonction_g_un (star1.get_mp()) ;
    fonction_g_un = 0.5*pow
	(sin((star1.get_mp().r-int_un)*M_PI/2./(lim_un-int_un)), 2.) ;
    fonction_g_un.std_spectral_base();
    
    Scalar fonction_f_deux (star2.get_mp()) ;
    fonction_f_deux = 0.5*pow(
	cos((star2.get_mp().r-int_deux)*M_PI/2./(lim_deux-int_deux)), 2.)+0.5 ;
    fonction_f_deux.std_spectral_base();
    
    Scalar fonction_g_deux (star2.get_mp()) ;
    fonction_g_deux = 0.5*pow(
	sin((star2.get_mp().r-int_deux)*M_PI/2./(lim_un-int_deux)), 2.) ;
    fonction_g_deux.std_spectral_base();
    
    // The functions total :
    Scalar decouple_un (star1.get_mp()) ;
    decouple_un.allocate_all() ;
    Scalar decouple_deux (star2.get_mp()) ;
    decouple_deux.allocate_all() ;
    
    Mtbl xabs_un (star1.get_mp().xa) ;
    Mtbl yabs_un (star1.get_mp().ya) ;
    Mtbl zabs_un (star1.get_mp().za) ;
	    
    Mtbl xabs_deux (star2.get_mp().xa) ;
    Mtbl yabs_deux (star2.get_mp().ya) ;
    Mtbl zabs_deux (star2.get_mp().za) ;
	    
    double xabs, yabs, zabs, air_un, air_deux, theta, phi ;
	    
    for (int l=0 ; l<nz_un ; l++) {
	int nr = star1.get_mp().get_mg()->get_nr (l) ;
		
	if (l==nz_un-1)
	    nr -- ;
		
	int np = star1.get_mp().get_mg()->get_np (l) ;
	int nt = star1.get_mp().get_mg()->get_nt (l) ;
		
	for (int k=0 ; k<np ; k++)
	    for (int j=0 ; j<nt ; j++)
		for (int i=0 ; i<nr ; i++) {
			    
		    xabs = xabs_un (l, k, j, i) ;
		    yabs = yabs_un (l, k, j, i) ;
		    zabs = zabs_un (l, k, j, i) ;
			    
		    // Coordinates of the point
		    star1.get_mp().convert_absolute 
			(xabs, yabs, zabs, air_un, theta, phi) ;
		    star2.get_mp().convert_absolute 
			(xabs, yabs, zabs, air_deux, theta, phi) ;
		
		    if (air_un <= lim_un)
			if (air_un < int_un)
			    decouple_un.set_grid_point(l, k, j, i) = 1 ;
			else
			// Close to star 1 :
			decouple_un.set_grid_point(l, k, j, i) = 
			    fonction_f_un.val_grid_point(l, k, j, i) ;
		    else 
			if (air_deux <= lim_deux)
			    if (air_deux < int_deux)
				decouple_un.set_grid_point(l, k, j, i) = 0 ;
			    else
			// Close to star 2 :
				decouple_un.set_grid_point(l, k, j, i) = 
		fonction_g_deux.val_point (air_deux, theta, phi) ;
		
			else
			    // Far from each star :
			    decouple_un.set_grid_point(l, k, j, i) = 0.5 ;
		}
			    
		// Case infinity :
		if (l==nz_un-1)
		    for (int k=0 ; k<np ; k++)
			for (int j=0 ; j<nt ; j++)
			    decouple_un.set_grid_point(nz_un-1, k, j, nr)=0.5 ;
	    }
    
    for (int l=0 ; l<nz_deux ; l++) {
	int nr = star2.get_mp().get_mg()->get_nr (l) ;
		
	if (l==nz_deux-1)
	    nr -- ;
		
	int np = star2.get_mp().get_mg()->get_np (l) ;
	int nt = star2.get_mp().get_mg()->get_nt (l) ;
		
	for (int k=0 ; k<np ; k++)
	    for (int j=0 ; j<nt ; j++)
		for (int i=0 ; i<nr ; i++) {
			    
		    xabs = xabs_deux (l, k, j, i) ;
		    yabs = yabs_deux (l, k, j, i) ;
		    zabs = zabs_deux (l, k, j, i) ;
			    
		    // coordinates of the point  :
		    star1.get_mp().convert_absolute 
			(xabs, yabs, zabs, air_un, theta, phi) ;
		    star2.get_mp().convert_absolute 
			(xabs, yabs, zabs, air_deux, theta, phi) ;
		    
		    if (air_deux <= lim_deux)
			if (air_deux < int_deux)
			    decouple_deux.set_grid_point(l, k, j, i) = 1 ;
			else
			// close to star two :
			decouple_deux.set_grid_point(l, k, j, i) = 
			    fonction_f_deux.val_grid_point(l, k, j, i) ;
		    else 
			if (air_un <= lim_un)
			    if (air_un < int_un)
				decouple_deux.set_grid_point(l, k, j, i) = 0 ;
			    else
			// close to star one :
				decouple_deux.set_grid_point(l, k, j, i) = 
			 fonction_g_un.val_point (air_un, theta, phi) ;
		
			else
			    // Far from each star :
			    decouple_deux.set_grid_point(l, k, j, i) = 0.5 ;
		}
			    
		// Case infinity :
		if (l==nz_deux-1)
		    for (int k=0 ; k<np ; k++)
			for (int j=0 ; j<nt ; j++)
			 decouple_deux.set_grid_point(nz_un-1, k, j, nr)=0.5 ;
   }
   
    decouple_un.std_spectral_base() ;
    decouple_deux.std_spectral_base() ;

    int nr = star1.get_mp().get_mg()->get_nr (1) ;
    int nt = star1.get_mp().get_mg()->get_nt (1) ;
    int np = star1.get_mp().get_mg()->get_np (1) ;
    cout << "decouple_un"  << endl << norme(decouple_un/(nr*nt*np)) << endl ;
    cout << "decouple_deux"  << endl << norme(decouple_deux/(nr*nt*np)) 
	 << endl ;

    star1.decouple = decouple_un ;
    star2.decouple = decouple_deux ;
}

void Binary::write_global(ostream& ost) const {

  using namespace Unites ;

  const Map& mp1 = star1.get_mp() ;
  const Mg3d* mg1 = mp1.get_mg() ;
  int nz1 = mg1->get_nzone() ; 

  ost.precision(5) ;
  ost << "# Grid 1 : " << nz1 << "x"
      << mg1->get_nr(0) << "x" << mg1->get_nt(0) << "x" << mg1->get_np(0) 
      << "  R_out(l) [km] : " ;
  for (int l=0; l<nz1; l++) {
    ost << " " << mp1.val_r(l, 1., M_PI/2, 0) / km ; 
  }
  ost << endl ; 

  ost << "#     VE(M)	   " << endl ;
  
  
  ost.setf(ios::scientific) ; 
  ost.width(14) ; 
  ost << virial() << endl ;
  
  ost << "#      d [km]         "  
      << "       d_G [km]       "
      << "     d/(a1 +a1')      "
      << "       f [Hz]         "
      << "    M_ADM [M_sol]     "     
      << "    M_ADM_vol [M_sol]     "     
      << "    M_Komar [M_sol]     "     
      << "    M_Komar_vol [M_sol]     "     
      << "   J [G M_sol^2/c]    "  << endl ;   
  
  ost.precision(14) ;
  ost.width(20) ; 
  ost << separation() / km ; ost.width(22) ;
  ost	<< ( star2.xa_barycenter() - star1.xa_barycenter() ) / km ; ost.width(22) ;
  ost	<< separation() / (star1.ray_eq() + star2.ray_eq()) ; ost.width(22) ;
  ost	<< omega / (2*M_PI)* f_unit ; ost.width(22) ;
  ost	<< mass_adm() / msol ; ost.width(22) ; 
  ost	<< mass_adm_vol() / msol ; ost.width(22) ; 
  ost	<< mass_kom() / msol ; ost.width(22) ; 
  ost	<< mass_kom_vol() / msol ; ost.width(22) ; 
  ost	<< angu_mom()(2)/ ( qpig / (4* M_PI) * msol*msol) << endl ; 
  
  ost << "#     H_c(1)[c^2]     "
      << "    e_c(1)[rho_nuc]   " 
      << "    M_B(1) [M_sol]    "
      << "     r_eq(1) [km]     "
      << "        a2/a1(1)	  " 
	    << "        a3/a1(1)	  " << endl ; 
  
  ost.width(20) ; 
  ost << star1.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
  ost	<< star1.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
  ost	<< star1.mass_b() / msol ; ost.width(22) ;	
  ost << star1.ray_eq() / km ; ost.width(22) ; 
  ost	<< star1.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
  ost	<< star1.ray_pole() / star1.ray_eq() << endl ;
  
  ost << "#     H_c(2)[c^2]     "
      << "    e_c(2)[rho_nuc]   " 
      << "    M_B(2) [M_sol]    "
      << "     r_eq(2) [km]     "
      << "        a2/a1(2)	  " 
      << "        a3/a1(2)	  " << endl ; 
  
  ost.width(20) ; 
  ost << star2.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
  ost	<< star2.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
  ost	<< star2.mass_b() / msol ; ost.width(22) ;	
  ost << star2.ray_eq() / km ; ost.width(22) ; 
  ost	<< star2.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
  ost	<< star2.ray_pole() / star1.ray_eq() << endl ;
  
  // Quantities in polytropic units if the EOS is a polytropic one
  // -------------------------------------------------------------
  const Eos* p_eos1 = &( star1.get_eos() ) ; 
  const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ; 	  
  
  if ((p_eos_poly != 0x0) && ( star1.get_eos() == star2.get_eos() )) {
      
      double kappa = p_eos_poly->get_kap() ; 
      double gamma = p_eos_poly->get_gam() ;  ; 
      double kap_ns2 = pow( kappa,  0.5 /(gamma-1.) ) ; 
      
      // Polytropic unit of length in terms of r_unit : 
      double r_poly = kap_ns2 / sqrt(ggrav) ; 
      
      // Polytropic unit of time in terms of t_unit :
      double t_poly = r_poly ; 
      
      // Polytropic unit of mass in terms of m_unit :
      double m_poly = r_poly / ggrav ; 
      
      // Polytropic unit of angular momentum in terms of j_unit :
      double j_poly = r_poly * r_poly / ggrav ; 
      
      ost << "#      d [poly]       "  
	  << "       d_G [poly]     "
	  << "     Omega [poly]     "
	  << "     M_ADM [poly]     "     
	  << "       J [poly]       "  
	  << "    M_B(1) [poly]     "
	  << "    M_B(2) [poly]     " << endl ; 
      
      ost.width(20) ; 
      ost << separation() / r_poly ; ost.width(22) ;
      ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly ; ost.width(22) ; 
      ost << omega * t_poly ; ost.width(22) ;
      ost << mass_adm() / m_poly ; ost.width(22) ;
      ost << angu_mom()(2) / j_poly ; ost.width(22) ;
      ost << star1.mass_b() / m_poly ; ost.width(22) ;
      ost << star2.mass_b() / m_poly << endl ; 
      
  }
  
}


   
		    //-------------------------------//
		    //		Miscellaneous	     //
		    //-------------------------------//

double Binary::separation() const {
    
    double dx = star1.mp.get_ori_x() - star2.mp.get_ori_x() ; 
    double dy = star1.mp.get_ori_y() - star2.mp.get_ori_y() ; 
    double dz = star1.mp.get_ori_z() - star2.mp.get_ori_z() ; 
    
    return sqrt( dx*dx + dy*dy + dz*dz ) ; 
    
}
}