LORENE
Lorene::Bhole Class Reference

Black hole. More...

#include <bhole.h>

Public Member Functions

 Bhole (Map_af &mapping)
 Standard constructor. More...
 
 Bhole (const Bhole &)
 Constructor by copy. More...
 
 Bhole (Map_af &, FILE *, bool old=false)
 Constructor from a Map_af and a file. More...
 
 ~Bhole ()
 Destructor. More...
 
void sauve (FILE *fich) const
 Write on a file. More...
 
void operator= (const Bhole &)
 Affectation. More...
 
const Map_afget_mp () const
 Returns the mapping (readonly). More...
 
Map_afset_mp ()
 Read/write of the mapping. More...
 
double get_rayon () const
 Returns the radius of the horizon. More...
 
void set_rayon (double ray)
 Sets the radius of the horizon to ray . More...
 
double get_omega () const
 Returns the angular velocity. More...
 
void set_omega (double ome)
 Sets the angular velocity to ome . More...
 
double get_omega_local () const
 Returns the local angular velocity. More...
 
void set_omega_local (double ome)
 Sets the local angular velocity to ome . More...
 
double get_rot_state () const
 Returns the state of rotation. More...
 
void set_rot_state (int rotation)
 Returns the state of rotation. More...
 
double * get_boost () const
 Returns the cartesian components of the boost with respect to the reference frame. More...
 
void set_boost (double vx, double vy, double vz)
 
double get_regul () const
 Returns the intensity of the regularisation on $\vec{N}$. More...
 
const Tenseurget_n_auto () const
 Returns the part of N generated by the hole. More...
 
Cmp return_n_auto () const
 Returns the part of N generated by the hole as a cmp. More...
 
const Tenseurget_n_comp () const
 Returns the part of N generated by the companion hole. More...
 
const Tenseurget_n_tot () const
 Returns the total N . More...
 
const Tenseurget_psi_auto () const
 Returns the part of $\Psi$ generated by the hole. More...
 
Cmp return_psi_auto () const
 Returns the part of $\Psi$ generated by the hole as a cmp. More...
 
const Tenseurget_psi_comp () const
 Returns the part of $\Psi$ generated by the companion hole. More...
 
const Tenseurget_psi_tot () const
 Returns the total $\Psi$. More...
 
const Tenseurget_grad_psi_tot () const
 Returns the gradient of $\Psi$. More...
 
const Tenseurget_grad_n_tot () const
 Returns the gradient of N . More...
 
const Tenseurget_shift_auto () const
 Returns the part of $\beta^i$ generated by the hole. More...
 
Cmp return_shift_auto (int i) const
 
const Tenseurget_taij_auto () const
 Returns the part of $A^{ij}$ generated by the hole. More...
 
const Tenseurget_taij_comp () const
 Returns the part of $A^{ij}$ generated by the companion hole. More...
 
const Tenseurget_taij_tot () const
 Returns the total $A^{ij}$. More...
 
const Tenseurget_tkij_tot () const
 Returns the total $K^{ij}$. More...
 
const Tenseurget_tkij_auto () const
 Returns the part of $K^{ij}$ generated by the hole. More...
 
const Cmp get_decouple () const
 Returns the function used to construct tkij_auto from tkij_tot . More...
 
Cmpset_n_auto ()
 
Cmpset_psi_auto ()
 
void set_lapse_hori (double xx)
 
void fait_n_comp (const Bhole &comp)
 Imports the part of N due to the companion hole comp . More...
 
void fait_psi_comp (const Bhole &comp)
 Imports the part of $\Psi$ due to the companion hole comp . More...
 
void fait_n_comp (const Et_bin_nsbh &comp)
 
Imports the part of N due to the companion ns comp . More...
 
void fait_psi_comp (const Et_bin_nsbh &comp)
 Imports the part of $\Psi$ due to the companion ns comp . More...
 
void fait_taij_auto ()
 Calculates the part of $A^{ij}$ generated by shift_auto . More...
 
void regularise_shift (const Tenseur &comp)
 Corrects shift_auto in such a way that the total $A^{ij}$ is equal to zero in the horizon, which should ensure the regularity of $K^{ij}$, using the stored values of the boost and the angular velocity. More...
 
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 $\Psi$ and N . More...
 
void init_bhole ()
 Sets the values of the fields to : More...
 
void init_kerr (double masse, double moment)
 Set the inital values to those of Kerr. More...
 
void solve_lapse_seul (double relax)
 
Solves the equation for N ~:

