LORENE
bhole.h
1 /*
2  * Definition of Lorene classes Bhole
3  * Bhole_binaire
4  *
5  */
6 
7 /*
8  * Copyright (c) 2000-2001 Philippe Grandclement
9  *
10  * This file is part of LORENE.
11  *
12  * LORENE is free software; you can redistribute it and/or modify
13  * it under the terms of the GNU General Public License as published by
14  * the Free Software Foundation; either version 2 of the License, or
15  * (at your option) any later version.
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 #ifndef __BHOLE_H_
30 #define __BHOLE_H_
31 
32 /*
33  * $Id: bhole.h,v 1.18 2014/10/13 08:52:32 j_novak Exp $
34  * $Log: bhole.h,v $
35  * Revision 1.18 2014/10/13 08:52:32 j_novak
36  * Lorene classes and functions now belong to the namespace Lorene.
37  *
38  * Revision 1.17 2007/05/15 12:44:18 p_grandclement
39  * Scalar version of reevaluate
40  *
41  * Revision 1.16 2007/04/24 20:15:30 f_limousin
42  * Implementation of Dirichlet and Neumann BC for the lapse
43  *
44  * Revision 1.15 2006/04/27 09:12:29 p_grandclement
45  * First try at irrotational black holes
46  *
47  * Revision 1.14 2005/08/29 15:10:12 p_grandclement
48  * Addition of things needed :
49  * 1) For BBH with different masses
50  * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
51  * WORKING YET !!!
52  *
53  * Revision 1.13 2004/12/02 09:33:04 p_grandclement
54  * *** empty log message ***
55  *
56  * Revision 1.12 2004/03/25 12:35:34 j_novak
57  * now using namespace Unites
58  *
59  * Revision 1.11 2004/03/25 09:14:06 j_novak
60  * *** empty log message ***
61  *
62  * Revision 1.10 2004/03/25 08:53:33 j_novak
63  * Error in the doc corrected.
64  *
65  * Revision 1.9 2004/03/24 17:14:04 j_novak
66  * Translation of comments to doxygen
67  *
68  * Revision 1.8 2003/11/25 07:12:58 k_taniguchi
69  * Change the argument of update_metric from the class Etoile_bin to Et_bin_nsbh.
70  *
71  * Revision 1.7 2003/11/13 13:43:53 p_grandclement
72  * Addition of things needed for Bhole::update_metric (const Etoile_bin&, double, double)
73  *
74  * Revision 1.6 2003/10/24 13:05:48 p_grandclement
75  * correction of the equations for Bin_ns_bh...
76  *
77  * Revision 1.5 2003/02/13 16:40:24 p_grandclement
78  * Addition of various things for the Bin_ns_bh project, non of them being
79  * completely tested
80  *
81  * Revision 1.4 2003/01/31 16:57:12 p_grandclement
82  * addition of the member Cmp decouple used to compute the K_ij auto, once
83  * the K_ij total is known
84  *
85  * Revision 1.3 2002/12/18 10:28:49 e_gourgoulhon
86  *
87  * Added the set_mp function.
88  *
89  * Revision 1.2 2002/09/13 09:17:31 j_novak
90  * Modif. commentaires
91  *
92  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
93  * LORENE
94  *
95  * Revision 2.45 2001/06/28 15:08:16 eric
96  * Ajout de la fonction area().
97  *
98  * Revision 2.44 2001/05/16 11:32:45 phil
99  * ajout de init_bhole_seul
100  *
101  * Revision 2.43 2001/05/07 12:35:28 phil
102  * mmodifi commentaires
103  *
104  * Revision 2.42 2001/05/07 11:43:43 phil
105  * *** empty log message ***
106  *
107  * Revision 2.41 2001/05/07 11:40:50 phil
108  * mise a jour des commentaires
109  *
110  * Revision 2.40 2001/05/07 09:30:37 phil
111  * *** empty log message ***
112  *
113  * Revision 2.39 2001/05/07 09:28:37 phil
114  * *** empty log message ***
115  *
116  * Revision 2.38 2001/05/07 09:11:31 phil
117  * *** empty log message ***
118  *
119  * Revision 2.37 2001/04/27 09:25:11 phil
120  * fait_decouple est devenu public
121  *
122  * Revision 2.36 2001/04/26 12:04:11 phil
123  * *** empty log message ***
124  *
125  * Revision 2.35 2001/04/05 13:42:40 phil
126  * *** empty log message ***
127  *
128  * Revision 2.34 2001/04/05 13:34:29 phil
129  * ajout resolution en utilisant phi
130  *
131  * Revision 2.33 2001/03/22 10:49:34 phil
132  * *** empty log message ***
133  *
134  * Revision 2.32 2001/03/22 10:41:25 phil
135  * pleins de modif
136  *
137  * Revision 2.31 2001/02/28 13:22:31 phil
138  * modif vire kk_auto
139  *
140  * Revision 2.30 2001/02/12 15:37:42 phil
141  * ajout calcul de J a linfini
142  *
143  * Revision 2.29 2001/01/29 14:29:56 phil
144  * ajout type rotation
145  *
146  * Revision 2.28 2001/01/24 10:57:08 phil
147  * ajout de set_rayonm
148  *
149  * Revision 2.27 2001/01/22 09:29:08 phil
150  * vire les procedures qui servent pas
151  *
152  * Revision 2.26 2001/01/10 09:31:05 phil
153  * modification de fait_kk_auto (membre de binaire maintenant)
154  *
155  * Revision 2.25 2001/01/09 13:38:15 phil
156  * ajout memebre kk_auto
157  *
158  * Revision 2.24 2000/12/21 10:47:07 phil
159  * retour eventuel a 2.21
160  *
161  * Revision 2.21 2000/12/20 09:07:56 phil
162  * ajout set_statiques
163  *
164  * Revision 2.20 2000/12/18 16:40:39 phil
165  * *** empty log message ***
166  *
167  * Revision 2.19 2000/12/18 16:38:04 phil
168  * ajout convergence vers une masse donnee
169  *
170  * Revision 2.18 2000/12/15 16:41:28 phil
171  * ajout calcul de la separation
172  *
173  * Revision 2.17 2000/12/14 10:44:28 phil
174  * ATTENTION : PASSAGE DE PHI A PSI
175  *
176  * Revision 2.16 2000/12/13 15:35:18 phil
177  * ajout calcul bare_masse
178  *
179  * Revision 2.15 2000/12/04 14:29:37 phil
180  * ajout de grad_n_tot
181  *
182  * Revision 2.14 2000/12/01 16:12:52 phil
183  * *** empty log message ***
184  *
185  * Revision 2.13 2000/12/01 14:16:20 phil
186  * *** empty log message ***
187  *
188  * Revision 2.12 2000/11/24 15:14:44 phil
189  * ajout de find_horizon
190  *
191  * Revision 2.11 2000/11/24 09:57:18 phil
192  * ajout calcul masse et moment pour systeme
193  *
194  * Revision 2.10 2000/11/17 10:03:13 phil
195  * ajout de coal
196  *
197  * Revision 2.9 2000/11/15 18:26:29 phil
198  * simplifaction resolution du shift : on bosse a omega bloque
199  *
200  * Revision 2.8 2000/11/15 12:59:50 phil
201  * changement solve_shift_omega
202  *
203  * Revision 2.7 2000/11/15 09:41:03 phil
204  * *** empty log message ***
205  *
206  * Revision 2.6 2000/11/15 09:39:26 phil
207  * ajout viriel
208  *
209  * Revision 2.5 2000/11/03 12:56:41 phil
210  * ajout de const
211  *
212  * Revision 2.4 2000/10/26 08:20:45 phil
213  * ajout verifie_shift
214  *
215  * Revision 2.3 2000/10/23 09:15:12 phil
216  * modif commentaires
217  *
218  * Revision 2.2 2000/10/20 10:51:13 phil
219  * Modif commentaires (minimale !)
220  *
221  * Revision 2.1 2000/10/20 09:27:28 phil
222  * *** empty log message ***
223  *
224  * Revision 2.0 2000/10/20 09:22:04 phil
225  * *** empty log message ***
226  *
227  *
228  * $Header: /cvsroot/Lorene/C++/Include/bhole.h,v 1.18 2014/10/13 08:52:32 j_novak Exp $
229  *
230  */
231 
232 #include "tenseur.h"
233 
234 namespace Lorene {
235 
236 class Tenseur ;
237 class Tenseur_sym ;
238 class Map_af ;
239 class Bhole_binaire ;
240 class Et_bin_nsbh ;
241 
242 #define COROT 0
243 #define IRROT 1
244 
268 class Bhole {
269 
270  // Data :
271  protected:
272 
273  Map_af& mp ;
274  double rayon ;
275  double omega ;
276  double omega_local ;
277  int rot_state ;
278  double lapse_hori ;
279 
283  double* boost ;
284  double regul ;
285 
289 
293 
296 
298 
301 
306 
309 
319 
320  //Constructors :
321  public:
326  Bhole (Map_af& mapping) ;
327  Bhole (const Bhole&) ;
328  Bhole (Map_af&, FILE*, bool old=false) ;
329  ~Bhole() ;
330 
331  public:
332  //sauvegarde
333  void sauve (FILE* fich) const ;
334 
335  public:
336  void operator= (const Bhole&) ;
337 
339  const Map_af& get_mp() const {return mp;} ;
340 
342  Map_af& set_mp() {return mp; } ;
343 
347  double get_rayon() const {return rayon;} ;
348 
352  void set_rayon(double ray) {rayon = ray ;} ;
353 
357  double get_omega() const {return omega;} ;
361  void set_omega(double ome) {omega = ome ;} ;
365  double get_omega_local() const {return omega_local;} ;
369  void set_omega_local(double ome) {omega_local = ome ;} ;
373  double get_rot_state() const {return rot_state;} ;
377  void set_rot_state(int rotation) {rot_state = rotation;} ;
382  double* get_boost () const {return boost;} ;
383  void set_boost (double vx, double vy, double vz) {
384  boost[0] = vx ; boost[1] = vy ; boost[2] = vz ;
385  }
386 
390  double get_regul() const {return regul;} ;
391 
395  const Tenseur& get_n_auto() const {return n_auto ;} ;
399  Cmp return_n_auto() const {return n_auto() ;} ;
403  const Tenseur& get_n_comp() const {return n_comp ;} ;
407  const Tenseur& get_n_tot() const {return n_tot ;} ;
408 
412  const Tenseur& get_psi_auto() const {return psi_auto ;} ;
416  Cmp return_psi_auto() const {return psi_auto() ;} ;
420  const Tenseur& get_psi_comp() const {return psi_comp ;} ;
424  const Tenseur& get_psi_tot() const {return psi_tot ;} ;
425 
429  const Tenseur& get_grad_psi_tot() const {return grad_psi_tot ;} ;
430 
434  const Tenseur& get_grad_n_tot() const {return grad_n_tot ;} ;
435 
436 
440  const Tenseur& get_shift_auto() const {return shift_auto ;} ;
441 
442  Cmp return_shift_auto(int i) const {return shift_auto(i) ;} ;
443 
447  const Tenseur& get_taij_auto() const {return taij_auto ;} ;
451  const Tenseur& get_taij_comp() const {return taij_comp ;} ;
455  const Tenseur& get_taij_tot() const {return taij_tot ;}
456 
460  const Tenseur& get_tkij_tot() const {return tkij_tot ;}
464  const Tenseur& get_tkij_auto() const {return tkij_auto ;}
468  const Cmp get_decouple() const {return decouple ;}
469 
470  Cmp& set_n_auto() {return n_auto.set() ;}
471  Cmp& set_psi_auto() {return psi_auto.set() ;}
472 
473  void set_lapse_hori(double xx) {lapse_hori = xx;}
474  public:
475 
482  void fait_n_comp (const Bhole& comp) ;
483 
490  void fait_psi_comp (const Bhole& comp) ;
497  void fait_n_comp (const Et_bin_nsbh& comp) ;
498 
505  void fait_psi_comp (const Et_bin_nsbh& comp) ;
506 
510  void fait_taij_auto () ;
511 
512 
525  void regularise_shift (const Tenseur& comp) ;
526 
533  double viriel_seul () const ;
534 
544  void init_bhole () ;
545 
572  void init_kerr (double masse, double moment) ;
573 
585  void solve_lapse_seul (double relax) ;
586 
600  void solve_psi_seul (double relax) ;
601 
617  void solve_shift_seul (double precis, double relax) ;
618 
628  void regularise_seul () ;
641  void solve_lapse_with_ns (double relax, int bound_nn,
642  double lim_nn) ;
643 
657  void solve_psi_with_ns (double relax) ;
658 
675  void solve_shift_with_ns (const Et_bin_nsbh& ns,
676  double precis, double relax,
677  int bound_nn, double lim_nn) ;
678 
679  void equilibrium (const Et_bin_nsbh& ns, double precis, double relax,
680  int bound_nn, double lim_nn) ;
681  void update_metric (const Et_bin_nsbh& ns) ;
682 
683 
688  void fait_tkij() ;
689 
691  double area() const ;
692 
699  double masse_adm_seul () const ;
700 
707  double masse_komar_seul() const ;
708 
717  double moment_seul_inf() const ;
718 
727  double moment_seul_hor() const ;
728 
729  /*
730  * local angular momentum
731  */
732  double local_momentum() const ;
733 
737  void init_bhole_phi () ;
738 
745  void init_bhole_seul () ;
746 
747  friend class Bhole_binaire ;
748  friend class Bin_ns_bh ;
749 } ;
750 
758 
759  // data :
760  private:
761  // les deux trous noirs.
764 
765  // Tableau sur les deux trous.
766  Bhole* holes[2] ;
767 
768  double pos_axe ;
769  double omega ;
770 
771  public:
772  // constructeurs & destructeur
773  Bhole_binaire(Map_af& mp1, Map_af& mp2) ;
774  Bhole_binaire(const Bhole_binaire& ) ;
775  ~Bhole_binaire() ;
776 
777  public:
778  // trucs pour modifier
779  void operator=(const Bhole_binaire&) ;
780 
785  Bhole& set(int i)
786  { assert( (i==1) || (i==2) );
787  return *holes[i-1] ;} ;
791  void set_omega(double ome) {omega = ome ;
792  hole1.set_omega (ome) ;
793  hole2.set_omega (ome) ;} ;
794 
795  void set_pos_axe (double) ;
796 
797  public:
798  // trucs pour lire :
803  const Bhole& operator()(int i) const
804  { assert( (i==1) || (i==2) );
805  return *holes[i-1] ;} ;
806 
808  double get_omega() const {return omega; } ;
809 
817  void init_bhole_binaire() ;
818 
823  double viriel() const ;
837  void solve_lapse (double precis, double relax) ;
838 
854  void solve_psi (double precis, double relax) ;
855 
872  void solve_shift (double precis, double relax) ;
873 
883  void fait_tkij () ;
884 
891  void fait_decouple () ;
892 
893  public:
901  void set_statiques (double precis, double relax) ;
902 
913  void coal (double precis, double relax, int nbre_ome, double search_ome,
914  double m1, double m2, const int sortie = 0) ;
915 
916 
921  double adm_systeme() const ;
922 
927  double komar_systeme() const ;
928 
934  double moment_systeme_inf() ;
935 
943  double moment_systeme_hor() const ;
944 
950  double distance_propre(const int nr = 65) const ;
951 
952  Tbl linear_momentum_systeme_inf() const ;
953 
957  void solve_phi (double precision, double relax) ;
961  void init_phi() ;
962 } ;
963 
964 }
965 #endif
const Tenseur & get_grad_psi_tot() const
Returns the gradient of .
Definition: bhole.h:429
void fait_tkij()
Calculation af the extrinsic curvature tensor.
Definition: bhole_kij.C:91
void solve_phi(double precision, double relax)
Solve the equation for the logarithm of .
double moment_seul_hor() const
Calculates the angular momentum of the black hole using the formula on the horizon : where and H de...
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
const Tenseur & get_tkij_auto() const
Returns the part of generated by the hole.
Definition: bhole.h:464
Tenseur n_tot
Total N .
Definition: bhole.h:288
void init_kerr(double masse, double moment)
Set the inital values to those of Kerr.
Neutron star - black hole binary system.
Definition: bin_ns_bh.h:120
Tenseur grad_n_tot
Total gradient of N .
Definition: bhole.h:294
Bhole * holes[2]
Array on the black holes.
Definition: bhole.h:766
void solve_shift(double precis, double relax)
Solves the equation for the shift, using the Oohara-Nakarmure scheme~: using the current for the ...
Class intended to describe tensors with a symmetry on the two last indices *** DEPRECATED : use class...
Definition: tenseur.h:1256
~Bhole_binaire()
Destructor.
~Bhole()
Destructor.
Definition: bhole.C:216
Lorene prototypes.
Definition: app_hor.h:67
const Tenseur & get_psi_comp() const
Returns the part of generated by the companion hole.
Definition: bhole.h:420
void set_rayon(double ray)
Sets the radius of the horizon to ray .
Definition: bhole.h:352
Bhole(Map_af &mapping)
Standard constructor.
Definition: bhole.C:180
void init_bhole_seul()
Initiates for a single the black hole.
Definition: bhole.C:448
const Tenseur & get_psi_tot() const
Returns the total .
Definition: bhole.h:424
void solve_psi_seul(double relax)
Solves the equation for ~: with the condition that on the horizon, where f is the value of on the ...
void solve_lapse_seul(double relax)
Solves the equation for N ~: with the condition that N =0 on the horizon.
const Tenseur & get_n_auto() const
Returns the part of N generated by the hole.
Definition: bhole.h:395
double moment_systeme_inf()
Calculates the angular momentum of the black hole using the formula at infinity : where ...
Definition: bhole_glob.C:173
const Tenseur & get_tkij_tot() const
Returns the total .
Definition: bhole.h:460
double omega_local
local angular velocity
Definition: bhole.h:276
const Tenseur & get_n_comp() const
Returns the part of N generated by the companion hole.
Definition: bhole.h:403
Tenseur psi_auto
Part of generated by the hole.
Definition: bhole.h:290
const Tenseur & get_n_tot() const
Returns the total N .
Definition: bhole.h:407
Tenseur n_comp
Part of N generated by the companion hole.
Definition: bhole.h:287
void init_bhole_binaire()
Initialisation of the system.
Tenseur_sym taij_tot
Total , which must be zero on the horizon of the regularisation on the shift has been done...
Definition: bhole.h:305
void fait_taij_auto()
Calculates the part of generated by shift_auto .
Definition: bhole.C:382
Tenseur_sym taij_auto
Part of generated by the hole.
Definition: bhole.h:299
void fait_n_comp(const Bhole &comp)
Imports the part of N due to the companion hole comp .
Definition: bhole.C:257
double distance_propre(const int nr=65) const
Calculation of the proper distance between the two spheres of inversion, along the x axis...
Definition: bhole_glob.C:265
Black hole.
Definition: bhole.h:268
double get_omega() const
Returns the angular velocity.
Definition: bhole.h:808
double * get_boost() const
Returns the cartesian components of the boost with respect to the reference frame.
Definition: bhole.h:382
Tenseur shift_auto
Part of generated by the hole.
Definition: bhole.h:297
void regularise_seul()
Corrects the shift in the innermost shell, so that it remains and that equals zero on the horizon...
const Tenseur & get_psi_auto() const
Returns the part of generated by the hole.
Definition: bhole.h:412
Tenseur_sym taij_comp
Part of generated by the companion hole.
Definition: bhole.h:300
void operator=(const Bhole_binaire &)
Affectation operator.
const Tenseur & get_taij_comp() const
Returns the part of generated by the companion hole.
Definition: bhole.h:451
void coal(double precis, double relax, int nbre_ome, double search_ome, double m1, double m2, const int sortie=0)
Solves the equation for a particular angular velocity, the systeme being initialized to Misner-Lindqu...
Definition: bhole_coal.C:160
Bhole_binaire(Map_af &mp1, Map_af &mp2)
Standard constructor.
double viriel_seul() const
Computes the viriel error, that is the difference between the ADM and the Komar masses, calculated by the asymptotic behaviours of respectively and N .
void solve_shift_seul(double precis, double relax)
Solves the equation for ~: with on the horizon, being the boost and .
double get_regul() const
Returns the intensity of the regularisation on .
Definition: bhole.h:390
void solve_psi_with_ns(double relax)
Solves the equation for ~: with the condition that on the horizon, where f is the value of on th...
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
double omega
Position of the axis of rotation.
Definition: bhole.h:769
void fait_tkij()
Calculates the total .
void regularise_shift(const Tenseur &comp)
Corrects shift_auto in such a way that the total is equal to zero in the horizon, which should ensure the regularity of , using the stored values of the boost and the angular velocity.
Definition: bhole.C:406
Tenseur psi_comp
Part of generated by the companion hole.
Definition: bhole.h:291
void init_bhole_phi()
Initiates the black hole for a resolution with .
double area() const
Computes the area of the throat.
Definition: bhole.C:548
Map_af & mp
Affine mapping.
Definition: bhole.h:273
Bhole hole1
Black hole one.
Definition: bhole.h:762
void init_phi()
Initiates the system for a resolution using the logarithm of .
void set_omega_local(double ome)
Sets the local angular velocity to ome .
Definition: bhole.h:369
double get_rayon() const
Returns the radius of the horizon.
Definition: bhole.h:347
void solve_psi(double precis, double relax)
Solves the equation for the conformal factor~: The fields are the total values excpet those with s...
double viriel() const
Computes the viriel error, that is the difference between the ADM and the Komar masses, calculated by the asymptotic behaviours of respectively and N .
const Tenseur & get_taij_auto() const
Returns the part of generated by the hole.
Definition: bhole.h:447
void sauve(FILE *fich) const
Write on a file.
Definition: bhole.C:485
double rayon
Radius of the horizon in LORENE's units.
Definition: bhole.h:274
Tenseur_sym tkij_auto
Auto .
Definition: bhole.h:307
const Tenseur & get_grad_n_tot() const
Returns the gradient of N .
Definition: bhole.h:434
double * boost
Lapse on the horizon.
Definition: bhole.h:283
const Map_af & get_mp() const
Returns the mapping (readonly).
Definition: bhole.h:339
double get_rot_state() const
Returns the state of rotation.
Definition: bhole.h:373
void fait_decouple()
Calculates {tt decouple} which is used to obtain tkij_auto by the formula : tkij_auto = decouple * tk...
Definition: bhole_kij.C:317
void set_statiques(double precis, double relax)
Initialize the systeme to Misner Lindquist solution, that is solving for N and in the case ...
Definition: bhole_coal.C:121
Tenseur_sym tkij_tot
Total .
Definition: bhole.h:308
Tenseur grad_psi_tot
Total gradient of .
Definition: bhole.h:295
void set_omega(double ome)
Sets the orbital velocity to ome.
Definition: bhole.h:791
void solve_lapse(double precis, double relax)
Solves the equation for the lapse~: The fields are the total values excpet those with subscript ...
const Cmp get_decouple() const
Returns the function used to construct tkij_auto from tkij_tot .
Definition: bhole.h:468
void init_bhole()
Sets the values of the fields to :
Definition: bhole.C:412
void operator=(const Bhole &)
Affectation.
Definition: bhole.C:221
double get_omega() const
Returns the angular velocity.
Definition: bhole.h:357
Affine radial mapping.
Definition: map.h:2042
Map_af & set_mp()
Read/write of the mapping.
Definition: bhole.h:342
const Bhole & operator()(int i) const
Read only of a component of the system.
Definition: bhole.h:803
int rot_state
State of rotation.
Definition: bhole.h:277
Cmp decouple
Function used to construct the part of generated by the hole from the total .
Definition: bhole.h:318
double masse_komar_seul() const
Calculates the Komar mass of the black hole using : .
void solve_lapse_with_ns(double relax, int bound_nn, double lim_nn)
Solves the equation for N ~: with the condition that N =0 on the horizon.
Definition: bhole_with_ns.C:92
void set_omega(double ome)
Sets the angular velocity to ome .
Definition: bhole.h:361
Tenseur psi_tot
Total .
Definition: bhole.h:292
void solve_shift_with_ns(const Et_bin_nsbh &ns, double precis, double relax, int bound_nn, double lim_nn)
Solves the equation for ~: with on the horizon, being the boost and .
const Tenseur & get_taij_tot() const
Returns the total .
Definition: bhole.h:455
double moment_systeme_hor() const
Calculates the angular momentum of the black hole using the formula on the horizon : where and den...
Definition: bhole_glob.C:111
double adm_systeme() const
Calculates the ADM mass of the system using : .
Definition: bhole_glob.C:89
const Tenseur & get_shift_auto() const
Returns the part of generated by the hole.
Definition: bhole.h:440
double komar_systeme() const
Calculates the Komar mass of the system using : .
Definition: bhole_glob.C:100
Cmp return_psi_auto() const
Returns the part of generated by the hole as a cmp.
Definition: bhole.h:416
double omega
Angular velocity in LORENE's units.
Definition: bhole.h:275
double moment_seul_inf() const
Calculates the angular momentum of the black hole using the formula at infinity : where ...
Binary black holes system.
Definition: bhole.h:757
Tenseur n_auto
Part of N generated by the hole.
Definition: bhole.h:286
Basic array class.
Definition: tbl.h:164
void fait_psi_comp(const Bhole &comp)
Imports the part of due to the companion hole comp .
Definition: bhole.C:283
double masse_adm_seul() const
Calculates the ADM mass of the black hole using : .
void set_rot_state(int rotation)
Returns the state of rotation.
Definition: bhole.h:377
Cmp return_n_auto() const
Returns the part of N generated by the hole as a cmp.
Definition: bhole.h:399
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
double regul
Intensity of the correction on the shift vector.
Definition: bhole.h:284
double get_omega_local() const
Returns the local angular velocity.
Definition: bhole.h:365
Bhole hole2
Black hole two.
Definition: bhole.h:763