|
LORENE
|
Binary black holes system. More...
#include <bhole.h>
Public Member Functions | |
| Bhole_binaire (Map_af &mp1, Map_af &mp2) | |
| Standard constructor. | |
| Bhole_binaire (const Bhole_binaire &) | |
| Copy. | |
| ~Bhole_binaire () | |
| Destructor. | |
| void | operator= (const Bhole_binaire &) |
| Affectation operator. | |
| Bhole & | set (int i) |
| Read/write of a component of the system. | |
| void | set_omega (double ome) |
Sets the orbital velocity to ome. | |
| void | set_pos_axe (double) |
| const Bhole & | operator() (int i) const |
| Read only of a component of the system. | |
| double | get_omega () const |
| Returns the angular velocity. | |
| void | init_bhole_binaire () |
| Initialisation of the system. | |
| 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 ![]() | |
| void | solve_lapse (double precis, double relax) |
| Solves the equation for the lapse~: | |
| void | solve_psi (double precis, double relax) |
| Solves the equation for the conformal factor~: | |
| void | solve_shift (double precis, double relax) |
| Solves the equation for the shift, using the Oohara-Nakarmure scheme~: | |
| void | fait_tkij () |
| Calculation af the extrinsic curvature tensor. | |
| void | fait_decouple () |
Calculates {tt decouple} which is used to obtain tkij_auto by the formula : tkij_auto = decouple * tkij_tot . | |
| void | set_statiques (double precis, double relax) |
Initialize the systeme to Misner Lindquist solution, that is solving for N and ![]() ![]() | |
| 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-Lindquist solution. | |
| double | adm_systeme () const |
Calculates the ADM mass of the system using : ![]() | |
| double | komar_systeme () const |
Calculates the Komar mass of the system using : ![]() | |
| double | moment_systeme_inf () |
Calculates the angular momentum of the black hole using the formula at infinity : ![]() ![]() | |
| double | moment_systeme_hor () const |
Calculates the angular momentum of the black hole using the formula on the horizon : ![]() ![]() ![]() | |
| double | distance_propre (const int nr=65) const |
| Calculation of the proper distance between the two spheres of inversion, along the x axis. | |
| Tbl | linear_momentum_systeme_inf () const |
| void | solve_phi (double precision, double relax) |
Solve the equation for the logarithm of ![]() | |
| void | init_phi () |
Initiates the system for a resolution using the logarithm of ![]() | |
Private Attributes | |
| Bhole | hole1 |
| Black hole one. | |
| Bhole | hole2 |
| Black hole two. | |
| Bhole * | holes [2] |
| Array on the black holes. | |
| double | pos_axe |
| double | omega |
| Position of the axis of rotation. | |
Binary black holes system.
()
This class is intended for dealing with binary black holes configurations in the conformaly flat approximation.
Standard constructor.
Definition at line 211 of file bhole_binaire.C.
References Lorene::get_ori_x(), hole1, hole2, holes, and omega.
| Lorene::Bhole_binaire::Bhole_binaire | ( | const Bhole_binaire & | source | ) |
Copy.
Definition at line 219 of file bhole_binaire.C.
References Bhole_binaire(), hole1, hole2, holes, and omega.
| Lorene::Bhole_binaire::~Bhole_binaire | ( | ) |
Destructor.
Definition at line 226 of file bhole_binaire.C.
| double Lorene::Bhole_binaire::adm_systeme | ( | ) | const |
Calculates the ADM mass of the system using : 
Definition at line 86 of file bhole_glob.C.
| void Lorene::Bhole_binaire::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-Lindquist solution.
| precis | [input] : precision for the convergence (on ![]() |
| relax | [input] : relaxation parameter. |
| nbre_ome | [input] : number of intermediates velocities to go from 0 to omega , typically 10. |
| sortie | [input] : flag for the output on files (0 no output files). |
Definition at line 157 of file bhole_coal.C.
References Lorene::diffrelmax(), fait_tkij(), hole1, hole2, omega, Lorene::pow(), set_omega(), solve_lapse(), solve_psi(), solve_shift(), Lorene::sqrt(), and viriel().
| double Lorene::Bhole_binaire::distance_propre | ( | const int | nr = 65 | ) | const |
Calculation of the proper distance between the two spheres of inversion, along the x axis.
| nr | [input] : number of points used for the calculation. |
Definition at line 262 of file bhole_glob.C.
References Lorene::cos(), hole1, hole2, and Lorene::pow().
| void Lorene::Bhole_binaire::fait_decouple | ( | ) |
Calculates {tt decouple} which is used to obtain tkij_auto by the formula : tkij_auto = decouple * tkij_tot .
(see the membre {tt Cmp decouple for more precisions about its value).
Definition at line 314 of file bhole_kij.C.
References Lorene::Cmp::allocate_all(), Lorene::cos(), hole1, hole2, Lorene::phi, Lorene::pow(), Lorene::Cmp::set(), Lorene::sin(), Lorene::Cmp::std_base_scal(), and Lorene::Cmp::val_point().
| void Lorene::Bhole_binaire::fait_tkij | ( | ) |
Calculation af the extrinsic curvature tensor.
The regularisation of the shifts must have been done. All the contributions of 
Definition at line 88 of file bhole_kij.C.
References Lorene::Cmp::dec2_dzpuis(), Lorene::Tenseur::dec2_dzpuis(), hole1, hole2, Lorene::Cmp::inc2_dzpuis(), Lorene::phi, Lorene::Cmp::raccord(), Lorene::Cmp::set(), and Lorene::Cmp::val_point().
|
inline |
| void Lorene::Bhole_binaire::init_bhole_binaire | ( | ) |
Initialisation of the system.
Each hole is set close to a Schwarzschild one and the parts of the fields generated by the companion are calculated.
The angular velocity is set to zero.
Definition at line 239 of file bhole_binaire.C.
References fait_decouple(), hole1, hole2, and set_omega().
| void Lorene::Bhole_binaire::init_phi | ( | ) |
Initiates the system for a resolution using the logarithm of 
Definition at line 155 of file bhole_solve_phi.C.
References hole1, hole2, and set_omega().
| double Lorene::Bhole_binaire::komar_systeme | ( | ) | const |
Calculates the Komar mass of the system using : 
Definition at line 97 of file bhole_glob.C.
| Tbl Lorene::Bhole_binaire::linear_momentum_systeme_inf | ( | ) | const |
Definition at line 319 of file bhole_glob.C.
| double Lorene::Bhole_binaire::moment_systeme_hor | ( | ) | const |
Calculates the angular momentum of the black hole using the formula on the horizon : 


Definition at line 108 of file bhole_glob.C.
References Lorene::Tenseur::annule(), Lorene::Tenseur::change_triad(), hole1, hole2, omega, Lorene::pow(), Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), and Lorene::Cmp::std_base_scal().
| double Lorene::Bhole_binaire::moment_systeme_inf | ( | ) |
Calculates the angular momentum of the black hole using the formula at infinity : 

Definition at line 170 of file bhole_glob.C.
References Lorene::Tenseur::annule(), Lorene::Tenseur::change_triad(), Lorene::Tenseur::dec2_dzpuis(), Lorene::Tenseur::gradient(), hole1, hole2, Lorene::Cmp::import_asymy(), Lorene::Cmp::import_symy(), Lorene::Tenseur::inc2_dzpuis(), Lorene::Map_af::integrale_surface_infini(), Lorene::Valeur::mult_cp(), Lorene::Cmp::mult_r_zec(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), omega, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), and Lorene::Cmp::va.
|
inline |
| void Lorene::Bhole_binaire::operator= | ( | const Bhole_binaire & | source | ) |
Affectation operator.
Definition at line 230 of file bhole_binaire.C.
References Bhole_binaire(), hole1, hole2, and omega.
|
inline |
|
inline |
| void Lorene::Bhole_binaire::set_pos_axe | ( | double | new_pos | ) |
Definition at line 254 of file bhole_binaire.C.
| void Lorene::Bhole_binaire::set_statiques | ( | double | precis, |
| double | relax ) |
Initialize the systeme to Misner Lindquist solution, that is solving for N and 