\[ *\Delta N = -\frac{2}{\Psi}\nabla_i \Psi \nabla^i N + N \Psi^4 K_{ij} K^{ij} \]

with the condition that N =0 on the horizon. More...

 
void solve_psi_seul (double relax)
 Solves the equation for $\Psi$~:

\[ \Delta \Psi = - \frac{\Psi^5} {8} K_{ij}K^{ij} \]

with the condition that $\partial_r \Psi =-\frac{1}{2 {\tt rayon}} f\left(\theta, \phi\right)$ on the horizon, where f is the value of $\Psi$ on the horizon at the preceeding step. More...

 
void solve_shift_seul (double precis, double relax)
 
Solves the equation for $\vec{\beta}$~:

\[ *\Delta \beta^i +\frac{1}{3} \nabla^i \left(\nabla_j \beta^j\right) = -6A^{ij}\frac{\nabla_j \Psi}{\Psi} + 2 K^{ij}\nabla_j N \]

with $\vec{\beta} = -\Omega \vec{m} - \vec{V}$ on the horizon, $\vec{V}$ being the boost and $\vec{m} = \frac{\partial}{\partial \phi}$. More...

 
void regularise_seul ()
 Corrects the shift in the innermost shell, so that it remains $ {\mathcal{C}}^2$ and that $A^{ij}$ equals zero on the horizon. More...
 
void solve_lapse_with_ns (double relax, int bound_nn, double lim_nn)
 
Solves the equation for N ~:

\[ *\Delta N = -\frac{2}{\Psi}\nabla_i \Psi \nabla^i N + N \Psi^4 K_{ij} K^{ij} \]

with the condition that N =0 on the horizon. More...

 
void solve_psi_with_ns (double relax)
 
Solves the equation for $\Psi$~:

