LORENE
binaire.h
1 /*
2  * Definition of Lorene class Binaire
3  *
4  */
5 
6 /*
7  * Copyright (c) 2000-2003 Eric Gourgoulhon
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 #ifndef __BINAIRE_H_
29 #define __BINAIRE_H_
30 
31 /*
32  * $Id: binaire.h,v 1.7 2017/02/24 15:34:59 j_novak Exp $
33  * $Log: binaire.h,v $
34  * Revision 1.7 2017/02/24 15:34:59 j_novak
35  * Removal of spurious comments
36  *
37  * Revision 1.6 2014/10/13 08:52:32 j_novak
38  * Lorene classes and functions now belong to the namespace Lorene.
39  *
40  * Revision 1.5 2009/06/18 18:42:13 k_taniguchi
41  * Defined a slightly modified code to determine
42  * the orbital angular velocity.
43  *
44  * Revision 1.4 2003/09/15 15:09:47 e_gourgoulhon
45  * Added the member function write_global.
46  *
47  * Revision 1.3 2002/06/17 14:05:16 j_novak
48  * friend functions are now also declared outside the class definition
49  *
50  * Revision 1.2 2001/12/20 13:00:31 k_taniguchi
51  * Addition of the Komar mass, the virial error by Gourgoulhon and Bonazzola, and the virial error by Friedman, Uryu, and Shibata.
52  *
53  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
54  * LORENE
55  *
56  * Revision 2.11 2001/02/28 14:27:32 keisuke
57  * *** empty log message ***
58  *
59  * Revision 2.10 2000/07/07 14:10:04 eric
60  * AJout de display_poly.
61  *
62  * Revision 2.9 2000/03/17 15:26:53 eric
63  * Ajout de la fonction analytical_omega.
64  *
65  * Revision 2.8 2000/03/15 16:43:12 eric
66  * Ajout de la fonction analytical_shift().
67  *
68  * Revision 2.7 2000/03/13 14:24:50 eric
69  * Ajout des membres p_ham_constr et p_mom_constr ainsi que des
70  * fonctions de calcul associees (verification des equations de contrainte).
71  *
72  * Revision 2.6 2000/02/18 14:52:05 eric
73  * Ajout des membres p_virial et p_total_ener et des fonctions de calcul
74  * associees.
75  *
76  * Revision 2.5 2000/02/12 17:36:39 eric
77  * Ajout de la fonction separation().
78  *
79  * Revision 2.4 2000/02/12 17:09:02 eric
80  * Ajout de la fonction orbit.
81  *
82  * Revision 2.3 2000/02/04 17:15:05 eric
83  * Ajout du membre ref_triad.
84  *
85  * Revision 2.2 2000/02/02 10:13:08 eric
86  * Ajout des fonctions de lecture/ecriture omega, x_axe.
87  *
88  * Revision 2.1 2000/01/31 17:01:46 eric
89  * *** empty log message ***
90  *
91  * Revision 2.0 2000/01/31 15:57:39 eric
92  * *** empty log message ***
93  *
94  *
95  * $Header: /cvsroot/Lorene/C++/Include/binaire.h,v 1.7 2017/02/24 15:34:59 j_novak Exp $
96  *
97  */
98 
99 #include "etoile.h"
100 
101 namespace Lorene {
107 class Binaire {
108 
109  // Data :
110  // -----
111  private:
116 
119 
122 
127  Etoile_bin* et[2] ;
128 
132  double omega ;
133 
136  double x_axe ;
137 
138  // Derived data :
139  // ------------
140  private:
142  mutable double* p_mass_adm ;
143 
145  mutable double* p_mass_kom ;
146 
148  mutable Tbl* p_angu_mom ;
149 
151  mutable double* p_total_ener ;
152 
154  mutable double* p_virial ;
155 
157  mutable double* p_virial_gb ;
158 
160  mutable double* p_virial_fus ;
161 
163  mutable double* p_ham_constr ;
164 
166  mutable Tbl* p_mom_constr ;
167 
168 
169  // Constructors - Destructor
170  // -------------------------
171  public:
188  Binaire(Map& mp1, int nzet1, const Eos& eos1, int irrot1,
189  Map& mp2, int nzet2, const Eos& eos2, int irrot2,
190  int relat) ;
191 
192 
193  Binaire(const Binaire& ) ;
194 
204  Binaire(Map& mp1, const Eos& eos1, Map& mp2, const Eos& eos2,
205  FILE* fich) ;
206 
207  ~Binaire() ;
208 
209 
210  // Memory management
211  // -----------------
212  private:
214  void del_deriv() const ;
215 
217  void set_der_0x0() const ;
218 
219  // Mutators / assignment
220  // ---------------------
221  public:
223  void operator=(const Binaire&) ;
224 
226  Etoile_bin& set(int i)
227  { assert( (i==1) || (i==2) );
228  del_deriv() ;
229  return *et[i-1] ;} ;
230 
232  double& set_omega() {return omega; } ;
233 
235  double& set_x_axe() {return x_axe; } ;
236 
237 
238  // Accessors
239  // ---------
240  public:
242  const Etoile_bin& operator()(int i) const
243  { assert( (i==1) || (i==2) );
244  return *et[i-1] ;} ;
245 
247  double get_omega() const {return omega; } ;
248 
250  double get_x_axe() const {return x_axe; } ;
251 
255  double separation() const ;
256 
257  // Outputs
258  // -------
259  public:
260  void sauve(FILE *) const ;
261 
262  // Display
263  friend ostream& operator<<(ostream& , const Binaire& ) ;
264 
266  void display_poly(ostream& ) const ;
267 
271  void write_global(ostream& ) const ;
272 
273  private:
275  ostream& operator>>(ostream& ) const ;
276 
277 
278  // Computational routines
279  // ----------------------
280  public:
282  double mass_adm() const ;
283 
285  double mass_kom() const ;
286 
294  const Tbl& angu_mom() const ;
295 
304  double total_ener() const ;
305 
310  double virial() const ;
311 
317  double virial_gb() const ;
318 
326  double virial_fus() const ;
327 
335  double ham_constr() const ;
336 
343  const Tbl& mom_constr() const ;
344 
371  void orbit(double fact_omeg_min, double fact_omeg_max, double& xgg1,
372  double& xgg2) ;
373 
404  void orbit_eqmass(double fact_omeg_min, double fact_omeg_max,
405  double mass1, double mass2,
406  double& xgg1, double& xgg2) ;
407 
411  void analytical_omega() ;
412 
417  void analytical_shift() ;
418 
419 };
420 ostream& operator<<(ostream& , const Binaire& ) ;
421 
422 }
423 #endif
double total_ener() const
Total energy (excluding the rest mass energy).
const Etoile_bin & operator()(int i) const
Returns a reference to the star no. i.
Definition: binaire.h:242
Etoile_bin star1
First star of the system.
Definition: binaire.h:118
double mass_adm() const
Total ADM mass.
Lorene prototypes.
Definition: app_hor.h:67
Tbl * p_mom_constr
Relative error on the momentum constraint.
Definition: binaire.h:166
void display_poly(ostream &) const
Display in polytropic units.
Definition: binaire.C:287
Equation of state base class.
Definition: eos.h:209
double * p_ham_constr
Relative error on the Hamiltonian constraint.
Definition: binaire.h:163
double * p_virial
Virial theorem error.
Definition: binaire.h:154
Base class for coordinate mappings.
Definition: map.h:688
const Base_vect_cart ref_triad
Cartesian triad of the absolute reference frame.
Definition: binaire.h:115
double virial_gb() const
Estimates the relative error on the virial theorem calculated by E.Gourgoulhon and S...
ostream & operator>>(ostream &) const
Operator >> (function called by the operator <<).
Definition: binaire.C:260
Class for stars in binary system.
Definition: etoile.h:817
double x_axe
Absolute X coordinate of the rotation axis.
Definition: binaire.h:136
void orbit_eqmass(double fact_omeg_min, double fact_omeg_max, double mass1, double mass2, double &xgg1, double &xgg2)
Computes the orbital angular velocity { omega} and the position of the rotation axis { x_axe}...
double * p_total_ener
Total energy of the system.
Definition: binaire.h:151
double & set_omega()
Sets the orbital angular velocity [{ f_unit}].
Definition: binaire.h:232
void analytical_shift()
Sets some analytical template for the shift vector (via the members { w_shift} and { khi_shift} of th...
const Tbl & mom_constr() const
Estimates the relative error on the momentum constraint equation by comparing ${}_j K^{ij}$ with {equ...
Etoile_bin * et[2]
Array of the two stars (to perform loops on the stars): { et[0]} contains the address of { star1} and...
Definition: binaire.h:127
double ham_constr() const
Estimates the relative error on the Hamiltonian constraint equation by comparing $ A$ with {equation}...
void operator=(const Binaire &)
Assignment to another { Binaire}.
Definition: binaire.C:220
double omega
Angular velocity with respect to an asymptotically inertial observer.
Definition: binaire.h:132
double separation() const
Returns the coordinate separation of the two stellar centers [{ r_unit}].
Definition: binaire.C:460
double * p_virial_fus
Virial theorem error by J.L.Friedman, K.Uryu, and M.Shibata.
Definition: binaire.h:160
double get_x_axe() const
Returns the absolute coordinate X of the rotation axis [{ r_unit}].
Definition: binaire.h:250
const Tbl & angu_mom() const
Total angular momentum.
void write_global(ostream &) const
Write global quantities in a formatted file.
Definition: binaire.C:337
Tbl * p_angu_mom
Total angular momentum of the system.
Definition: binaire.h:148
Binary systems.
Definition: binaire.h:107
Cartesian vectorial bases (triads).
Definition: base_vect.h:201
void orbit(double fact_omeg_min, double fact_omeg_max, double &xgg1, double &xgg2)
Computes the orbital angular velocity { omega} and the position of the rotation axis { x_axe}...
void analytical_omega()
Sets the orbital angular velocity to some 2-PN analytical value (Keplerian value in the Newtonian cas...
double get_omega() const
Returns the orbital angular velocity [{ f_unit}].
Definition: binaire.h:247
double mass_kom() const
Total Komar mass.
double virial_fus() const
Estimates the relative error on the virial theorem calculated by J.L.Friedman, K.Uryu, and M.Shibata (PRD accepted, gr-qc/0108070)
double & set_x_axe()
Sets the absolute coordinate X of the rotation axis [{ r_unit}].
Definition: binaire.h:235
friend ostream & operator<<(ostream &, const Binaire &)
Save in a file.
Definition: binaire.C:254
double * p_mass_adm
Total ADM mass of the system.
Definition: binaire.h:142
double * p_mass_kom
Total Komar mass of the system.
Definition: binaire.h:145
double * p_virial_gb
Virial theorem error by E.Gourgoulhon and S.Bonazzola.
Definition: binaire.h:157
void del_deriv() const
Destructor.
Definition: binaire.C:180
double virial() const
Estimates the relative error on the virial theorem (for a relativistic one, it returns $|1 - M_{ Koma...
void set_der_0x0() const
Sets to { 0x0} all the pointers on derived quantities.
Definition: binaire.C:198
Basic array class.
Definition: tbl.h:164
Binaire(Map &mp1, int nzet1, const Eos &eos1, int irrot1, Map &mp2, int nzet2, const Eos &eos2, int irrot2, int relat)
Standard constructor.
Definition: binaire.C:112
Etoile_bin star2
Second star of the system.
Definition: binaire.h:121