LORENE
binary.C
1/*
2 * Methods of class Binary
3 *
4 */
5
6/*
7 * Copyright (c) 2004 Francois Limousin
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 as published by
13 * the Free Software Foundation; either version 2 of the License, or
14 * (at your option) any later version.
15 *
16 * LORENE is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 * GNU General Public License for more details.
20 *
21 * You should have received a copy of the GNU General Public License
22 * along with LORENE; if not, write to the Free Software
23 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
24 *
25 */
26
27
28/*
29 * $Id: binary.C,v 1.18 2018/11/16 14:34:35 j_novak Exp $
30 * $Log: binary.C,v $
31 * Revision 1.18 2018/11/16 14:34:35 j_novak
32 * Changed minor points to avoid some compilation warnings.
33 *
34 * Revision 1.17 2016/12/05 16:17:47 j_novak
35 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36 *
37 * Revision 1.16 2014/10/13 08:52:44 j_novak
38 * Lorene classes and functions now belong to the namespace Lorene.
39 *
40 * Revision 1.15 2014/10/06 15:12:59 j_novak
41 * Modified #include directives to use c++ syntax.
42 *
43 * Revision 1.14 2005/09/15 14:39:14 e_gourgoulhon
44 * Added printing of angular momentum in display_poly.
45 *
46 * Revision 1.13 2005/09/13 19:38:31 f_limousin
47 * Reintroduction of the resolution of the equations in cartesian coordinates.
48 *
49 * Revision 1.12 2005/02/24 17:31:27 f_limousin
50 * Update of the function decouple().
51 *
52 * Revision 1.11 2005/02/18 13:14:06 j_novak
53 * Changing of malloc/free to new/delete + suppression of some unused variables
54 * (trying to avoid compilation warnings).
55 *
56 * Revision 1.10 2004/07/21 11:47:10 f_limousin
57 * Add / Delete comments.
58 *
59 * Revision 1.9 2004/05/25 14:25:12 f_limousin
60 * Add the virial theorem for conformally flat configurations.
61 *
62 * Revision 1.8 2004/03/25 10:29:01 j_novak
63 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
64 *
65 * Revision 1.7 2004/03/23 10:00:47 f_limousin
66 * Minor changes.
67 *
68 * Revision 1.6 2004/02/27 09:59:33 f_limousin
69 * Modification in the routine decouple().
70 *
71 * Revision 1.5 2004/02/21 17:05:12 e_gourgoulhon
72 * Method Scalar::point renamed Scalar::val_grid_point.
73 * Method Scalar::set_point renamed Scalar::set_grid_point.
74 *
75 * Revision 1.4 2004/01/20 15:21:12 f_limousin
76 * First version
77 *
78 *
79 * $Header: /cvsroot/Lorene/C++/Source/Binary/binary.C,v 1.18 2018/11/16 14:34:35 j_novak Exp $
80 *
81 */
82
83// Headers C
84#include <cmath>
85
86// Headers Lorene
87#include "binary.h"
88#include "eos.h"
89#include "utilitaires.h"
90#include "graphique.h"
91#include "param.h"
92#include "unites.h"
93
94 //--------------//
95 // Constructors //
96 //--------------//
97
98// Standard constructor
99// --------------------
100
101namespace Lorene {
102Binary::Binary(Map& mp1, int nzet1, const Eos& eos1, int irrot1,
103 Map& mp2, int nzet2, const Eos& eos2, int irrot2,
104 int conf_flat)
105 : star1(mp1, nzet1, eos1, irrot1, conf_flat),
106 star2(mp2, nzet2, eos2, irrot2, conf_flat)
107{
108
109 et[0] = &star1 ;
110 et[1] = &star2 ;
111
112 omega = 0 ;
113 x_axe = 0 ;
114
115 // Pointers of derived quantities initialized to zero :
116 set_der_0x0() ;
117}
118
119// Copy constructor
120// ----------------
122 : star1(bibi.star1),
123 star2(bibi.star2),
124 omega(bibi.omega),
125 x_axe(bibi.x_axe)
126{
127 et[0] = &star1 ;
128 et[1] = &star2 ;
129
130 // Pointers of derived quantities initialized to zero :
131 set_der_0x0() ;
132}
133
134// Constructor from a file
135// -----------------------
136Binary::Binary(Map& mp1, const Eos& eos1, Map& mp2, const Eos& eos2,
137 FILE* fich)
138 : star1(mp1, eos1, fich),
139 star2(mp2, eos2, fich)
140{
141 et[0] = &star1 ;
142 et[1] = &star2 ;
143
144 // omega and x_axe are read in the file:
145 fread_be(&omega, sizeof(double), 1, fich) ;
146 fread_be(&x_axe, sizeof(double), 1, fich) ;
147
148 // Pointers of derived quantities initialized to zero :
149 set_der_0x0() ;
150
151}
152
153 //------------//
154 // Destructor //
155 //------------//
156
158
159 del_deriv() ;
160
161}
162
163 //----------------------------------//
164 // Management of derived quantities //
165 //----------------------------------//
166
167void Binary::del_deriv() const {
168
169 if (p_mass_adm != 0x0) delete p_mass_adm ;
170 if (p_mass_kom != 0x0) delete p_mass_kom ;
171 if (p_angu_mom != 0x0) delete p_angu_mom ;
172 if (p_total_ener != 0x0) delete p_total_ener ;
173 if (p_virial != 0x0) delete p_virial ;
174 if (p_ham_constr != 0x0) delete p_ham_constr ;
175 if (p_mom_constr != 0x0) delete p_mom_constr ;
176
177 set_der_0x0() ;
178}
179
180
181
182
184
185 p_mass_adm = 0x0 ;
186 p_mass_kom = 0x0 ;
187 p_angu_mom = 0x0 ;
188 p_total_ener = 0x0 ;
189 p_virial = 0x0 ;
190 p_ham_constr = 0x0 ;
191 p_mom_constr = 0x0 ;
192
193}
194
195
196 //--------------//
197 // Assignment //
198 //--------------//
199
200// Assignment to another Binary
201// --------------------------------
202
203void Binary::operator=(const Binary& bibi) {
204
205 star1 = bibi.star1 ;
206 star2 = bibi.star2 ;
207
208 omega = bibi.omega ;
209 x_axe = bibi.x_axe ;
210
211 del_deriv() ; // Deletes all derived quantities
212
213}
214
215 //--------------//
216 // Outputs //
217 //--------------//
218
219// Save in a file
220// --------------
221void Binary::sauve(FILE* fich) const {
222
223 star1.sauve(fich) ;
224 star2.sauve(fich) ;
225
226 fwrite_be(&omega, sizeof(double), 1, fich) ;
227 fwrite_be(&x_axe, sizeof(double), 1, fich) ;
228
229}
230
231// Printing
232// --------
233ostream& operator<<(ostream& ost, const Binary& bibi) {
234 bibi >> ost ;
235 return ost ;
236}
237
238
239ostream& Binary::operator>>(ostream& ost) const {
240
241 using namespace Unites ;
242
243 ost << endl ;
244 ost << "Binary neutron stars" << endl ;
245 ost << "=============" << endl ;
246 ost << endl <<
247 "Orbital angular velocity : " << omega * f_unit << " rad/s" << endl ;
248 ost << endl <<
249 "Coordinate separation between the two stellar centers : "
250 << separation() / km << " km" << endl ;
251 ost <<
252 "Absolute coordinate X of the rotation axis : " << x_axe / km
253 << " km" << endl ;
254 ost << endl << "Star 1 : " << endl ;
255 ost << "====== " << endl ;
256 ost << star1 << endl ;
257 ost << "Star 2 : " << endl ;
258 ost << "====== " << endl ;
259 ost << star2 << endl ;
260 return ost ;
261}
262
263// Display in polytropic units
264// ---------------------------
265
266void Binary::display_poly(ostream& ost) const {
267
268 using namespace Unites ;
269
270 const Eos* p_eos1 = &( star1.get_eos() ) ;
271 const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;
272
273 if (p_eos_poly != 0x0) {
274
275 assert( star1.get_eos() == star2.get_eos() ) ;
276
277 double kappa = p_eos_poly->get_kap() ;
278 double gamma = p_eos_poly->get_gam() ; ;
279 double kap_ns2 = pow( kappa, 0.5 /(gamma-1) ) ;
280
281 // Polytropic unit of length in terms of r_unit :
282 double r_poly = kap_ns2 / sqrt(ggrav) ;
283
284 // Polytropic unit of time in terms of t_unit :
285 double t_poly = r_poly ;
286
287 // Polytropic unit of mass in terms of m_unit :
288 double m_poly = r_poly / ggrav ;
289
290 // Polytropic unit of angular momentum in terms of j_unit :
291 double j_poly = r_poly * r_poly / ggrav ;
292
293 ost.precision(10) ;
294 ost << endl << "Quantities in polytropic units : " << endl ;
295 ost << "==============================" << endl ;
296 ost << " ( r_poly = " << r_poly / km << " km )" << endl ;
297 ost << " d_e_max : " << separation() / r_poly << endl ;
298 ost << " d_G : "
299 << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly
300 << endl ;
301 ost << " Omega : " << omega * t_poly << endl ;
302 ost << " J : " << angu_mom()(2) / j_poly << endl ;
303 ost << " M_ADM : " << mass_adm() / m_poly << endl ;
304 ost << " M_Komar : " << mass_kom() / m_poly << endl ;
305 ost << " M_bar(star 1) : " << star1.mass_b() / m_poly << endl ;
306 ost << " M_bar(star 2) : " << star2.mass_b() / m_poly << endl ;
307 ost << " R_0(star 1) : " <<
308 0.5 * ( star1.ray_eq() + star1.ray_eq_pi() ) / r_poly << endl ;
309 ost << " R_0(star 2) : " <<
310 0.5 * ( star2.ray_eq() + star2.ray_eq_pi() ) / r_poly << endl ;
311
312 }
313
314
315}
316
317
319
320 int nz_un = star1.get_mp().get_mg()->get_nzone() ;
321 int nz_deux = star2.get_mp().get_mg()->get_nzone() ;
322
323 // We determine R_limite :
324 double distance = star2.get_mp().get_ori_x() - star1.get_mp().get_ori_x() ;
325 cout << "distance = " << distance << endl ;
326 double lim_un = distance/2. ;
327 double lim_deux = distance/2. ;
328 double int_un = distance/6. ;
329 double int_deux = distance/6. ;
330
331 // The functions used.
332 Scalar fonction_f_un (star1.get_mp()) ;
333 fonction_f_un = 0.5*pow(
334 cos((star1.get_mp().r-int_un)*M_PI/2./(lim_un-int_un)), 2.)+0.5 ;
335 fonction_f_un.std_spectral_base();
336
337 Scalar fonction_g_un (star1.get_mp()) ;
338 fonction_g_un = 0.5*pow
339 (sin((star1.get_mp().r-int_un)*M_PI/2./(lim_un-int_un)), 2.) ;
340 fonction_g_un.std_spectral_base();
341
342 Scalar fonction_f_deux (star2.get_mp()) ;
343 fonction_f_deux = 0.5*pow(
344 cos((star2.get_mp().r-int_deux)*M_PI/2./(lim_deux-int_deux)), 2.)+0.5 ;
345 fonction_f_deux.std_spectral_base();
346
347 Scalar fonction_g_deux (star2.get_mp()) ;
348 fonction_g_deux = 0.5*pow(
349 sin((star2.get_mp().r-int_deux)*M_PI/2./(lim_un-int_deux)), 2.) ;
350 fonction_g_deux.std_spectral_base();
351
352 // The functions total :
353 Scalar decouple_un (star1.get_mp()) ;
354 decouple_un.allocate_all() ;
355 Scalar decouple_deux (star2.get_mp()) ;
356 decouple_deux.allocate_all() ;
357
358 Mtbl xabs_un (star1.get_mp().xa) ;
359 Mtbl yabs_un (star1.get_mp().ya) ;
360 Mtbl zabs_un (star1.get_mp().za) ;
361
362 Mtbl xabs_deux (star2.get_mp().xa) ;
363 Mtbl yabs_deux (star2.get_mp().ya) ;
364 Mtbl zabs_deux (star2.get_mp().za) ;
365
366 double xabs, yabs, zabs, air_un, air_deux, theta, phi ;
367
368 for (int l=0 ; l<nz_un ; l++) {
369 int nr = star1.get_mp().get_mg()->get_nr (l) ;
370
371 if (l==nz_un-1)
372 nr -- ;
373
374 int np = star1.get_mp().get_mg()->get_np (l) ;
375 int nt = star1.get_mp().get_mg()->get_nt (l) ;
376
377 for (int k=0 ; k<np ; k++)
378 for (int j=0 ; j<nt ; j++)
379 for (int i=0 ; i<nr ; i++) {
380
381 xabs = xabs_un (l, k, j, i) ;
382 yabs = yabs_un (l, k, j, i) ;
383 zabs = zabs_un (l, k, j, i) ;
384
385 // Coordinates of the point
386 star1.get_mp().convert_absolute
387 (xabs, yabs, zabs, air_un, theta, phi) ;
388 star2.get_mp().convert_absolute
389 (xabs, yabs, zabs, air_deux, theta, phi) ;
390
391 if (air_un <= lim_un)
392 if (air_un < int_un)
393 decouple_un.set_grid_point(l, k, j, i) = 1 ;
394 else
395 // Close to star 1 :
396 decouple_un.set_grid_point(l, k, j, i) =
397 fonction_f_un.val_grid_point(l, k, j, i) ;
398 else
399 if (air_deux <= lim_deux)
400 if (air_deux < int_deux)
401 decouple_un.set_grid_point(l, k, j, i) = 0 ;
402 else
403 // Close to star 2 :
404 decouple_un.set_grid_point(l, k, j, i) =
405 fonction_g_deux.val_point (air_deux, theta, phi) ;
406
407 else
408 // Far from each star :
409 decouple_un.set_grid_point(l, k, j, i) = 0.5 ;
410 }
411
412 // Case infinity :
413 if (l==nz_un-1)
414 for (int k=0 ; k<np ; k++)
415 for (int j=0 ; j<nt ; j++)
416 decouple_un.set_grid_point(nz_un-1, k, j, nr)=0.5 ;
417 }
418
419 for (int l=0 ; l<nz_deux ; l++) {
420 int nr = star2.get_mp().get_mg()->get_nr (l) ;
421
422 if (l==nz_deux-1)
423 nr -- ;
424
425 int np = star2.get_mp().get_mg()->get_np (l) ;
426 int nt = star2.get_mp().get_mg()->get_nt (l) ;
427
428 for (int k=0 ; k<np ; k++)
429 for (int j=0 ; j<nt ; j++)
430 for (int i=0 ; i<nr ; i++) {
431
432 xabs = xabs_deux (l, k, j, i) ;
433 yabs = yabs_deux (l, k, j, i) ;
434 zabs = zabs_deux (l, k, j, i) ;
435
436 // coordinates of the point :
437 star1.get_mp().convert_absolute
438 (xabs, yabs, zabs, air_un, theta, phi) ;
439 star2.get_mp().convert_absolute
440 (xabs, yabs, zabs, air_deux, theta, phi) ;
441
442 if (air_deux <= lim_deux)
443 if (air_deux < int_deux)
444 decouple_deux.set_grid_point(l, k, j, i) = 1 ;
445 else
446 // close to star two :
447 decouple_deux.set_grid_point(l, k, j, i) =
448 fonction_f_deux.val_grid_point(l, k, j, i) ;
449 else
450 if (air_un <= lim_un)
451 if (air_un < int_un)
452 decouple_deux.set_grid_point(l, k, j, i) = 0 ;
453 else
454 // close to star one :
455 decouple_deux.set_grid_point(l, k, j, i) =
456 fonction_g_un.val_point (air_un, theta, phi) ;
457
458 else
459 // Far from each star :
460 decouple_deux.set_grid_point(l, k, j, i) = 0.5 ;
461 }
462
463 // Case infinity :
464 if (l==nz_deux-1)
465 for (int k=0 ; k<np ; k++)
466 for (int j=0 ; j<nt ; j++)
467 decouple_deux.set_grid_point(nz_un-1, k, j, nr)=0.5 ;
468 }
469
470 decouple_un.std_spectral_base() ;
471 decouple_deux.std_spectral_base() ;
472
473 int nr = star1.get_mp().get_mg()->get_nr (1) ;
474 int nt = star1.get_mp().get_mg()->get_nt (1) ;
475 int np = star1.get_mp().get_mg()->get_np (1) ;
476 cout << "decouple_un" << endl << norme(decouple_un/(nr*nt*np)) << endl ;
477 cout << "decouple_deux" << endl << norme(decouple_deux/(nr*nt*np))
478 << endl ;
479
480 star1.decouple = decouple_un ;
481 star2.decouple = decouple_deux ;
482}
483
484void Binary::write_global(ostream& ost) const {
485
486 using namespace Unites ;
487
488 const Map& mp1 = star1.get_mp() ;
489 const Mg3d* mg1 = mp1.get_mg() ;
490 int nz1 = mg1->get_nzone() ;
491
492 ost.precision(5) ;
493 ost << "# Grid 1 : " << nz1 << "x"
494 << mg1->get_nr(0) << "x" << mg1->get_nt(0) << "x" << mg1->get_np(0)
495 << " R_out(l) [km] : " ;
496 for (int l=0; l<nz1; l++) {
497 ost << " " << mp1.val_r(l, 1., M_PI/2, 0) / km ;
498 }
499 ost << endl ;
500
501 ost << "# VE(M) " << endl ;
502
503
504 ost.setf(ios::scientific) ;
505 ost.width(14) ;
506 ost << virial() << endl ;
507
508 ost << "# d [km] "
509 << " d_G [km] "
510 << " d/(a1 +a1') "
511 << " f [Hz] "
512 << " M_ADM [M_sol] "
513 << " M_ADM_vol [M_sol] "
514 << " M_Komar [M_sol] "
515 << " M_Komar_vol [M_sol] "
516 << " J [G M_sol^2/c] " << endl ;
517
518 ost.precision(14) ;
519 ost.width(20) ;
520 ost << separation() / km ; ost.width(22) ;
521 ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / km ; ost.width(22) ;
522 ost << separation() / (star1.ray_eq() + star2.ray_eq()) ; ost.width(22) ;
523 ost << omega / (2*M_PI)* f_unit ; ost.width(22) ;
524 ost << mass_adm() / msol ; ost.width(22) ;
525 ost << mass_adm_vol() / msol ; ost.width(22) ;
526 ost << mass_kom() / msol ; ost.width(22) ;
527 ost << mass_kom_vol() / msol ; ost.width(22) ;
528 ost << angu_mom()(2)/ ( qpig / (4* M_PI) * msol*msol) << endl ;
529
530 ost << "# H_c(1)[c^2] "
531 << " e_c(1)[rho_nuc] "
532 << " M_B(1) [M_sol] "
533 << " r_eq(1) [km] "
534 << " a2/a1(1) "
535 << " a3/a1(1) " << endl ;
536
537 ost.width(20) ;
538 ost << star1.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
539 ost << star1.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
540 ost << star1.mass_b() / msol ; ost.width(22) ;
541 ost << star1.ray_eq() / km ; ost.width(22) ;
542 ost << star1.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
543 ost << star1.ray_pole() / star1.ray_eq() << endl ;
544
545 ost << "# H_c(2)[c^2] "
546 << " e_c(2)[rho_nuc] "
547 << " M_B(2) [M_sol] "
548 << " r_eq(2) [km] "
549 << " a2/a1(2) "
550 << " a3/a1(2) " << endl ;
551
552 ost.width(20) ;
553 ost << star2.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
554 ost << star2.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
555 ost << star2.mass_b() / msol ; ost.width(22) ;
556 ost << star2.ray_eq() / km ; ost.width(22) ;
557 ost << star2.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
558 ost << star2.ray_pole() / star1.ray_eq() << endl ;
559
560 // Quantities in polytropic units if the EOS is a polytropic one
561 // -------------------------------------------------------------
562 const Eos* p_eos1 = &( star1.get_eos() ) ;
563 const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;
564
565 if ((p_eos_poly != 0x0) && ( star1.get_eos() == star2.get_eos() )) {
566
567 double kappa = p_eos_poly->get_kap() ;
568 double gamma = p_eos_poly->get_gam() ; ;
569 double kap_ns2 = pow( kappa, 0.5 /(gamma-1.) ) ;
570
571 // Polytropic unit of length in terms of r_unit :
572 double r_poly = kap_ns2 / sqrt(ggrav) ;
573
574 // Polytropic unit of time in terms of t_unit :
575 double t_poly = r_poly ;
576
577 // Polytropic unit of mass in terms of m_unit :
578 double m_poly = r_poly / ggrav ;
579
580 // Polytropic unit of angular momentum in terms of j_unit :
581 double j_poly = r_poly * r_poly / ggrav ;
582
583 ost << "# d [poly] "
584 << " d_G [poly] "
585 << " Omega [poly] "
586 << " M_ADM [poly] "
587 << " J [poly] "
588 << " M_B(1) [poly] "
589 << " M_B(2) [poly] " << endl ;
590
591 ost.width(20) ;
592 ost << separation() / r_poly ; ost.width(22) ;
593 ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly ; ost.width(22) ;
594 ost << omega * t_poly ; ost.width(22) ;
595 ost << mass_adm() / m_poly ; ost.width(22) ;
596 ost << angu_mom()(2) / j_poly ; ost.width(22) ;
597 ost << star1.mass_b() / m_poly ; ost.width(22) ;
598 ost << star2.mass_b() / m_poly << endl ;
599
600 }
601
602}
603
604
605
606 //-------------------------------//
607 // Miscellaneous //
608 //-------------------------------//
609
610double Binary::separation() const {
611
612 double dx = star1.mp.get_ori_x() - star2.mp.get_ori_x() ;
613 double dy = star1.mp.get_ori_y() - star2.mp.get_ori_y() ;
614 double dz = star1.mp.get_ori_z() - star2.mp.get_ori_z() ;
615
616 return sqrt( dx*dx + dy*dy + dz*dz ) ;
617
618}
619}
double mass_adm_vol() const
Total ADM mass (computed by a volume integral).
void sauve(FILE *) const
Save in a file.
Definition binary.C:221
void display_poly(ostream &) const
Display in polytropic units.
Definition binary.C:266
double mass_kom_vol() const
Total Komar mass (computed by a volume integral).
void write_global(ostream &) const
Write global quantities in a formatted file.
Definition binary.C:484
void del_deriv() const
Deletes all the derived quantities.
Definition binary.C:167
Star_bin * 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.h:89
double omega
Angular velocity with respect to an asymptotically inertial observer.
Definition binary.h:94
double * p_total_ener
Total energy of the system.
Definition binary.h:114
Tbl * p_angu_mom
Total angular momentum of the system.
Definition binary.h:111
Star_bin star2
Second star of the system.
Definition binary.h:83
const Tbl & angu_mom() const
Total angular momentum.
friend ostream & operator<<(ostream &, const Binary &)
Display.
Definition binary.C:233
void fait_decouple()
Calculates decouple which is used to obtain qq_auto by the formula : qq_auto = decouple * qq.
Definition binary.C:318
Star_bin star1
First star of the system.
Definition binary.h:80
double * p_mass_adm
Total ADM mass of the system.
Definition binary.h:105
double * p_virial
Virial theorem error.
Definition binary.h:117
double mass_adm() const
Total ADM mass.
~Binary()
Destructor.
Definition binary.C:157
void operator=(const Binary &)
Assignment to another Binary.
Definition binary.C:203
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition binary.C:183
double virial() const
Estimates the relative error on the virial theorem.
double x_axe
Absolute X coordinate of the rotation axis.
Definition binary.h:98
Tbl * p_mom_constr
Relative error on the momentum constraint.
Definition binary.h:123
Binary(Map &mp1, int nzet1, const Eos &eos1, int irrot1, Map &mp2, int nzet2, const Eos &eos2, int irrot2, int conf_flat)
Standard constructor.
Definition binary.C:102
double * p_mass_kom
Total Komar mass of the system.
Definition binary.h:108
double mass_kom() const
Total Komar mass.
double separation() const
Returns the coordinate separation of the two stellar centers [r_unit].
Definition binary.C:610
double * p_ham_constr
Relative error on the Hamiltonian constraint.
Definition binary.h:120
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<).
Definition binary.C:239
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
Multi-domain array.
Definition mtbl.h:118
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
Definition scalar.C:790
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition scalar.h:643
virtual void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition scalar.C:371
double & set_grid_point(int l, int k, int j, int i)
Setting the value of the field at a given grid point.
Definition scalar.h:690
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
Definition scalar.C:896
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Cmp sin(const Cmp &)
Sine.
Definition cmp_math.C:72
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition cmp_math.C:484
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
Cmp cos(const Cmp &)
Cosine.
Definition cmp_math.C:97
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
Coord phi
coordinate centered on the grid
Definition map.h:732
Standard units of space, time and mass.