LORENE
|
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_af & | get_mp () const |
Returns the mapping (readonly). More... | |
Map_af & | set_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 . More... | |
const Tenseur & | get_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 Tenseur & | get_n_comp () const |
Returns the part of N generated by the companion hole. More... | |
const Tenseur & | get_n_tot () const |
Returns the total N . More... | |
const Tenseur & | get_psi_auto () const |
Returns the part of generated by the hole. More... | |
Cmp | return_psi_auto () const |
Returns the part of generated by the hole as a cmp. More... | |
const Tenseur & | get_psi_comp () const |
Returns the part of generated by the companion hole. More... | |
const Tenseur & | get_psi_tot () const |
Returns the total . More... | |
const Tenseur & | get_grad_psi_tot () const |
Returns the gradient of . More... | |
const Tenseur & | get_grad_n_tot () const |
Returns the gradient of N . More... | |
const Tenseur & | get_shift_auto () const |
Returns the part of generated by the hole. More... | |
Cmp | return_shift_auto (int i) const |
const Tenseur & | get_taij_auto () const |
Returns the part of generated by the hole. More... | |
const Tenseur & | get_taij_comp () const |
Returns the part of generated by the companion hole. More... | |
const Tenseur & | get_taij_tot () const |
Returns the total . More... | |
const Tenseur & | get_tkij_tot () const |
Returns the total . More... | |
const Tenseur & | get_tkij_auto () const |
Returns the part of generated by the hole. More... | |
const Cmp | get_decouple () const |
Returns the function used to construct tkij_auto from tkij_tot . More... | |
Cmp & | set_n_auto () |
Cmp & | set_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 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 due to the companion ns comp . More... | |
void | fait_taij_auto () |
Calculates the part of generated by shift_auto . More... | |
void | regularise_shift (const Tenseur &comp) |
Corrects shift_auto in such a way that the total is equal to zero in the horizon, which should ensure the regularity of , 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 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 ~:
with the condition that | |
void | solve_psi_seul (double relax) |
Solves the equation for ~:
with the condition that on the horizon, where | |
void | solve_shift_seul (double precis, double relax) |
Solves the equation for ~:
with on the horizon, being the boost and . More... | |
void | regularise_seul () |
Corrects the shift in the innermost shell, so that it remains and that equals zero on the horizon. More... | |
void | solve_lapse_with_ns (double relax, int bound_nn, double lim_nn) |
Solves the equation for N ~:
with the condition that | |
void | solve_psi_with_ns (double relax) |
Solves the equation for ~:
with the condition that on the horizon, where | |
void | solve_shift_with_ns (const Et_bin_nsbh &ns, double precis, double relax, int bound_nn, double lim_nn) |
Solves the equation for ~:
with on the horizon, being the boost and . 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 . 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 : . More... | |
double | masse_komar_seul () const |
Calculates the Komar mass of the black hole using : . More... | |
double | moment_seul_inf () const |
Calculates the angular momentum of the black hole using the formula at infinity : where . More... | |
double | moment_seul_hor () const |
Calculates the angular momentum of the black hole using the formula on the horizon : where and H denotes the horizon. More... | |
double | local_momentum () const |
void | init_bhole_phi () |
Initiates the black hole for a resolution with . More... | |
void | init_bhole_seul () |
Initiates for a single the black hole. More... | |
Protected Attributes | |
Map_af & | mp |
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 generated by the hole. More... | |
Tenseur | psi_comp |
Part of generated by the companion hole. More... | |
Tenseur | psi_tot |
Total . More... | |
Tenseur | grad_n_tot |
Total gradient of N . More... | |
Tenseur | grad_psi_tot |
Total gradient of . More... | |
Tenseur | shift_auto |
Part of generated by the hole. More... | |
Tenseur_sym | taij_auto |
Part of generated by the hole. More... | |
Tenseur_sym | taij_comp |
Part of generated by the companion hole. More... | |
Tenseur_sym | taij_tot |
Total , which must be zero on the horizon of the regularisation on the shift has been done. More... | |
Tenseur_sym | tkij_auto |
Auto . More... | |
Tenseur_sym | tkij_tot |
Total . More... | |
Cmp | decouple |
Function used to construct the part of generated by the hole from the total . More... | |
Friends | |
class | Bhole_binaire |
Binary black hole system. More... | |
class | Bin_ns_bh |
Binary NS-BH. More... | |
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 is the the extrinsic curavature tensor with the conformal factor extracted.
The tensor denotes , that is .
Lorene::Bhole::Bhole | ( | Map_af & | mapping | ) |
Lorene::Bhole::Bhole | ( | const Bhole & | source | ) |
Lorene::Bhole::Bhole | ( | Map_af & | mpi, |
FILE * | fich, | ||
bool | old = false |
||
) |
Constructor from a Map_af
and a file.
Definition at line 499 of file bhole.C.
References boost, decouple, Lorene::fread_be(), Lorene::Map::get_bvect_cart(), grad_n_tot, grad_psi_tot, mp, n_auto, omega, omega_local, psi_auto, regul, rot_state, Lorene::Tenseur::set_etat_qcq(), Lorene::Tenseur::set_etat_zero(), Lorene::Cmp::set_etat_zero(), shift_auto, taij_auto, taij_comp, taij_tot, tkij_auto, and tkij_tot.
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().
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 .
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().
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 .
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().
void Lorene::Bhole::fait_psi_comp | ( | const Bhole & | comp | ) |
Imports the part of due to the companion hole comp
.
The total is then calculated.
It also imports the gradient of and construct the total .
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().
void Lorene::Bhole::fait_psi_comp | ( | const Et_bin_nsbh & | comp | ) |
Imports the part of due to the companion ns comp
.
The total is then calculated.
It also imports the gradient of and construct the total .
Definition at line 347 of file bhole.C.
References Lorene::Tenseur::dec2_dzpuis(), Lorene::Map::get_bvect_cart(), Lorene::Et_bin_nsbh::get_confpsi_auto(), Lorene::Et_bin_nsbh::get_d_confpsi_auto(), grad_psi_tot, Lorene::Tenseur::gradient(), Lorene::Cmp::import(), Lorene::Cmp::import_symy(), mp, psi_auto, psi_comp, psi_tot, regul, Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), and Lorene::Tenseur::set_std_base().
void Lorene::Bhole::fait_taij_auto | ( | ) |
Calculates the part of generated by shift_auto
.
Definition at line 382 of file bhole.C.
References Lorene::Tenseur::get_etat(), and shift_auto.
void Lorene::Bhole::fait_tkij | ( | ) |
Calculates the total .
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.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Returns the part of generated by the hole.
Definition at line 440 of file bhole.h.
References shift_auto.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
void Lorene::Bhole::init_bhole | ( | ) |
Sets the values of the fields to :
n_auto
n_comp
psi_auto
psi_comp
a being the radius of the hole, the other fields being set to zero.
Definition at line 412 of file bhole.C.
References Lorene::Cmp::annule(), decouple, grad_n_tot, grad_psi_tot, Lorene::Tenseur::gradient(), mp, n_auto, n_comp, n_tot, psi_auto, psi_comp, psi_tot, Lorene::Map::r, Lorene::Cmp::raccord(), rayon, Lorene::Tenseur::set(), Lorene::Cmp::set_dzpuis(), Lorene::Tenseur::set_etat_zero(), Lorene::Cmp::set_etat_zero(), Lorene::Tenseur::set_std_base(), shift_auto, taij_auto, taij_comp, taij_tot, tkij_auto, and tkij_tot.
void Lorene::Bhole::init_bhole_phi | ( | ) |
Initiates the black hole for a resolution with .
Definition at line 78 of file bhole_solve_phi.C.
References Lorene::Cmp::annule(), decouple, grad_n_tot, grad_psi_tot, Lorene::Tenseur::gradient(), Lorene::log(), mp, n_auto, n_comp, n_tot, psi_auto, psi_comp, psi_tot, Lorene::Map::r, Lorene::Cmp::raccord(), rayon, Lorene::Tenseur::set(), Lorene::Cmp::set_dzpuis(), Lorene::Tenseur::set_etat_zero(), Lorene::Cmp::set_etat_zero(), Lorene::Tenseur::set_std_base(), shift_auto, taij_auto, taij_comp, taij_tot, tkij_auto, and tkij_tot.
void Lorene::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 448 of file bhole.C.
References Lorene::Cmp::annule(), decouple, grad_n_tot, grad_psi_tot, Lorene::Tenseur::gradient(), mp, n_auto, n_comp, n_tot, psi_auto, psi_comp, psi_tot, Lorene::Map::r, Lorene::Cmp::raccord(), rayon, Lorene::Tenseur::set(), Lorene::Cmp::set_dzpuis(), Lorene::Tenseur::set_etat_zero(), Lorene::Cmp::set_etat_zero(), Lorene::Tenseur::set_std_base(), Lorene::Cmp::set_val_inf(), shift_auto, taij_auto, taij_comp, taij_tot, tkij_auto, and tkij_tot.
void Lorene::Bhole::init_kerr | ( | double | masse, |
double | moment | ||
) |
Set the inital values to those of Kerr.
masse | [input] : ADM mass in LORENE's units. |
moment | [input] : , in LORENE's units. |
The radius of *this
is supposed to be equal to .
The fields have then the following values~:
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.
double Lorene::Bhole::masse_adm_seul | ( | ) | const |
Calculates the ADM mass of the black hole using : .
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.
double Lorene::Bhole::masse_komar_seul | ( | ) | const |
Calculates the Komar mass of the black hole using : .
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.
double Lorene::Bhole::moment_seul_hor | ( | ) | const |
Calculates the angular momentum of the black hole using the formula on the horizon : where 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.
double Lorene::Bhole::moment_seul_inf | ( | ) | const |
Calculates the angular momentum of the black hole using the formula at infinity : where .
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.
void Lorene::Bhole::operator= | ( | const Bhole & | source | ) |
void Lorene::Bhole::regularise_seul | ( | ) |
Corrects the shift in the innermost shell, so that it remains and that 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.
void Lorene::Bhole::regularise_shift | ( | const Tenseur & | comp | ) |
Corrects shift_auto
in such a way that the total is equal to zero in the horizon, which should ensure the regularity of , 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
.
comp | [input] : the part of generated by the companion hole. |
Definition at line 406 of file bhole.C.
References omega, omega_local, regul, and shift_auto.
|
inline |
|
inline |
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.
|
inline |
|
inline |
|
inline |
Sets the local angular velocity to ome
.
Definition at line 369 of file bhole.h.
References omega_local.
|
inline |
|
inline |
void Lorene::Bhole::solve_lapse_seul | ( | double | relax | ) |
Solves the equation for N ~:
with the condition that N
=0 on the horizon.
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.
void Lorene::Bhole::solve_lapse_with_ns | ( | double | relax, |
int | bound_nn, | ||
double | lim_nn | ||
) |
Solves the equation for N ~:
with the condition that N
=0 on the horizon.
relax | [input] : the relaxation parameter. |
WARNING this should only be used for BH in a Bin_ns_bh
.
Definition at line 92 of file bhole_with_ns.C.
References Lorene::Valeur::annule_hard(), Lorene::flat_scalar_prod(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), grad_n_tot, Lorene::Tenseur::gradient(), mp, n_auto, n_comp, n_tot, Lorene::Cmp::poisson_dirichlet(), Lorene::Cmp::poisson_neumann(), Lorene::pow(), psi_auto, psi_tot, Lorene::Cmp::raccord(), Lorene::Tenseur::set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Cmp::std_base_scal(), tkij_auto, and tkij_tot.
void Lorene::Bhole::solve_psi_seul | ( | double | relax | ) |
Solves the equation for ~:
with the condition that on the horizon, where f
is the value of on the horizon at the preceeding step.
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.
void Lorene::Bhole::solve_psi_with_ns | ( | double | relax | ) |
Solves the equation for ~:
with the condition that on the horizon, where f
is the value of on the horizon at the preceeding step.
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.
void Lorene::Bhole::solve_shift_seul | ( | double | precis, |
double | relax | ||
) |
Solves the equation for ~:
with on the horizon, being the boost and .
The solution is solved using Oohara-scheme and an iteration.
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.
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 ~:
with on the horizon, being the boost and .
The solution is solved using Oohara-scheme and an iteration.
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.
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 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.
|
friend |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |