LORENE
Lorene::Bhole_binaire Class Reference

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...
 
Bholeset (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 Bholeoperator() (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 $\Psi$ and N . More...
 
void 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} \]

The fields are the total values excpet those with subscript $_a$, 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~:

\[ *\Delta \Psi_a = -\frac {\Psi^5}{8} K_{ij}K_a^{ij} \]

The fields are the total values excpet those with subscript $_a$, 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~:

\[ *\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 \]

using the current $\Omega$ for the boudary conditions~: $\vec{N} = -\Omega \vec{m}$ 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 $Psi$ in the case $\Omega = 0$. 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 : $M = -\frac{1}{2 \pi} \oint_{\infty} \nabla^i \Psi {\mathrm d} S_i$. More...
 
double komar_systeme () const
 Calculates the Komar mass of the system using : $M = \frac{1}{4 \pi} \oint_{\infty} \nabla^i N {\mathrm d} S_i$. More...
 
double moment_systeme_inf ()
 Calculates the angular momentum of the black hole using the formula at infinity : $ J = -\frac{1}{8 \pi} \oint_{\infty} K_j^i m^j {\mathrm d}S_i$ where $\vec{m} = \frac{\partial}{\partial \phi}$. More...
 
double moment_systeme_hor () const
 Calculates the angular momentum of the black hole using the formula on the horizon : $ J = -\sum_{a= 1, 2} \frac{1}{8 \pi} \oint_{H_a} \Psi^6 K_j^i m^j {\mathrm d}S_i$ where $\vec{m} = \frac{\partial}{\partial \phi}$ and $H_a$ 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 $\Psi$. More...
 
void init_phi ()
 Initiates the system for a resolution using the logarithm of $\Psi$. More...
 

Private Attributes

Bhole hole1
 Black hole one. More...
 
Bhole hole2
 Black hole two. More...
 
Bholeholes [2]
 Array on the black holes. More...
 
double pos_axe
 
double omega
 Position of the axis of rotation. More...
 

Detailed Description

Binary black holes system.

()

This class is intended for dealing with binary black holes configurations in the conformaly flat approximation.

Definition at line 757 of file bhole.h.

Constructor & Destructor Documentation

◆ Bhole_binaire() [1/2]

Lorene::Bhole_binaire::Bhole_binaire ( Map_af mp1,
Map_af mp2 
)

Standard constructor.

Definition at line 214 of file bhole_binaire.C.

References hole1, hole2, and holes.

◆ Bhole_binaire() [2/2]

Lorene::Bhole_binaire::Bhole_binaire ( const Bhole_binaire source)

Copy.

Definition at line 222 of file bhole_binaire.C.

References hole1, hole2, and holes.

◆ ~Bhole_binaire()

Lorene::Bhole_binaire::~Bhole_binaire ( )

Destructor.

Definition at line 229 of file bhole_binaire.C.

Member Function Documentation

◆ adm_systeme()

double Lorene::Bhole_binaire::adm_systeme ( ) const

Calculates the ADM mass of the system using : $M = -\frac{1}{2 \pi} \oint_{\infty} \nabla^i \Psi {\mathrm d} S_i$.

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.

◆ coal()

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.

Parameters
precis[input] : precision for the convergence (on $\beta$).
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).
Returns
: the virial error.

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().

◆ distance_propre()

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.

Parameters
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.

◆ fait_decouple()

void Lorene::Bhole_binaire::fait_decouple ( )

◆ fait_tkij()

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 $A^{ij}$ 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.

◆ get_omega()

double Lorene::Bhole_binaire::get_omega ( ) const
inline

Returns the angular velocity.

Definition at line 808 of file bhole.h.

References omega.

◆ init_bhole_binaire()

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().

◆ init_phi()

void Lorene::Bhole_binaire::init_phi ( )

Initiates the system for a resolution using the logarithm of $\Psi$.

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().

◆ komar_systeme()

double Lorene::Bhole_binaire::komar_systeme ( ) const

Calculates the Komar mass of the system using : $M = \frac{1}{4 \pi} \oint_{\infty} \nabla^i N {\mathrm d} S_i$.

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.

◆ moment_systeme_hor()

double Lorene::Bhole_binaire::moment_systeme_hor ( ) const

◆ moment_systeme_inf()

◆ operator()()

const Bhole& Lorene::Bhole_binaire::operator() ( int  i) const
inline

Read only of a component of the system.

i must be equal to 1 or 2.

Definition at line 803 of file bhole.h.

References holes.

◆ operator=()

void Lorene::Bhole_binaire::operator= ( const Bhole_binaire source)

Affectation operator.

Definition at line 233 of file bhole_binaire.C.

References hole1, hole2, and omega.

◆ set()

Bhole& Lorene::Bhole_binaire::set ( int  i)
inline

Read/write of a component of the system.

i must be equal to 1 or 2.

Definition at line 785 of file bhole.h.

References holes.

◆ set_omega()

void Lorene::Bhole_binaire::set_omega ( double  ome)
inline

Sets the orbital velocity to ome.

Definition at line 791 of file bhole.h.

References hole1, hole2, omega, and Lorene::Bhole::set_omega().

◆ set_statiques()

void Lorene::Bhole_binaire::set_statiques ( double  precis,
double  relax 
)

Initialize the systeme to Misner Lindquist solution, that is solving for N and $Psi$ in the case $\Omega = 0$.

Parameters
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().

◆ solve_lapse()

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} \]

The fields are the total values excpet those with subscript $_a$, 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.

Parameters
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.

◆ solve_phi()

◆ solve_psi()

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} \]

The fields are the total values excpet those with subscript $_a$, which are the fields generated by each holes (a = 1, 2).

The boundary conditions are such that $\partial_r \Psi = -\frac{1}{2 r_0} f\left(\theta, \phi\right)$ on both horizons, $r_0$ being the radii of those horizons and f the value of $\Psi$ on the horizons at the preceeding step. The companion contributions are then obtained.

Parameters
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.

◆ solve_shift()

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 \]

using the current $\Omega$ for the boudary conditions~: $\vec{N} = -\Omega \vec{m}$ on both horizons.

The fields are the total values excpet those with subscript $_a$, which are the fields generated by each holes (a = 1, 2). The companion contributions are then obtained.

Parameters
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.

◆ viriel()

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 $\Psi$ 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.

Member Data Documentation

◆ hole1

Bhole Lorene::Bhole_binaire::hole1
private

Black hole one.

Definition at line 762 of file bhole.h.

◆ hole2

Bhole Lorene::Bhole_binaire::hole2
private

Black hole two.

Definition at line 763 of file bhole.h.

◆ holes

Bhole* Lorene::Bhole_binaire::holes[2]
private

Array on the black holes.

Definition at line 766 of file bhole.h.

◆ omega

double Lorene::Bhole_binaire::omega
private

Position of the axis of rotation.

Angular velocity

Definition at line 769 of file bhole.h.


The documentation for this class was generated from the following files: