Bhole Class Reference
[Stars and black holes]

Black hole. More...

#include <bhole.h>

List of all members.

Public Member Functions

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

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.

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

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

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.

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

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}$.
double area () const
 Computes the area of the throat.
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$.
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$.
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}$.
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.
double local_momentum () const
void init_bhole_phi ()
 Initiates the black hole for a resolution with $\Phi = \log \Psi$.
void init_bhole_seul ()
 Initiates for a single the black hole.

Protected Attributes

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

Friends

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

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 :

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 263 of file bhole.h.


Constructor & Destructor Documentation

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 173 of file bhole.C.

References boost.

Bhole::Bhole ( const Bhole source  ) 

Constructor by copy.

Definition at line 193 of file bhole.C.

References boost.

Bhole::Bhole ( Map_af mpi,
FILE *  fich,
bool  old = false 
)
Bhole::~Bhole (  ) 

Destructor.

Definition at line 209 of file bhole.C.

References boost.


Member Function Documentation

double Bhole::area (  )  const

Computes the area of the throat.

Definition at line 541 of file bhole.C.

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

void 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 303 of file bhole.C.

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

void 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 250 of file bhole.C.

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

void Bhole::fait_psi_comp ( const Et_bin_nsbh comp  ) 

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

The total $\Psi$ is then calculated.

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

Definition at line 340 of file bhole.C.

References Map::get_bvect_cart(), Et_bin_nsbh::get_confpsi_auto(), Et_bin_nsbh::get_d_confpsi_auto(), grad_psi_tot, Tenseur::gradient(), Cmp::import(), Cmp::import_symy(), mp, psi_auto, psi_comp, psi_tot, regul, Tenseur::set(), Tenseur::set_etat_qcq(), and Tenseur::set_std_base().

void 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 276 of file bhole.C.

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

void Bhole::fait_taij_auto (  ) 

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

Definition at line 375 of file bhole.C.

References Tenseur::get_etat(), Tenseur::gradient(), mp, Cmp::raccord(), Tenseur::set(), Tenseur::set_etat_qcq(), Tenseur::set_etat_zero(), Tenseur::set_std_base(), shift_auto, and taij_auto.

void Bhole::fait_tkij (  ) 

Calculates the total $K^{ij}$.

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

Definition at line 322 of file bhole_pseudo_kerr.C.

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

double* Bhole::get_boost (  )  const [inline]

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

Definition at line 377 of file bhole.h.

References boost.

const Cmp Bhole::get_decouple (  )  const [inline]

Returns the function used to construct tkij_auto from tkij_tot .

Definition at line 463 of file bhole.h.

References decouple.

const Tenseur& Bhole::get_grad_n_tot (  )  const [inline]

Returns the gradient of N .

Definition at line 429 of file bhole.h.

References grad_n_tot.

const Tenseur& Bhole::get_grad_psi_tot (  )  const [inline]

Returns the gradient of $\Psi$.

Definition at line 424 of file bhole.h.

References grad_psi_tot.

const Map_af& Bhole::get_mp (  )  const [inline]

Returns the mapping (readonly).

Definition at line 334 of file bhole.h.

References mp.

const Tenseur& Bhole::get_n_auto (  )  const [inline]

Returns the part of N generated by the hole.

Definition at line 390 of file bhole.h.

References n_auto.

const Tenseur& Bhole::get_n_comp (  )  const [inline]

Returns the part of N generated by the companion hole.

Definition at line 398 of file bhole.h.

References n_comp.

const Tenseur& Bhole::get_n_tot (  )  const [inline]

Returns the total N .

Definition at line 402 of file bhole.h.

References n_tot.

double Bhole::get_omega (  )  const [inline]

Returns the angular velocity.

Definition at line 352 of file bhole.h.

References omega.

double Bhole::get_omega_local (  )  const [inline]

Returns the local angular velocity.

Definition at line 360 of file bhole.h.

References omega_local.

const Tenseur& Bhole::get_psi_auto (  )  const [inline]

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

Definition at line 407 of file bhole.h.

References psi_auto.

const Tenseur& Bhole::get_psi_comp (  )  const [inline]

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

Definition at line 415 of file bhole.h.

References psi_comp.

const Tenseur& Bhole::get_psi_tot (  )  const [inline]

Returns the total $\Psi$.

Definition at line 419 of file bhole.h.

References psi_tot.

double Bhole::get_rayon (  )  const [inline]

Returns the radius of the horizon.

Definition at line 342 of file bhole.h.

References rayon.

double Bhole::get_regul (  )  const [inline]

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

Definition at line 385 of file bhole.h.

References regul.

double Bhole::get_rot_state (  )  const [inline]

Returns the state of rotation.

Definition at line 368 of file bhole.h.

References rot_state.

const Tenseur& Bhole::get_shift_auto (  )  const [inline]

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

Definition at line 435 of file bhole.h.

References shift_auto.

const Tenseur& Bhole::get_taij_auto (  )  const [inline]

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

Definition at line 442 of file bhole.h.

References taij_auto.

const Tenseur& Bhole::get_taij_comp (  )  const [inline]

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

Definition at line 446 of file bhole.h.

References taij_comp.

const Tenseur& Bhole::get_taij_tot (  )  const [inline]

Returns the total $A^{ij}$.

Definition at line 450 of file bhole.h.

References taij_tot.

const Tenseur& Bhole::get_tkij_auto (  )  const [inline]

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

Definition at line 459 of file bhole.h.

References tkij_auto.

const Tenseur& Bhole::get_tkij_tot (  )  const [inline]

Returns the total $K^{ij}$.

Definition at line 455 of file bhole.h.

References tkij_tot.

void Bhole::init_bhole (  ) 

Sets the values of the fields to :.

  • n_auto $= \frac{1}{2}-2\frac{a}{r}$
  • n_comp $= \frac{1}{2}$
  • psi_auto $= \frac{1}{2}+\frac{a}{r}$
  • psi_comp $= \frac{1}{2}$

a being the radius of the hole, the other fields being set to zero.

Definition at line 405 of file bhole.C.

References Cmp::annule(), decouple, grad_n_tot, grad_psi_tot, Tenseur::gradient(), mp, n_auto, n_comp, n_tot, psi_auto, psi_comp, psi_tot, Map::r, Cmp::raccord(), rayon, Tenseur::set(), Cmp::set_dzpuis(), Cmp::set_etat_zero(), Tenseur::set_etat_zero(), Tenseur::set_std_base(), shift_auto, taij_auto, taij_comp, taij_tot, tkij_auto, and tkij_tot.

void Bhole::init_bhole_phi (  ) 
void Bhole::init_bhole_seul (  ) 

Initiates for a single the black hole.

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

Definition at line 441 of file bhole.C.

References Cmp::annule(), decouple, grad_n_tot, grad_psi_tot, Tenseur::gradient(), mp, n_auto, n_comp, n_tot, psi_auto, psi_comp, psi_tot, Map::r, Cmp::raccord(), rayon, Tenseur::set(), Cmp::set_dzpuis(), Cmp::set_etat_zero(), Tenseur::set_etat_zero(), Tenseur::set_std_base(), Cmp::set_val_inf(), shift_auto, taij_auto, taij_comp, taij_tot, tkij_auto, and tkij_tot.

void 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 57 of file bhole_init_kerr.C.

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

double 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 340 of file bhole_pseudo_kerr.C.

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

double 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 348 of file bhole_pseudo_kerr.C.

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

double 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 406 of file bhole_pseudo_kerr.C.

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

double 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 359 of file bhole_pseudo_kerr.C.

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

void Bhole::operator= ( const Bhole source  ) 
void 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 250 of file bhole_pseudo_kerr.C.

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

void 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 399 of file bhole.C.

References omega, omega_local, regul, and shift_auto.

Cmp Bhole::return_n_auto (  )  const [inline]

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

Definition at line 394 of file bhole.h.

References n_auto.

Cmp Bhole::return_psi_auto (  )  const [inline]

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

Definition at line 411 of file bhole.h.

References psi_auto.

void Bhole::sauve ( FILE *  fich  )  const

Write on a file.

Definition at line 478 of file bhole.C.

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

Map_af& Bhole::set_mp (  )  [inline]

Read/write of the mapping.

Definition at line 337 of file bhole.h.

References mp.

void Bhole::set_omega ( double  ome  )  [inline]

Sets the angular velocity to ome .

Definition at line 356 of file bhole.h.

References omega.

void Bhole::set_omega_local ( double  ome  )  [inline]

Sets the local angular velocity to ome .

Definition at line 364 of file bhole.h.

References omega_local.

void Bhole::set_rayon ( double  ray  )  [inline]

Sets the radius of the horizon to ray .

Definition at line 347 of file bhole.h.

References rayon.

void Bhole::set_rot_state ( int  rotation  )  [inline]