\[ (\Delta \Psi = - \frac{\Psi^5} {8} K_{ij}K^{ij} \]

with the condition that $\partial_r \Psi=-\frac{1}{2 {\tt rayon}} f\left(\theta, \phi\right)$ on the horizon, where f is the value of $\Psi$ on the horizon at the preceeding step. More...

 
void solve_shift_with_ns (const Et_bin_nsbh &ns, double precis, double relax, int bound_nn, double lim_nn)
 
Solves the equation for $\vec{\beta}$~:

\[ *\Delta \beta^i +\frac{1}{3} \nabla^i \left(\nabla_j \beta^j\right) = -6A^{ij}\frac{\nabla_j \Psi}{\Psi} + 2 K^{ij}\nabla_j N \]

with $\vec{\beta} = -\Omega \vec{m} - \vec{V}$ on the horizon, $\vec{V}$ being the boost and $\vec{m} = \frac{\partial}{\partial \phi}$. More...

 
void equilibrium (const Et_bin_nsbh &ns, double precis, double relax, int bound_nn, double lim_nn)
 
void update_metric (const Et_bin_nsbh &ns)
 
void fait_tkij ()
 Calculates the total $K^{ij}$. More...
 
double area () const
 Computes the area of the throat. More...
 
double masse_adm_seul () const
 Calculates the ADM mass of the black hole using : $M = -\frac{1}{2 \pi} \oint_{\infty} \nabla^i \Psi {\mathrm d} S_i$. More...
 
double masse_komar_seul () const
 Calculates the Komar mass of the black hole using : $M = \frac{1}{4 \pi} \oint_{\infty} \nabla^i N {\mathrm d} S_i$. More...
 
double moment_seul_inf () const
 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_seul_hor () const
 Calculates the angular momentum of the black hole using the formula on the horizon : $ J = -\frac{1}{8 \pi} \oint_{H} \Psi^6 K_j^i m^j {\mathrm d}S_i$ where $\vec{m} = \frac{\partial}{\partial \phi}$ and H denotes the horizon. More...
 
double local_momentum () const
 
void init_bhole_phi ()
 Initiates the black hole for a resolution with $\Phi = \log \Psi$. More...
 
void init_bhole_seul ()
 Initiates for a single the black hole. More...
 

Protected Attributes

Map_afmp
 Affine mapping. More...
 
double rayon
 Radius of the horizon in LORENE's units. More...
 
double omega
 Angular velocity in LORENE's units. More...
 
double omega_local
 local angular velocity More...
 
int rot_state
 State of rotation. More...
 
double lapse_hori
 
double * boost
 Lapse on the horizon. More...
 
double regul
 Intensity of the correction on the shift vector. More...
 
Tenseur n_auto
 Part of N generated by the hole. More...
 
Tenseur n_comp
 Part of N generated by the companion hole. More...
 
Tenseur n_tot
 Total N . More...
 
Tenseur psi_auto
 Part of $\Psi$ generated by the hole. More...
 
Tenseur psi_comp
 Part of $\Psi$ generated by the companion hole. More...
 
Tenseur psi_tot
 Total $\Psi$. More...
 
Tenseur grad_n_tot
 Total gradient of N . More...
 
Tenseur grad_psi_tot
 Total gradient of $\Psi$. More...
 
Tenseur shift_auto
 Part of $\beta^i$ generated by the hole. More...
 
Tenseur_sym taij_auto
 Part of $A^{ij}$ generated by the hole. More...
 
Tenseur_sym taij_comp
 Part of $A^{ij}$ generated by the companion hole. More...
 
Tenseur_sym taij_tot
 Total $A^{ij}$, which must be zero on the horizon of the regularisation on the shift has been done. More...
 
Tenseur_sym tkij_auto
 Auto $K^{ij}$. More...
 
Tenseur_sym tkij_tot
 Total $K^{ij}$. More...
 
Cmp decouple
 Function used to construct the part of $K^{ij}$ generated by the hole from the total $K^{ij}$. More...
 

Friends

class Bhole_binaire
 Binary black hole system. More...
 
class Bin_ns_bh
 Binary NS-BH. More...
 

Detailed Description

Black hole.

()

This class represents a black hole in a comformaly flat approximation. It is defined with an affine mapping with a nucleus, not used in calculations, and must have, at least two shells. The horizon (i.e. the surface where N=0 ) is located at the inner boudary of the innermost shell (i.e. the domain number 1).

It can described either :

  • an isolated black hole . The fields of the companion and the total fields being zero.
  • a black hole in a binary system . The boost having no direct meaning in this case.

The tensor $K^{ij}$ is the the extrinsic curavature tensor with the conformal factor extracted.

The tensor $A^{ij}$ denotes $\nabla^i N^j + \nabla^j N^i-\frac{2}{3} \nabla_k N^k \delta^{ij}$, that is $2NK^{ij}$.

Definition at line 268 of file bhole.h.

Constructor & Destructor Documentation

◆ Bhole() [1/3]

Lorene::Bhole::Bhole ( Map_af mapping)

Standard constructor.

All the fields are set to zero, except rayon , which value is deducted from the mapping

Definition at line 180 of file bhole.C.

References boost.

◆ Bhole() [2/3]

Lorene::Bhole::Bhole ( const Bhole source)

Constructor by copy.

Definition at line 200 of file bhole.C.

References boost.

◆ Bhole() [3/3]

◆ ~Bhole()

Lorene::Bhole::~Bhole ( )

Destructor.

Definition at line 216 of file bhole.C.

References boost.

Member Function Documentation

◆ area()

double Lorene::Bhole::area ( ) const

Computes the area of the throat.

Definition at line 548 of file bhole.C.

References Lorene::Map_af::integrale_surface(), mp, Lorene::pow(), psi_tot, Lorene::Cmp::raccord(), rayon, and Lorene::Cmp::std_base_scal().

◆ fait_n_comp() [1/2]

void Lorene::Bhole::fait_n_comp ( const Bhole comp)

Imports the part of N due to the companion hole comp .

The total N is then calculated.

It also imports the gradient of N and construct the total $\nabla N$.

Definition at line 257 of file bhole.C.

References Lorene::Tenseur::dec2_dzpuis(), Lorene::Map::get_bvect_cart(), grad_n_tot, Lorene::Tenseur::gradient(), Lorene::Cmp::import_symy(), mp, n_auto, n_comp, n_tot, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), and Lorene::Tenseur::set_std_base().

◆ fait_n_comp() [2/2]

void Lorene::Bhole::fait_n_comp ( const Et_bin_nsbh comp)


Imports the part of N due to the companion ns comp .

The total N is then calculated.

It also imports the gradient of N and construct the total $\nabla N$.

Definition at line 310 of file bhole.C.

References Lorene::Tenseur::dec2_dzpuis(), Lorene::Map::get_bvect_cart(), Lorene::Et_bin_nsbh::get_d_n_auto(), Lorene::Et_bin_nsbh::get_n_auto(), grad_n_tot, Lorene::Tenseur::gradient(), Lorene::Cmp::import(), Lorene::Cmp::import_symy(), mp, n_auto, n_comp, n_tot, regul, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), and Lorene::Tenseur::set_std_base().

◆ fait_psi_comp() [1/2]

void Lorene::Bhole::fait_psi_comp ( const Bhole comp)

Imports the part of $\Psi$ due to the companion hole comp .