| precis | [input] : precision for the convergence (on N ). |
| relax | [input] : relaxation parameter. |
Definition at line 118 of file bhole_coal.C.
References Lorene::diffrelmax(), hole1, hole2, init_bhole_binaire(), set_omega(), solve_lapse(), and solve_psi().
| void Lorene::Bhole_binaire::solve_lapse | ( | double | precis, |
| double | relax ) |
Solves the equation for the lapse~:
![\[\Delta N_a = -\frac{2}{\Psi}\nabla_i \Psi_a \nabla^i N + N \Psi^4 K_{ij}K_a^{ij}
\]](form_78.png)
The fields are the total values excpet those with subscript 
| precis | [input] : precision, for the boudary conditions are obtained by iteration. |
| relax | [input] : relaxation parameter. |
Definition at line 84 of file bhole_equations_bin.C.
References Lorene::flat_scalar_prod(), hole1, hole2, Lorene::pow(), and Lorene::Cmp::std_base_scal().
| void Lorene::Bhole_binaire::solve_phi | ( | double | precision, |
| double | relax ) |
Solve the equation for the logarithm of 
Definition at line 108 of file bhole_solve_phi.C.
References Lorene::diffrelmax(), Lorene::flat_scalar_prod(), hole1, hole2, Lorene::Cmp::std_base_scal(), and Lorene::Valeur::std_base_scal().
| void Lorene::Bhole_binaire::solve_psi | ( | double | precis, |
| double | relax ) |
Solves the equation for the conformal factor~:
![\[*\Delta \Psi_a = -\frac {\Psi^5}{8} K_{ij}K_a^{ij}
\]](form_80.png)
The fields are the total values excpet those with subscript 



