Black hole. More...
#include <bhole.h>
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_af & | get_mp () const |
Returns the mapping (readonly). | |
Map_af & | set_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 ![]() | |
const Tenseur & | get_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 Tenseur & | get_n_comp () const |
Returns the part of N generated by the companion hole. | |
const Tenseur & | get_n_tot () const |
Returns the total N . | |
const Tenseur & | get_psi_auto () const |
Returns the part of ![]() | |
Cmp | return_psi_auto () const |
Returns the part of ![]() | |
const Tenseur & | get_psi_comp () const |
Returns the part of ![]() | |
const Tenseur & | get_psi_tot () const |
Returns the total ![]() | |
const Tenseur & | get_grad_psi_tot () const |
Returns the gradient of ![]() | |
const Tenseur & | get_grad_n_tot () const |
Returns the gradient of N . | |
const Tenseur & | get_shift_auto () const |
Returns the part of ![]() | |
Cmp | return_shift_auto (int i) const |
const Tenseur & | get_taij_auto () const |
Returns the part of ![]() | |
const Tenseur & | get_taij_comp () const |
Returns the part of ![]() | |
const Tenseur & | get_taij_tot () const |
Returns the total ![]() | |
const Tenseur & | get_tkij_tot () const |
Returns the total ![]() | |
const Tenseur & | get_tkij_auto () const |
Returns the part of ![]() | |
const Cmp | get_decouple () const |
Returns the function used to construct tkij_auto from tkij_tot . | |
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 . | |
void | fait_psi_comp (const Bhole &comp) |
Imports the part of ![]() 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 ![]() comp . | |
void | fait_taij_auto () |
Calculates the part of ![]() shift_auto . | |
void | regularise_shift (const Tenseur &comp) |
Corrects shift_auto in such a way that the total ![]() ![]() | |
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 ![]() | |
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 ~:
| |
void | solve_psi_seul (double relax) |
Solves the equation for ![]()
with the condition that | |
void | solve_shift_seul (double precis, double relax) |
Solves the equation for ![]()
| |
void | regularise_seul () |
Corrects the shift in the innermost shell, so that it remains ![]() ![]() | |
void | solve_lapse_with_ns (double relax, int bound_nn, double lim_nn) |
Solves the equation for N ~:
| |
void | solve_psi_with_ns (double relax) |
Solves the equation for ![]()
with the condition that | |
void | solve_shift_with_ns (const Et_bin_nsbh &ns, double precis, double relax, int bound_nn, double lim_nn) |
Solves the equation for ![]()
| |
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 ![]() | |
double | area () const |
Computes the area of the throat. | |
double | masse_adm_seul () const |
Calculates the ADM mass of the black hole using : ![]() | |
double | masse_komar_seul () const |
Calculates the Komar mass of the black hole using : ![]() | |
double | moment_seul_inf () const |
Calculates the angular momentum of the black hole using the formula at infinity : ![]() ![]() | |
double | moment_seul_hor () const |
Calculates the angular momentum of the black hole using the formula on the horizon : ![]() ![]() | |
double | local_momentum () const |
void | init_bhole_phi () |
Initiates the black hole for a resolution with ![]() | |
void | init_bhole_seul () |
Initiates for a single the black hole. | |
Protected Attributes | |
Map_af & | mp |
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 ![]() | |
Tenseur | psi_comp |
Part of ![]() | |
Tenseur | psi_tot |
Total ![]() | |
Tenseur | grad_n_tot |
Total gradient of N . | |
Tenseur | grad_psi_tot |
Total gradient of ![]() | |
Tenseur | shift_auto |
Part of ![]() | |
Tenseur_sym | taij_auto |
Part of ![]() | |
Tenseur_sym | taij_comp |
Part of ![]() | |
Tenseur_sym | taij_tot |
Total ![]() | |
Tenseur_sym | tkij_auto |
Auto ![]() | |
Tenseur_sym | tkij_tot |
Total ![]() | |
Cmp | decouple |
Function used to construct the part of ![]() ![]() | |
Friends | |
class | Bhole_binaire |
Binary black hole system. | |
class | Bin_ns_bh |
Binary NS-BH. |
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
.
Definition at line 263 of file bhole.h.
Bhole::Bhole | ( | Map_af & | mapping | ) |
Bhole::Bhole | ( | const Bhole & | source | ) |
Bhole::Bhole | ( | Map_af & | mpi, | |
FILE * | fich, | |||
bool | old = false | |||
) |
Constructor from a Map_af
and a file.
Definition at line 492 of file bhole.C.
References boost, decouple, fread_be(), Map::get_bvect_cart(), grad_n_tot, grad_psi_tot, mp, n_auto, omega, omega_local, psi_auto, regul, rot_state, Tenseur::set_etat_qcq(), Cmp::set_etat_zero(), Tenseur::set_etat_zero(), shift_auto, taij_auto, taij_comp, taij_tot, tkij_auto, and tkij_tot.
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 .
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 .
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 due to the companion ns
comp
.
The total is then calculated.
It also imports the gradient of and construct the total
.
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 due to the companion hole
comp
.
The total is then calculated.
It also imports the gradient of and construct the total
.
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 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 .
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] |
const Cmp Bhole::get_decouple | ( | ) | const [inline] |
const Tenseur& Bhole::get_grad_n_tot | ( | ) | const [inline] |
const Tenseur& Bhole::get_grad_psi_tot | ( | ) | const [inline] |
const Map_af& Bhole::get_mp | ( | ) | const [inline] |
const Tenseur& Bhole::get_n_auto | ( | ) | const [inline] |
const Tenseur& Bhole::get_n_comp | ( | ) | const [inline] |
const Tenseur& Bhole::get_n_tot | ( | ) | const [inline] |
double Bhole::get_omega | ( | ) | const [inline] |
double Bhole::get_omega_local | ( | ) | const [inline] |
const Tenseur& Bhole::get_psi_auto | ( | ) | const [inline] |
const Tenseur& Bhole::get_psi_comp | ( | ) | const [inline] |
const Tenseur& Bhole::get_psi_tot | ( | ) | const [inline] |
double Bhole::get_rayon | ( | ) | const [inline] |
double Bhole::get_regul | ( | ) | const [inline] |
double Bhole::get_rot_state | ( | ) | const [inline] |
const Tenseur& Bhole::get_shift_auto | ( | ) | const [inline] |
Returns the part of generated by the hole.
Definition at line 435 of file bhole.h.
References shift_auto.
const Tenseur& Bhole::get_taij_auto | ( | ) | const [inline] |
const Tenseur& Bhole::get_taij_comp | ( | ) | const [inline] |
const Tenseur& Bhole::get_taij_tot | ( | ) | const [inline] |
const Tenseur& Bhole::get_tkij_auto | ( | ) | const [inline] |
const Tenseur& Bhole::get_tkij_tot | ( | ) | const [inline] |
void 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 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 | ( | ) |
Initiates the black hole for a resolution with .
Definition at line 71 of file bhole_solve_phi.C.
References Cmp::annule(), decouple, grad_n_tot, grad_psi_tot, Tenseur::gradient(), log(), 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_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.
masse | [input] : ADM mass in LORENE's units. | |
moment | [input] : ![]() |
The radius of
*this
is supposed to be equal to .
The fields have then the following values~:
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 : .
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 : .
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 : 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 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 : where
.
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 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 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 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 ![]() |
Definition at line 399 of file bhole.C.
References omega, omega_local, regul, and shift_auto.
Cmp Bhole::return_n_auto | ( | ) | const [inline] |
Cmp Bhole::return_psi_auto | ( | ) | const [inline] |
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] |
void Bhole::set_omega | ( | double | ome | ) | [inline] |
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] |
void Bhole::set_rot_state | ( | int | rotation | ) | [inline] |
void 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 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 ~:
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 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 ~:
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 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 ~:
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 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 ~:
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 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 ~:
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 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 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.
friend class Bhole_binaire [friend] |
double* Bhole::boost [protected] |
Cmp Bhole::decouple [protected] |
Tenseur Bhole::grad_n_tot [protected] |
Tenseur Bhole::grad_psi_tot [protected] |
Tenseur Bhole::n_auto [protected] |
Tenseur Bhole::n_comp [protected] |
Tenseur Bhole::n_tot [protected] |
double Bhole::omega [protected] |
double Bhole::omega_local [protected] |
Tenseur Bhole::psi_auto [protected] |
Tenseur Bhole::psi_comp [protected] |
Tenseur Bhole::psi_tot [protected] |
double Bhole::rayon [protected] |
double Bhole::regul [protected] |
int Bhole::rot_state [protected] |
Tenseur Bhole::shift_auto [protected] |
Tenseur_sym Bhole::taij_auto [protected] |
Tenseur_sym Bhole::taij_comp [protected] |
Tenseur_sym Bhole::taij_tot [protected] |
Tenseur_sym Bhole::tkij_auto [protected] |
Tenseur_sym Bhole::tkij_tot [protected] |