The total $\Psi$ is then calculated.

It also imports the gradient of $\Psi$ and construct the total $\nabla \Psi$.

Definition at line 283 of file bhole.C.

References Lorene::Tenseur::dec2_dzpuis(), Lorene::Map::get_bvect_cart(), grad_psi_tot, Lorene::Tenseur::gradient(), Lorene::Cmp::import_symy(), mp, psi_auto, psi_comp, psi_tot, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), and Lorene::Tenseur::set_std_base().

◆ fait_psi_comp() [2/2]

void Lorene::Bhole::fait_psi_comp ( const Et_bin_nsbh comp)

◆ fait_taij_auto()

void Lorene::Bhole::fait_taij_auto ( )

Calculates the part of $A^{ij}$ generated by shift_auto .

Definition at line 382 of file bhole.C.

References Lorene::Tenseur::get_etat(), and shift_auto.

◆ fait_tkij()

void Lorene::Bhole::fait_tkij ( )

Calculates the total $K^{ij}$.

The regularisation of the shift must be done before to ensure regularity.

Definition at line 329 of file bhole_pseudo_kerr.C.

References fait_taij_auto(), mp, n_auto, Lorene::Cmp::raccord(), Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), taij_auto, and tkij_auto.

◆ get_boost()

double* Lorene::Bhole::get_boost ( ) const
inline

Returns the cartesian components of the boost with respect to the reference frame.

Definition at line 382 of file bhole.h.

References boost.

◆ get_decouple()

const Cmp Lorene::Bhole::get_decouple ( ) const
inline

Returns the function used to construct tkij_auto from tkij_tot .

Definition at line 468 of file bhole.h.

References decouple.

◆ get_grad_n_tot()

const Tenseur& Lorene::Bhole::get_grad_n_tot ( ) const
inline

Returns the gradient of N .

Definition at line 434 of file bhole.h.

References grad_n_tot.

◆ get_grad_psi_tot()

const Tenseur& Lorene::Bhole::get_grad_psi_tot ( ) const
inline

Returns the gradient of $\Psi$.

Definition at line 429 of file bhole.h.

References grad_psi_tot.

◆ get_mp()

const Map_af& Lorene::Bhole::get_mp ( ) const
inline

Returns the mapping (readonly).

Definition at line 339 of file bhole.h.

References mp.

◆ get_n_auto()

const Tenseur& Lorene::Bhole::get_n_auto ( ) const
inline

Returns the part of N generated by the hole.

Definition at line 395 of file bhole.h.

References n_auto.

◆ get_n_comp()

const Tenseur& Lorene::Bhole::get_n_comp ( ) const
inline

Returns the part of N generated by the companion hole.

Definition at line 403 of file bhole.h.

References n_comp.

◆ get_n_tot()

const Tenseur& Lorene::Bhole::get_n_tot ( ) const
inline

Returns the total N .

Definition at line 407 of file bhole.h.

References n_tot.

◆ get_omega()

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

Returns the angular velocity.

Definition at line 357 of file bhole.h.

References omega.

◆ get_omega_local()

double Lorene::Bhole::get_omega_local ( ) const
inline

Returns the local angular velocity.

Definition at line 365 of file bhole.h.

References omega_local.

◆ get_psi_auto()

const Tenseur& Lorene::Bhole::get_psi_auto ( ) const
inline

Returns the part of $\Psi$ generated by the hole.

Definition at line 412 of file bhole.h.

References psi_auto.

◆ get_psi_comp()

const Tenseur& Lorene::Bhole::get_psi_comp ( ) const
inline

Returns the part of $\Psi$ generated by the companion hole.

Definition at line 420 of file bhole.h.

References psi_comp.

◆ get_psi_tot()

const Tenseur& Lorene::Bhole::get_psi_tot ( ) const
inline

Returns the total $\Psi$.

Definition at line 424 of file bhole.h.

References psi_tot.

◆ get_rayon()

double Lorene::Bhole::get_rayon ( ) const
inline

Returns the radius of the horizon.

Definition at line 347 of file bhole.h.

References rayon.

◆ get_regul()

double Lorene::Bhole::get_regul ( ) const
inline

Returns the intensity of the regularisation on $\vec{N}$.

Definition at line 390 of file bhole.h.

References regul.

◆ get_rot_state()

double Lorene::Bhole::get_rot_state ( ) const
inline

Returns the state of rotation.

Definition at line 373 of file bhole.h.

References rot_state.

◆ get_shift_auto()

const Tenseur& Lorene::Bhole::get_shift_auto ( ) const
inline

