LORENE
|
Binary black holes system. More...
#include <bhole.h>
Public Member Functions | |
Bhole_binaire (Map_af &mp1, Map_af &mp2) | |
Standard constructor. More... | |
Bhole_binaire (const Bhole_binaire &) | |
Copy. More... | |
~Bhole_binaire () | |
Destructor. More... | |
void | operator= (const Bhole_binaire &) |
Affectation operator. More... | |
Bhole & | set (int i) |
Read/write of a component of the system. More... | |
void | set_omega (double ome) |
Sets the orbital velocity to ome . More... | |
void | set_pos_axe (double) |
const Bhole & | operator() (int i) const |
Read only of a component of the system. More... | |
double | get_omega () const |
Returns the angular velocity. More... | |
void | init_bhole_binaire () |
Initialisation of the system. More... | |
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 . More... | |
void | solve_lapse (double precis, double relax) |
Solves the equation for the lapse~:
The fields are the total values excpet those with subscript , which are the fields generated by each holes (a = 1, 2). More... | |
void | solve_psi (double precis, double relax) |
Solves the equation for the conformal factor~:
The fields are the total values excpet those with subscript , which are the fields generated by each holes (a = 1, 2). More... | |
void | solve_shift (double precis, double relax) |
Solves the equation for the shift, using the Oohara-Nakarmure scheme~:
using the current for the boudary conditions~: on both horizons. More... | |
void | fait_tkij () |
Calculation af the extrinsic curvature tensor. More... | |
void | fait_decouple () |
Calculates {tt decouple} which is used to obtain tkij_auto by the formula : tkij_auto = decouple * tkij_tot . More... | |
void | set_statiques (double precis, double relax) |
Initialize the systeme to Misner Lindquist solution, that is solving for N and in the case . More... | |
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. More... | |
double | adm_systeme () const |
Calculates the ADM mass of the system using : . More... | |
double | komar_systeme () const |
Calculates the Komar mass of the system using : . More... | |
double | moment_systeme_inf () |
Calculates the angular momentum of the black hole using the formula at infinity : where . More... | |
double | moment_systeme_hor () const |
Calculates the angular momentum of the black hole using the formula on the horizon : where and denotes the horizon of the hole a . More... | |
double | distance_propre (const int nr=65) const |
Calculation of the proper distance between the two spheres of inversion, along the x axis. More... | |
Tbl | linear_momentum_systeme_inf () const |
void | solve_phi (double precision, double relax) |
Solve the equation for the logarithm of . More... | |
void | init_phi () |
Initiates the system for a resolution using the logarithm of . More... | |
Private Attributes | |
Bhole | hole1 |
Black hole one. More... | |
Bhole | hole2 |
Black hole two. More... | |
Bhole * | holes [2] |
Array on the black holes. More... | |
double | pos_axe |
double | omega |
Position of the axis of rotation. More... | |
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 214 of file bhole_binaire.C.
Lorene::Bhole_binaire::Bhole_binaire | ( | const Bhole_binaire & | source | ) |
Lorene::Bhole_binaire::~Bhole_binaire | ( | ) |
Destructor.
Definition at line 229 of file bhole_binaire.C.
double Lorene::Bhole_binaire::adm_systeme | ( | ) | const |
Calculates the ADM mass of the system using : .
Definition at line 89 of file bhole_glob.C.
References hole1, hole2, Lorene::Map_af::integrale_surface_infini(), Lorene::Bhole::mp, and Lorene::Bhole::psi_auto.
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 160 of file bhole_coal.C.
References Lorene::Bhole::area(), Lorene::diffrelmax(), fait_tkij(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), Lorene::Bhole::get_regul(), Lorene::Bhole::get_shift_auto(), hole1, hole2, Lorene::Map_af::homothetie_interne(), Lorene::Bhole::mp, omega, Lorene::pow(), Lorene::Bhole::rayon, Lorene::Bhole::regul, set_omega(), solve_lapse(), solve_psi(), solve_shift(), Lorene::sqrt(), Lorene::Map_af::val_r(), 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 265 of file bhole_glob.C.
References Lorene::Map::convert_absolute(), Lorene::cos(), Lorene::Map::get_ori_x(), hole1, hole2, Lorene::Bhole::mp, Lorene::pow(), Lorene::Bhole::psi_auto, and Lorene::Bhole::rayon.
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 317 of file bhole_kij.C.
References Lorene::Cmp::allocate_all(), Lorene::Map::convert_absolute(), Lorene::cos(), Lorene::Bhole::decouple, Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), hole1, hole2, Lorene::Bhole::mp, Lorene::pow(), Lorene::Map::r, Lorene::Cmp::set(), Lorene::sin(), Lorene::Cmp::std_base_scal(), Lorene::Cmp::val_point(), Lorene::Map::xa, Lorene::Map::ya, and Lorene::Map::za.
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 are then calculated and the total tensor must be zero on both horizons. The computation is done to avoid every singularity close to the horizons (division by N =0) for it is done in the coefficient space in the two regions surrounding the holes.
Definition at line 91 of file bhole_kij.C.
References Lorene::Tenseur::dec2_dzpuis(), Lorene::Bhole::fait_taij_auto(), Lorene::Tenseur::get_etat(), Lorene::Tenseur::get_mp(), Lorene::Map::get_rot_phi(), hole1, hole2, Lorene::Tenseur::set_etat_qcq(), Lorene::Bhole::taij_auto, Lorene::Bhole::taij_comp, and Lorene::Bhole::tkij_tot.
|
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 242 of file bhole_binaire.C.
References fait_decouple(), Lorene::Bhole::fait_n_comp(), Lorene::Bhole::fait_psi_comp(), hole1, hole2, Lorene::Bhole::init_bhole(), and set_omega().
void Lorene::Bhole_binaire::init_phi | ( | ) |
Initiates the system for a resolution using the logarithm of .
Definition at line 158 of file bhole_solve_phi.C.
References Lorene::Bhole::fait_psi_comp(), hole1, hole2, Lorene::Bhole::init_bhole_phi(), and set_omega().
double Lorene::Bhole_binaire::komar_systeme | ( | ) | const |
Calculates the Komar mass of the system using : .
Definition at line 100 of file bhole_glob.C.
References hole1, hole2, Lorene::Map_af::integrale_surface_infini(), Lorene::Bhole::mp, and Lorene::Bhole::n_auto.
double Lorene::Bhole_binaire::moment_systeme_hor | ( | ) | const |
Calculates the angular momentum of the black hole using the formula on the horizon : where and denotes the horizon of the hole a .
Definition at line 111 of file bhole_glob.C.
References Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_rot_phi(), hole1, hole2, Lorene::Map_af::integrale_surface(), Lorene::Bhole::mp, omega, Lorene::pow(), Lorene::Bhole::psi_tot, Lorene::Bhole::rayon, Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), Lorene::Cmp::std_base_scal(), Lorene::Bhole::tkij_tot, Lorene::Map::xa, and Lorene::Map::ya.
double Lorene::Bhole_binaire::moment_systeme_inf | ( | ) |
Calculates the angular momentum of the black hole using the formula at infinity : where .
Definition at line 173 of file bhole_glob.C.
References Lorene::Tenseur::annule(), Lorene::Tenseur::dec2_dzpuis(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), Lorene::Map::get_rot_phi(), Lorene::Tenseur::gradient(), hole1, hole2, Lorene::Tenseur::inc2_dzpuis(), Lorene::Bhole::mp, Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), omega, Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), Lorene::Bhole::shift_auto, and Lorene::Cmp::va.
|
inline |
void Lorene::Bhole_binaire::operator= | ( | const Bhole_binaire & | source | ) |
Affectation operator.
Definition at line 233 of file bhole_binaire.C.
|
inline |
|
inline |
void Lorene::Bhole_binaire::set_statiques | ( | double | precis, |
double | relax | ||
) |
Initialize the systeme to Misner Lindquist solution, that is solving for N and in the case .
precis | [input] : precision for the convergence (on N ). |
relax | [input] : relaxation parameter. |
Definition at line 121 of file bhole_coal.C.
References Lorene::diffrelmax(), Lorene::Map::get_mg(), Lorene::Bhole::get_n_auto(), Lorene::Mg3d::get_nzone(), hole1, hole2, init_bhole_binaire(), Lorene::Bhole::mp, set_omega(), solve_lapse(), and solve_psi().
void Lorene::Bhole_binaire::solve_lapse | ( | double | precis, |
double | relax | ||
) |
Solves the equation for the lapse~:
The fields are the total values excpet those with subscript , which are the fields generated by each holes (a = 1, 2).
The boundary conditions are such that N =0 on both horizons. The companion contributions are then obtained.
precis | [input] : precision, for the boudary conditions are obtained by iteration. |
relax | [input] : relaxation parameter. |
Definition at line 87 of file bhole_equations_bin.C.
References Lorene::Bhole::fait_n_comp(), Lorene::flat_scalar_prod(), Lorene::Bhole::grad_n_tot, Lorene::Tenseur::gradient(), hole1, hole2, Lorene::Bhole::n_auto, Lorene::Bhole::n_tot, Lorene::pow(), Lorene::Bhole::psi_auto, Lorene::Bhole::psi_tot, Lorene::Cmp::raccord(), Lorene::Tenseur::set(), Lorene::Cmp::std_base_scal(), Lorene::Bhole::tkij_auto, and Lorene::Bhole::tkij_tot.
void Lorene::Bhole_binaire::solve_phi | ( | double | precision, |
double | relax | ||
) |
Solve the equation for the logarithm of .
Definition at line 111 of file bhole_solve_phi.C.
References Lorene::diffrelmax(), Lorene::Bhole::fait_psi_comp(), Lorene::flat_scalar_prod(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Bhole::grad_psi_tot, Lorene::Tenseur::gradient(), hole1, hole2, Lorene::Bhole::mp, Lorene::Bhole::psi_auto, Lorene::Cmp::raccord(), Lorene::Bhole::rayon, Lorene::Tenseur::set(), and Lorene::Cmp::std_base_scal().
void Lorene::Bhole_binaire::solve_psi | ( | double | precis, |
double | relax | ||
) |
Solves the equation for the conformal factor~:
The fields are the total values excpet those with subscript , which are the fields generated by each holes (a = 1, 2).
The boundary conditions are such that on both horizons, being the radii of those horizons and f the value of on the horizons at the preceeding step. The companion contributions are then obtained.
precis | [input] : precision, for the boudary conditions are being obtained by iteration. |
relax | [input] : relaxation parameter. |
Definition at line 137 of file bhole_equations_bin.C.
References Lorene::Bhole::fait_psi_comp(), Lorene::flat_scalar_prod(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), hole1, hole2, Lorene::Bhole::mp, Lorene::pow(), Lorene::Bhole::psi_auto, Lorene::Bhole::psi_tot, Lorene::Cmp::raccord(), Lorene::Bhole::rayon, Lorene::Valeur::set(), Lorene::Tenseur::set(), Lorene::Cmp::std_base_scal(), Lorene::Bhole::tkij_auto, and Lorene::Bhole::tkij_tot.
void Lorene::Bhole_binaire::solve_shift | ( | double | precis, |
double | relax | ||
) |
Solves the equation for the shift, using the Oohara-Nakarmure scheme~:
using the current for the boudary conditions~: on both horizons.
The fields are the total values excpet those with subscript , which are the fields generated by each holes (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 199 of file bhole_equations_bin.C.
References Lorene::flat_scalar_prod(), Lorene::Tenseur::get_etat(), Lorene::Tenseur::gradient(), hole1, Lorene::Bhole::n_auto, omega, Lorene::Bhole::psi_auto, Lorene::Bhole::psi_tot, Lorene::Bhole::taij_tot, and Lorene::Bhole::tkij_tot.
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 and N .
Definition at line 96 of file bhole_pseudo_viriel.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), hole1, hole2, Lorene::Bhole::mp, Lorene::Bhole::n_auto, and Lorene::Bhole::psi_auto.
|
private |
|
private |