| precis | [input] : precision, for the boudary conditions are being obtained by iteration. |
| relax | [input] : relaxation parameter. |
Definition at line 134 of file bhole_equations_bin.C.
References Lorene::flat_scalar_prod(), hole1, hole2, Lorene::pow(), Lorene::Valeur::set(), Lorene::Cmp::std_base_scal(), and Lorene::Valeur::std_base_scal().
| void Lorene::Bhole_binaire::solve_shift | ( | double | precis, |
| double | relax ) |
Solves the equation for the shift, using the Oohara-Nakarmure scheme~:
![\[*\Delta \beta_a^i +\frac{1}{3} \nabla^i \left(\nabla_j \beta_a^j\right) =
-6A^{ij}\frac{\nabla_i \Psi_a}{\Psi} + 2K^{ij}\nabla_j N_a
\]](form_82.png)
using the current 


a = 1, 2). The companion contributions are then obtained.
| precis | [input] : precision for the solver, the boundary conditions and the inversion of the operator being obtained by iteration. |
| relax | [input] : relaxation parameter. |
Definition at line 196 of file bhole_equations_bin.C.
References Lorene::Valeur::base, Lorene::Cmp::filtre(), Lorene::flat_scalar_prod(), Lorene::Tenseur::get_etat(), Lorene::Tenseur::get_mp(), hole1, hole2, omega, Lorene::Tenseur::set(), Lorene::Valeur::set(), Lorene::Valeur::set_etat_c_qcq(), Lorene::Mtbl::set_etat_qcq(), Lorene::Tenseur::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Tenseur::set_etat_zero(), and Lorene::Tenseur::set_std_base().
| double Lorene::Bhole_binaire::viriel | ( | ) | const |
Computes the viriel error, that is the difference between the ADM and the Komar masses, calculated by the asymptotic behaviours of respectively 
Definition at line 93 of file bhole_pseudo_viriel.C.
|
private |
|
private |