Returns the part of $\beta^i$ generated by the hole.

Definition at line 440 of file bhole.h.

References shift_auto.

◆ get_taij_auto()

const Tenseur& Lorene::Bhole::get_taij_auto ( ) const
inline

Returns the part of $A^{ij}$ generated by the hole.

Definition at line 447 of file bhole.h.

References taij_auto.

◆ get_taij_comp()

const Tenseur& Lorene::Bhole::get_taij_comp ( ) const
inline

Returns the part of $A^{ij}$ generated by the companion hole.

Definition at line 451 of file bhole.h.

References taij_comp.

◆ get_taij_tot()

const Tenseur& Lorene::Bhole::get_taij_tot ( ) const
inline

Returns the total $A^{ij}$.

Definition at line 455 of file bhole.h.

References taij_tot.

◆ get_tkij_auto()

const Tenseur& Lorene::Bhole::get_tkij_auto ( ) const
inline

Returns the part of $K^{ij}$ generated by the hole.

Definition at line 464 of file bhole.h.

References tkij_auto.

◆ get_tkij_tot()

const Tenseur& Lorene::Bhole::get_tkij_tot ( ) const
inline

Returns the total $K^{ij}$.

Definition at line 460 of file bhole.h.

References tkij_tot.

◆ init_bhole()

void Lorene::Bhole::init_bhole ( )

◆ init_bhole_phi()

◆ init_bhole_seul()

void Lorene::Bhole::init_bhole_seul ( )

◆ init_kerr()

void Lorene::Bhole::init_kerr ( double  masse,
double  moment 
)

Set the inital values to those of Kerr.

Parameters
masse[input] : ADM mass in LORENE's units.
moment[input] : $\frac{J}{M}=a$, in LORENE's units.

The radius $r_0$ of *this is supposed to be equal to $\frac{1}{2}\sqrt{M^2-a^2}$.

The fields have then the following values~:

  • $\Omega = \frac{a}{2M\sqrt{M^2-a^2}}$
  • $\Psi^4 = 1+\frac{2M}{r}+\frac{3M^2+a^2\cos^2\theta}{2r^2}+\frac{2Mr_0^2}{r^3}+ \frac{r_0^4}{r^4}$.
  • $N = \left(1-\frac{2MR}{\Sigma}+\frac{4a^2M^2R^2\sin^2\theta} {\Sigma^2\left(R^2+a^2\right)+2a^2\Sigma R\sin^2\theta}\right)^{\frac{1}{2}}$
  • $N^\phi = \frac{2aMR}{\Sigma\left(R^2+a^2\right)+2a^2MR\sin^2\theta}$ {itemize} where~:
  • $R = r+\frac{M^2-a^2}{4r}+M$.
  • $\Sigma = R^2+a^2-2MR$
  • $r_0$ is equal to rayon .

Definition at line 64 of file bhole_init_kerr.C.

References Lorene::Cmp::annule(), Lorene::Map::cost, Lorene::Map::get_mg(), Lorene::Tenseur::inc_dzpuis(), mp, Lorene::Valeur::mult_cp(), Lorene::Cmp::mult_r(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), n_auto, omega, Lorene::pow(), psi_auto, Lorene::Map::r, Lorene::Cmp::raccord(), rayon, Lorene::Tenseur::set(), Lorene::Cmp::set_dzpuis(), Lorene::Tenseur::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Tenseur::set_std_base(), Lorene::Cmp::set_val_hor(), Lorene::Cmp::set_val_inf(), shift_auto, Lorene::Map::sint, Lorene::sqrt(), Lorene::Cmp::std_base_scal(), and Lorene::Cmp::va.

◆ masse_adm_seul()

double Lorene::Bhole::masse_adm_seul ( ) const

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

WARNING this should only be used for an isolated black hole.

Definition at line 347 of file bhole_pseudo_kerr.C.

References Lorene::Map_af::integrale_surface_infini(), mp, and psi_auto.

◆ masse_komar_seul()

double Lorene::Bhole::masse_komar_seul ( ) const

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

WARNING this should only be used for an isolated black hole.

Definition at line 355 of file bhole_pseudo_kerr.C.

References Lorene::Map_af::integrale_surface_infini(), mp, and n_auto.

◆ moment_seul_hor()

double Lorene::Bhole::moment_seul_hor ( ) const

Calculates the angular momentum of the black hole using the formula on the horizon : $ J = -\frac{1}{8 \pi} \oint_{H} \Psi^6 K_j^i m^j {\mathrm d}S_i$ where $\vec{m} = \frac{\partial}{\partial \phi}$ and H denotes the horizon.