Returns the state of rotation.

Definition at line 372 of file bhole.h.

References rot_state.

void 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 93 of file bhole_pseudo_kerr.C.

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

void Bhole::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.

Parameters:
relax [input] : the relaxation parameter.

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

Definition at line 85 of file bhole_with_ns.C.

References flat_scalar_prod(), Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), grad_n_tot, Tenseur::gradient(), mp, n_auto, n_comp, n_tot, Cmp::poisson_dirichlet(), Cmp::poisson_neumann(), pow(), psi_auto, psi_tot, Cmp::raccord(), Tenseur::set(), Tenseur::set_etat_qcq(), Cmp::std_base_scal(), tkij_auto, and tkij_tot.

void 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 131 of file bhole_pseudo_kerr.C.

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

void 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 149 of file bhole_with_ns.C.

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

void 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 172 of file bhole_pseudo_kerr.C.

References boost, Cmp::filtre(), flat_scalar_prod(), Mg3d::get_angu(), Tenseur::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Tenseur::gradient(), mp, n_auto, omega, psi_auto, Tenseur::set(), Tenseur::set_etat_zero(), Tenseur::set_std_base(), shift_auto, Mg3d::std_base_vect_cart(), taij_auto, tkij_auto, Map::x, and Map::y.

void 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 218 of file bhole_with_ns.C.

References Map::convert_absolute(), cos(), Cmp::filtre(), flat_scalar_prod(), Mg3d::get_angu(), Tenseur::get_etat(), Map::get_mg(), Etoile::get_mp(), Mg3d::get_np(), Mg3d::get_nt(), Etoile_bin::get_shift_auto(), Tenseur::gradient(), mp, n_auto, n_tot, omega, omega_local, Map::phi, psi_auto, psi_tot, Cmp::raccord(), regul, Tenseur::set(), Tenseur::set_etat_zero(), Tenseur::set_std_base(), shift_auto, sin(), Mg3d::std_base_vect_cart(), taij_tot, Map::tet, tkij_tot, Map::x, Map::xa, Map::y, Map::ya, and Map::za.

double 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 75 of file bhole_pseudo_viriel.C.

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


Friends And Related Function Documentation

friend class Bhole_binaire [friend]

Binary black hole system.

Definition at line 742 of file bhole.h.

friend class Bin_ns_bh [friend]

Binary NS-BH.

Definition at line 743 of file bhole.h.


Member Data Documentation

double* Bhole::boost [protected]

Lapse on the horizon.

Cartesian components of the boost in the reference frame.

Definition at line 278 of file bhole.h.

Cmp 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 313 of file bhole.h.

Total gradient of N .

Definition at line 289 of file bhole.h.

Total gradient of $\Psi$.

Definition at line 290 of file bhole.h.

Map_af& Bhole::mp [protected]

Affine mapping.

Definition at line 268 of file bhole.h.

Tenseur Bhole::n_auto [protected]

Part of N generated by the hole.

Definition at line 281 of file bhole.h.

Tenseur Bhole::n_comp [protected]

Part of N generated by the companion hole.

Definition at line 282 of file bhole.h.

Tenseur Bhole::n_tot [protected]

Total N .

Definition at line 283 of file bhole.h.

double Bhole::omega [protected]

Angular velocity in LORENE's units.

Definition at line 270 of file bhole.h.

double Bhole::omega_local [protected]

local angular velocity

Definition at line 271 of file bhole.h.

Tenseur Bhole::psi_auto [protected]

Part of $\Psi$ generated by the hole.

Definition at line 285 of file bhole.h.

Tenseur Bhole::psi_comp [protected]

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

Definition at line 286 of file bhole.h.

Tenseur Bhole::psi_tot [protected]

Total $\Psi$.

Definition at line 287 of file bhole.h.

double Bhole::rayon [protected]

Radius of the horizon in LORENE's units.

Definition at line 269 of file bhole.h.

double Bhole::regul [protected]

Intensity of the correction on the shift vector.

Definition at line 279 of file bhole.h.

int Bhole::rot_state [protected]

State of rotation.

Definition at line 272 of file bhole.h.

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

Definition at line 292 of file bhole.h.

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

Definition at line 294 of file bhole.h.

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

Definition at line 295 of file bhole.h.

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

Definition at line 300 of file bhole.h.

Auto $K^{ij}$.

Definition at line 302 of file bhole.h.

Total $K^{ij}$.

Definition at line 303 of file bhole.h.


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

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1