WARNING It supposes that the boost is zero and should only be used for an isolated black hole..

Definition at line 413 of file bhole_pseudo_kerr.C.

References boost, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map_af::integrale_surface(), mp, omega, Lorene::pow(), psi_auto, rayon, Lorene::Tenseur::set_etat_qcq(), Lorene::Cmp::std_base_scal(), tkij_auto, Lorene::Map::xa, and Lorene::Map::ya.

◆ moment_seul_inf()

double Lorene::Bhole::moment_seul_inf ( ) const

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}$.

WARNING It supposes that the boost is zero and should only be used for an isolated black hole..

Definition at line 366 of file bhole_pseudo_kerr.C.

References boost, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map_af::integrale_surface_infini(), mp, Lorene::Valeur::mult_cp(), Lorene::Cmp::mult_r_zec(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), omega, Lorene::Tenseur::set_etat_qcq(), tkij_auto, and Lorene::Cmp::va.

◆ operator=()

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

◆ regularise_seul()

void Lorene::Bhole::regularise_seul ( )

Corrects the shift in the innermost shell, so that it remains $ {\mathcal{C}}^2$ and that $A^{ij}$ equals zero on the horizon.

regul is then, the relative difference between the shift before and after the regularisation.

WARNING this should only be used for an isolated black hole.

Definition at line 257 of file bhole_pseudo_kerr.C.

References Lorene::Cmp::annule(), boost, Lorene::Tenseur::c, Lorene::Valeur::coef_i(), Lorene::diffrelmax(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Tenseur::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), mp, Lorene::norme(), omega, Lorene::pow(), Lorene::Map::r, regul, Lorene::Tenseur::set(), Lorene::Valeur::set_etat_c_qcq(), Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_std_base(), shift_auto, Lorene::Cmp::va, Lorene::Map_af::val_r(), Lorene::Map::x, and Lorene::Map::y.

◆ regularise_shift()

void Lorene::Bhole::regularise_shift ( const Tenseur comp)

Corrects shift_auto in such a way that the total $A^{ij}$ is equal to zero in the horizon, which should ensure the regularity of $K^{ij}$, using the stored values of the boost and the angular velocity.

WARNING : this should only be used for a black hole in a binary system Bhole_binaire or Bin_ns_bh .

Parameters
comp[input] : the part of $\beta^i$ generated by the companion hole.

Definition at line 406 of file bhole.C.

References omega, omega_local, regul, and shift_auto.

◆ return_n_auto()

Cmp Lorene::Bhole::return_n_auto ( ) const
inline

Returns the part of N generated by the hole as a cmp.

Definition at line 399 of file bhole.h.

References n_auto.

◆ return_psi_auto()

Cmp Lorene::Bhole::return_psi_auto ( ) const
inline

Returns the part of $\Psi$ generated by the hole as a cmp.

Definition at line 416 of file bhole.h.

References psi_auto.

◆ sauve()

void Lorene::Bhole::sauve ( FILE *  fich) const

Write on a file.

Definition at line 485 of file bhole.C.

References boost, Lorene::fwrite_be(), n_auto, omega, omega_local, psi_auto, regul, rot_state, Lorene::Tenseur::sauve(), and shift_auto.

◆ set_mp()

Map_af& Lorene::Bhole::set_mp ( )
inline

Read/write of the mapping.

Definition at line 342 of file bhole.h.

References mp.

◆ set_omega()

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

Sets the angular velocity to ome .

Definition at line 361 of file bhole.h.

References omega.

◆ set_omega_local()

void Lorene::Bhole::set_omega_local ( double  ome)
inline

Sets the local angular velocity to ome .

Definition at line 369 of file bhole.h.

References omega_local.

◆ set_rayon()

void Lorene::Bhole::set_rayon ( double  ray)
inline

Sets the radius of the horizon to ray .

Definition at line 352 of file bhole.h.

References rayon.

◆ set_rot_state()

void Lorene::Bhole::set_rot_state ( int  rotation)
inline

Returns the state of rotation.

Definition at line 377 of file bhole.h.

References rot_state.

◆ solve_lapse_seul()

void Lorene::Bhole::solve_lapse_seul ( double  relax)


Solves the equation for N ~:

\[ *\Delta N = -\frac{2}{\Psi}\nabla_i \Psi \nabla^i N + N \Psi^4 K_{ij} K^{ij} \]

with the condition that N =0 on the horizon.

Parameters
relax[input] : the relaxation parameter.

WARNING this should only be used for an isolated black hole.

Definition at line 100 of file bhole_pseudo_kerr.C.

References Lorene::flat_scalar_prod(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Tenseur::gradient(), mp, n_auto, Lorene::pow(), psi_auto, Lorene::Cmp::raccord(), Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Valeur::std_base_scal(), Lorene::Cmp::std_base_scal(), and tkij_auto.

◆ solve_lapse_with_ns()

void Lorene::Bhole::solve_lapse_with_ns ( double  relax,
int  bound_nn,
double  lim_nn 
)

◆ solve_psi_seul()

void Lorene::Bhole::solve_psi_seul ( double  relax)

Solves the equation for $\Psi$~:

\[ \Delta \Psi = - \frac{\Psi^5} {8} K_{ij}K^{ij} \]

with the condition that $\partial_r \Psi =-\frac{1}{2 {\tt rayon}} f\left(\theta, \phi\right)$ on the horizon, where f is the value of $\Psi$ on the horizon at the preceeding step.

Parameters
relax[input] : the relaxation parameter.

WARNING this should only be used for an isolated black hole.

Definition at line 138 of file bhole_pseudo_kerr.C.

References Lorene::flat_scalar_prod(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Cmp::poisson_neumann(), Lorene::pow(), psi_auto, Lorene::Cmp::raccord(), rayon, Lorene::Valeur::set(), Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Cmp::std_base_scal(), and tkij_auto.

◆ solve_psi_with_ns()

void Lorene::Bhole::solve_psi_with_ns ( double  relax)


Solves the equation for $\Psi$~:

\[ (\Delta \Psi = - \frac{\Psi^5} {8} K_{ij}K^{ij} \]

with the condition that $\partial_r \Psi=-\frac{1}{2 {\tt rayon}} f\left(\theta, \phi\right)$ on the horizon, where f is the value of $\Psi$ on the horizon at the preceeding step.

Parameters
relax[input] : the relaxation parameter.

WARNING this should only be used for BH in a Bin_ns_bh .

Definition at line 156 of file bhole_with_ns.C.

References Lorene::cos(), Lorene::flat_scalar_prod(), Lorene::Mg3d::get_angu(), Lorene::Tenseur::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Map::phi, Lorene::pow(), psi_auto, psi_tot, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::sin(), Lorene::Cmp::std_base_scal(), Lorene::Map::tet, tkij_auto, and tkij_tot.

◆ solve_shift_seul()

void Lorene::Bhole::solve_shift_seul ( double  precis,
double  relax 
)


Solves the equation for $\vec{\beta}$~:

\[ *\Delta \beta^i +\frac{1}{3} \nabla^i \left(\nabla_j \beta^j\right) = -6A^{ij}\frac{\nabla_j \Psi}{\Psi} + 2 K^{ij}\nabla_j N \]

with $\vec{\beta} = -\Omega \vec{m} - \vec{V}$ on the horizon, $\vec{V}$ being the boost and $\vec{m} = \frac{\partial}{\partial \phi}$.

The solution is solved using Oohara-scheme and an iteration.

Parameters
precis[input] : parameter for the Oohara-solver which is an iterative scheme.
relax[input] : the relaxation parameter.

WARNING this should only be used for an isolated black hole.

Definition at line 179 of file bhole_pseudo_kerr.C.

References Lorene::flat_scalar_prod(), Lorene::Tenseur::get_etat(), Lorene::Tenseur::gradient(), n_auto, psi_auto, Lorene::Tenseur::set_std_base(), shift_auto, taij_auto, and tkij_auto.

◆ solve_shift_with_ns()

void Lorene::Bhole::solve_shift_with_ns ( const Et_bin_nsbh ns,
double  precis,
double  relax,
int  bound_nn,
double  lim_nn 
)


Solves the equation for $\vec{\beta}$~:

\[ *\Delta \beta^i +\frac{1}{3} \nabla^i \left(\nabla_j \beta^j\right) = -6A^{ij}\frac{\nabla_j \Psi}{\Psi} + 2 K^{ij}\nabla_j N \]

with $\vec{\beta} = -\Omega \vec{m} - \vec{V}$ on the horizon, $\vec{V}$ being the boost and $\vec{m} = \frac{\partial}{\partial \phi}$.

The solution is solved using Oohara-scheme and an iteration.

Parameters
ns[input] : the companion.
precis[input] : parameter for the Oohara-solver which is an iterative scheme.
relax[input] : the relaxation parameter.

WARNING this should only be used for BH in a Bin_ns_bh .

Definition at line 225 of file bhole_with_ns.C.

References Lorene::flat_scalar_prod(), Lorene::Tenseur::get_etat(), Lorene::Tenseur::gradient(), n_auto, psi_auto, psi_tot, Lorene::Tenseur::set_std_base(), shift_auto, taij_tot, and tkij_tot.

◆ viriel_seul()

double Lorene::Bhole::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 $\Psi$ and N .

WARNING this should only be used for an isolated black hole.

Definition at line 82 of file bhole_pseudo_viriel.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), mp, n_auto, and psi_auto.

Friends And Related Function Documentation

◆ Bhole_binaire

friend class Bhole_binaire
friend

Binary black hole system.

Definition at line 747 of file bhole.h.

◆ Bin_ns_bh

friend class Bin_ns_bh
friend

Binary NS-BH.

Definition at line 748 of file bhole.h.

Member Data Documentation

◆ boost

double* Lorene::Bhole::boost
protected

Lapse on the horizon.

Cartesian components of the boost in the reference frame.

Definition at line 283 of file bhole.h.

◆ decouple

Cmp Lorene::Bhole::decouple
protected

Function used to construct the part of $K^{ij}$ generated by the hole from the total $K^{ij}$.

Only used for a binary system.

Mainly this Cmp is 1 around the hole and 0 around the companion and the sum of decouple for the hole and his companion is 1 everywhere.

Definition at line 318 of file bhole.h.

◆ grad_n_tot

Tenseur Lorene::Bhole::grad_n_tot
protected

Total gradient of N .

Definition at line 294 of file bhole.h.

◆ grad_psi_tot

Tenseur Lorene::Bhole::grad_psi_tot
protected

Total gradient of $\Psi$.

Definition at line 295 of file bhole.h.

◆ mp

Map_af& Lorene::Bhole::mp
protected

Affine mapping.

Definition at line 273 of file bhole.h.

◆ n_auto

Tenseur Lorene::Bhole::n_auto
protected

Part of N generated by the hole.

Definition at line 286 of file bhole.h.

◆ n_comp

Tenseur Lorene::Bhole::n_comp
protected

Part of N generated by the companion hole.

Definition at line 287 of file bhole.h.

◆ n_tot

Tenseur Lorene::Bhole::n_tot
protected

Total N .

Definition at line 288 of file bhole.h.

◆ omega

double Lorene::Bhole::omega
protected

Angular velocity in LORENE's units.

Definition at line 275 of file bhole.h.

◆ omega_local

double Lorene::Bhole::omega_local
protected

local angular velocity

Definition at line 276 of file bhole.h.

◆ psi_auto

Tenseur Lorene::Bhole::psi_auto
protected

Part of $\Psi$ generated by the hole.

Definition at line 290 of file bhole.h.

◆ psi_comp

Tenseur Lorene::Bhole::psi_comp
protected

Part of $\Psi$ generated by the companion hole.

Definition at line 291 of file bhole.h.

◆ psi_tot

Tenseur Lorene::Bhole::psi_tot
protected

Total $\Psi$.

Definition at line 292 of file bhole.h.

◆ rayon

double Lorene::Bhole::rayon
protected

Radius of the horizon in LORENE's units.

Definition at line 274 of file bhole.h.

◆ regul

double Lorene::Bhole::regul
protected

Intensity of the correction on the shift vector.

Definition at line 284 of file bhole.h.

◆ rot_state

int Lorene::Bhole::rot_state
protected

State of rotation.

Definition at line 277 of file bhole.h.

◆ shift_auto

Tenseur Lorene::Bhole::shift_auto
protected

Part of $\beta^i$ generated by the hole.

Definition at line 297 of file bhole.h.

◆ taij_auto

Tenseur_sym Lorene::Bhole::taij_auto
protected

Part of $A^{ij}$ generated by the hole.

Definition at line 299 of file bhole.h.

◆ taij_comp

Tenseur_sym Lorene::Bhole::taij_comp
protected

Part of $A^{ij}$ generated by the companion hole.

Definition at line 300 of file bhole.h.

◆ taij_tot

Tenseur_sym Lorene::Bhole::taij_tot
protected

Total $A^{ij}$, which must be zero on the horizon of the regularisation on the shift has been done.

Definition at line 305 of file bhole.h.

◆ tkij_auto

Tenseur_sym Lorene::Bhole::tkij_auto
protected

Auto $K^{ij}$.

Definition at line 307 of file bhole.h.

◆ tkij_tot

Tenseur_sym Lorene::Bhole::tkij_tot
protected

Total $K^{ij}$.

Definition at line 308 of file bhole.h.


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