Base class for black holes. More...
#include <blackhole.h>
Public Member Functions | |
Black_hole (Map &mp_i, bool Kerr_schild, double massbh) | |
Standard constructor. | |
Black_hole (const Black_hole &) | |
Copy constructor. | |
Black_hole (Map &mp_i, FILE *fich) | |
Constructor from a file (see sauve(FILE*) ). | |
virtual | ~Black_hole () |
Destructor. | |
void | operator= (const Black_hole &) |
Assignment to another Black_hole . | |
Map & | set_mp () |
Read/write of the mapping. | |
double & | set_mass_bh () |
Read/write of the gravitational mass of BH [{ m}]. | |
const Map & | get_mp () const |
Returns the mapping. | |
bool | is_kerrschild () const |
Returns true for a Kerr-Schild background, false for a Conformally flat one. | |
double | get_mass_bh () const |
Returns the gravitational mass of BH [{ m}]. | |
const Scalar & | get_lapconf () const |
Returns lapconf generated by the black hole. | |
const Scalar & | get_lapconf_rs () const |
Returns the part of lapconf from the numerical computation. | |
const Scalar & | get_lapse () const |
Returns the lapse function generated by the black hole. | |
const Vector & | get_shift () const |
Returns the shift vector generated by the black hole. | |
const Vector & | get_shift_rs () const |
Returns the part of the shift vector from the numerical computation. | |
const Scalar & | get_confo () const |
Returns the conformal factor generated by the black hole. | |
virtual void | sauve (FILE *) const |
Save in a file. | |
virtual double | mass_irr () const |
Irreducible mass of the black hole. | |
virtual double | mass_adm () const |
ADM mass. | |
virtual double | mass_kom () const |
Komar mass. | |
virtual double | rad_ah () const |
Radius of the apparent horizon. | |
double | spin_am_bh (bool bclapconf_nd, bool bclapconf_fs, const Tbl &xi_i, const double &phi_i, const double &theta_i, const int &nrk_phi, const int &nrk_theta) const |
Spin angular momentum. | |
const Tbl & | angu_mom_bh () const |
Total angular momentum. | |
const Valeur | bc_lapconf (bool neumann, bool first) const |
Boundary condition on the apparent horizon of the black hole for the lapse function: 2-D Valeur . | |
const Valeur | bc_shift_x (double omega_r) const |
Boundary condition on the apparent horizon of the black hole for the shift vector of the direction: 2-D Valeur . | |
const Valeur | bc_shift_y (double omega_r) const |
Boundary condition on the apparent horizon of the black hole for the shift vector of the direction: 2-D Valeur . | |
const Valeur | bc_shift_z () const |
Boundary condition on the apparent horizon of the black hole for the shift vector of the direction: 2-D Valeur . | |
const Valeur | bc_confo () const |
Boundary condition on the apparent horizon of the black hole for the conformal factor: 2-D Valeur . | |
void | extr_curv_bh () |
Computes taij , taij_quad from shift , lapse , confo . | |
void | equilibrium_spher (bool neumann, bool first, double spin_omega, double precis=1.e-14, double precis_shift=1.e-8) |
Computes a spherical, static black-hole by giving boundary conditions on the apparent horizon. | |
void | static_bh (bool neumann, bool first) |
Sets the metric quantities to a spherical, static blach-hole analytic solution. | |
double | rah_iso (bool neumann, bool first) const |
Computes a radius of apparent horizon (excised surface) in isotropic coordinates. | |
const Scalar | r_coord (bool neumann, bool first) const |
Expresses the areal radial coordinate by that in spatially isotropic coordinates. | |
Tbl | runge_kutta_phi_bh (const Tbl &xi_i, const double &phi_i, const int &nrk) const |
Compute a forth-order Runge-Kutta integration to the phi direction for the solution of the Killing vectors on the equator. | |
Tbl | runge_kutta_theta_bh (const Tbl &xi_i, const double &theta_i, const double &phi, const int &nrk) const |
Compute a forth-order Runge-Kutta integration to the theta direction for the solution of the Killing vectors. | |
Vector | killing_vect_bh (const Tbl &xi_i, const double &phi_i, const double &theta_i, const int &nrk_phi, const int &nrk_theta) const |
Compute the Killing vector of a black hole normalized so that its affine length is 2 M_PI. | |
Protected Member Functions | |
virtual void | del_deriv () const |
Deletes all the derived quantities. | |
void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. | |
virtual ostream & | operator>> (ostream &) const |
Operator >> (virtual function called by the operator <<). | |
Protected Attributes | |
Map & | mp |
Mapping associated with the black hole. | |
bool | kerrschild |
true for a Kerr-Schild background, false for a conformally flat background | |
double | mass_bh |
Gravitational mass of BH. | |
Scalar | lapconf |
A function (lapse function * conformal factor) lapconf generated by the black hole. | |
Scalar | lapconf_rs |
Part of lapconf from the numerical computation. | |
Scalar | lapconf_bh |
Part of lapconf from the analytic background. | |
Scalar | lapse |
Lapse function generated by the black hole. | |
Vector | shift |
Shift vector generated by the black hole. | |
Vector | shift_rs |
Part of the shift vector from the numerical computation. | |
Vector | shift_bh |
Part of the shift vector from the analytic background. | |
Scalar | confo |
Conformal factor generated by the black hole. | |
Sym_tensor | taij |
Trace of the extrinsic curvature. | |
Sym_tensor | taij_rs |
Part of the extrinsic curvature tensor. | |
Scalar | taij_quad |
Part of the scalar ![]() ![]() | |
Scalar | taij_quad_rs |
Part of the scalar. | |
Metric_flat | flat |
Flat metric defined on the mapping (Spherical components with respect to the mapping of the black hole ). | |
double * | p_mass_irr |
Conformal metric ![]() | |
double * | p_mass_adm |
Irreducible mass of the black hole. | |
double * | p_mass_kom |
ADM mass. | |
double * | p_rad_ah |
Komar mass. | |
double * | p_spin_am_bh |
Radius of the apparent horizon. | |
Tbl * | p_angu_mom_bh |
Spin angular momentum. | |
Friends | |
ostream & | operator<< (ostream &, const Black_hole &) |
Display. |
Base class for black holes.
()
According to the 3+1 formalism, the spacetime metric is written
where is the spatial metric.
Definition at line 70 of file blackhole.h.
Black_hole::Black_hole | ( | Map & | mp_i, | |
bool | Kerr_schild, | |||
double | massbh | |||
) |
Standard constructor.
mp_i | Mapping on which the black hole will be defined |
Definition at line 67 of file blackhole.C.
References confo, lapconf, lapconf_bh, lapconf_rs, lapse, set_der_0x0(), Tensor::set_etat_zero(), shift, shift_bh, shift_rs, Scalar::std_spectral_base(), taij, taij_quad, taij_quad_rs, and taij_rs.
Black_hole::Black_hole | ( | const Black_hole & | bh | ) |
Black_hole::Black_hole | ( | Map & | mp_i, | |
FILE * | fich | |||
) |
Constructor from a file (see sauve(FILE*)
).
mp_i | Mapping on which the black hole will be defined | |
fich | input file (must have been created by the function sauve ) |
Definition at line 138 of file blackhole.C.
References confo, kerrschild, lapconf, lapconf_bh, lapconf_rs, lapse, mass_bh, set_der_0x0(), Tensor::set_etat_zero(), shift, shift_bh, shift_rs, Vector::std_spectral_base(), Scalar::std_spectral_base(), taij, taij_quad, taij_quad_rs, and taij_rs.
Black_hole::~Black_hole | ( | ) | [virtual] |
const Tbl & Black_hole::angu_mom_bh | ( | ) | const |
Total angular momentum.
Definition at line 630 of file blackhole_global.C.
References Tensor::annule_domain(), Tbl::annule_hard(), Map::cosp, Map::cost, Scalar::dec_dzpuis(), Map::get_bvect_cart(), Map_af::integrale_surface(), mp, p_angu_mom_bh, Scalar::raccord(), Tbl::set(), Tensor::set(), Map::sinp, Map::sint, Scalar::std_spectral_base(), taij, Map::val_r(), Map::x, Map::y, and Map::z.
const Valeur Black_hole::bc_confo | ( | ) | const |
Boundary condition on the apparent horizon of the black hole for the conformal factor: 2-D Valeur
.
Definition at line 406 of file blackhole_bc.C.
References confo, Map::cosp, Map::cost, Scalar::dsdr(), Mg3d::get_angu(), Scalar::get_dzpuis(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), kerrschild, lapconf, mass_bh, mp, pow(), Map::r, Valeur::set(), Scalar::set_dzpuis(), shift, Map::sinp, Map::sint, Valeur::std_base_scal(), Scalar::std_spectral_base(), and Scalar::val_grid_point().
const Valeur Black_hole::bc_lapconf | ( | bool | neumann, | |
bool | first | |||
) | const |
Boundary condition on the apparent horizon of the black hole for the lapse function: 2-D Valeur
.
Definition at line 67 of file blackhole_bc.C.
References Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), kerrschild, lapconf, mp, Map::r, Valeur::set(), Valeur::std_base_scal(), Scalar::std_spectral_base(), and Scalar::val_grid_point().
const Valeur Black_hole::bc_shift_x | ( | double | omega_r | ) | const |
Boundary condition on the apparent horizon of the black hole for the shift vector of the direction: 2-D Valeur
.
Definition at line 207 of file blackhole_bc.C.
References Valeur::base, confo, Map::cosp, Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), kerrschild, lapconf, mp, pow(), Map::r, Valeur::set(), Map::sint, Mg3d::std_base_vect_cart(), Scalar::std_spectral_base(), Scalar::val_grid_point(), and Map::y.
const Valeur Black_hole::bc_shift_y | ( | double | omega_r | ) | const |
Boundary condition on the apparent horizon of the black hole for the shift vector of the direction: 2-D Valeur
.
Definition at line 277 of file blackhole_bc.C.
References Valeur::base, confo, Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), kerrschild, lapconf, mp, pow(), Map::r, Valeur::set(), Map::sinp, Map::sint, Mg3d::std_base_vect_cart(), Scalar::std_spectral_base(), Scalar::val_grid_point(), and Map::x.
const Valeur Black_hole::bc_shift_z | ( | ) | const |
Boundary condition on the apparent horizon of the black hole for the shift vector of the direction: 2-D Valeur
.
Reimplemented in Hole_bhns.
Definition at line 346 of file blackhole_bc.C.
References Valeur::base, confo, Map::cost, Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), kerrschild, lapconf, mp, pow(), Map::r, Valeur::set(), Mg3d::std_base_vect_cart(), Scalar::std_spectral_base(), and Scalar::val_grid_point().
void Black_hole::del_deriv | ( | ) | const [protected, virtual] |
Deletes all the derived quantities.
Reimplemented in Hole_bhns.
Definition at line 201 of file blackhole.C.
References p_angu_mom_bh, p_mass_adm, p_mass_irr, p_mass_kom, p_rad_ah, p_spin_am_bh, and set_der_0x0().
void Black_hole::equilibrium_spher | ( | bool | neumann, | |
bool | first, | |||
double | spin_omega, | |||
double | precis = 1.e-14 , |
|||
double | precis_shift = 1.e-8 | |||
) |
Computes a spherical, static black-hole by giving boundary conditions on the apparent horizon.
neumann | [input] true for the neumann bc, false for the dirichlet one for the lapse | |
first | [input] true for the first type of bc, false for the second type | |
spin_omega | [input] spin parameter of the black hole | |
precis | [input] threshold in the relative difference of a metric quantity between two succesive steps to stop the iterative procedure (default value: 1.e-14) | |
precis_shift | [input] threshold in the relative difference of the shift vector between two succesive steps to stop the iterative procedure (default value: 1.e-8) |
Definition at line 66 of file blackhole_eq_spher.C.
References Tensor::annule_domain(), bc_confo(), bc_lapconf(), bc_shift_x(), bc_shift_y(), bc_shift_z(), confo, contract(), Map::cosp, Map::cost, del_deriv(), Scalar::deriv(), diffrel(), extr_curv_bh(), Mg3d::get_angu(), Map::get_bvect_cart(), Scalar::get_dzpuis(), Map::get_mg(), Mg3d::get_nzone(), Scalar::inc_dzpuis(), kerrschild, lapconf, lapconf_bh, lapconf_rs, lapse, mass_adm(), mass_bh, mass_irr(), mass_kom(), mp, Scalar::poisson_dirichlet(), Scalar::poisson_neumann(), pow(), Map::r, r_coord(), Scalar::raccord(), rad_ah(), Tbl::set(), Vector::set(), Tensor::set(), Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Tensor::set_etat_zero(), shift, shift_bh, shift_rs, Map::sinp, Map::sint, sqrt(), Vector::std_spectral_base(), Scalar::std_spectral_base(), taij, and taij_quad.
void Black_hole::extr_curv_bh | ( | ) |
Computes taij
, taij_quad
from shift
, lapse
, confo
.
Definition at line 58 of file blackhole_extr_curv.C.
References Tensor::annule_domain(), Metric_flat::con(), confo, Map::cosp, Map::cost, Metric_flat::cov(), flat, Map::get_bvect_cart(), Scalar::inc_dzpuis(), kerrschild, lapconf, mass_bh, mp, pow(), Map::r, Scalar::raccord(), Tensor::set(), Tensor::set_etat_qcq(), shift, shift_rs, Map::sinp, Map::sint, sqrt(), Tensor::std_spectral_base(), Scalar::std_spectral_base(), taij, taij_quad, taij_quad_rs, and taij_rs.
const Scalar& Black_hole::get_confo | ( | ) | const [inline] |
Returns the conformal factor generated by the black hole.
Definition at line 237 of file blackhole.h.
References confo.
const Scalar& Black_hole::get_lapconf | ( | ) | const [inline] |
Returns lapconf generated by the black hole.
Definition at line 220 of file blackhole.h.
References lapconf.
const Scalar& Black_hole::get_lapconf_rs | ( | ) | const [inline] |
Returns the part of lapconf from the numerical computation.
Definition at line 223 of file blackhole.h.
References lapconf_rs.
const Scalar& Black_hole::get_lapse | ( | ) | const [inline] |
Returns the lapse function generated by the black hole.
Definition at line 226 of file blackhole.h.
References lapse.
double Black_hole::get_mass_bh | ( | ) | const [inline] |
Returns the gravitational mass of BH [{ m}].
Definition at line 217 of file blackhole.h.
References mass_bh.
const Map& Black_hole::get_mp | ( | ) | const [inline] |
const Vector& Black_hole::get_shift | ( | ) | const [inline] |
Returns the shift vector generated by the black hole.
Definition at line 229 of file blackhole.h.
References shift.
const Vector& Black_hole::get_shift_rs | ( | ) | const [inline] |
Returns the part of the shift vector from the numerical computation.
Definition at line 234 of file blackhole.h.
References shift_rs.
bool Black_hole::is_kerrschild | ( | ) | const [inline] |
Returns true
for a Kerr-Schild background, false
for a Conformally flat one.
Definition at line 214 of file blackhole.h.
References kerrschild.
Vector Black_hole::killing_vect_bh | ( | const Tbl & | xi_i, | |
const double & | phi_i, | |||
const double & | theta_i, | |||
const int & | nrk_phi, | |||
const int & | nrk_theta | |||
) | const |
Compute the Killing vector of a black hole normalized so that its affine length is 2 M_PI.
xi_i | [input] initial set of the Killing vectors at the starting point on the equator | |
phi_i | [input] initial phi coordinate at which the integration starts | |
theta_i | [input] initial theta coordinate at which the integration starts | |
nrk_phi | [input] number of the Runge-Kutta integration between two successive collocation points for the phi direction | |
nrk_theta | [input] number of the Runge-Kutta integration between two successive collocation points for the theta direction |
Definition at line 62 of file blackhole_killing.C.
References confo, Map::get_bvect_spher(), Tbl::get_dim(), Map::get_mg(), Tbl::get_ndim(), Mg3d::get_nr(), kerrschild, mp, pow(), Map::r, rad_ah(), runge_kutta_phi_bh(), runge_kutta_theta_bh(), Tbl::set(), Tbl::set_etat_qcq(), Scalar::set_grid_point(), Map::sint, Scalar::std_spectral_base(), Scalar::val_grid_point(), and Scalar::val_point().
double Black_hole::mass_adm | ( | ) | const [virtual] |
ADM mass.
Definition at line 99 of file blackhole_global.C.
References Tensor::annule_domain(), confo, Map::cosp, Map::cost, Scalar::dec_dzpuis(), Scalar::dsdr(), Map::get_bvect_cart(), Scalar::integrale(), Map_af::integrale_surface(), Map_af::integrale_surface_infini(), kerrschild, lapconf, mp, p_mass_adm, pow(), Map::r, Scalar::raccord(), Tensor::set(), Scalar::set_etat_qcq(), shift, Map::sinp, Map::sint, Scalar::std_spectral_base(), taij_quad, and Map::val_r().
double Black_hole::mass_irr | ( | ) | const [virtual] |
Irreducible mass of the black hole.
Definition at line 65 of file blackhole_global.C.
References Tensor::annule_domain(), confo, Map_af::integrale_surface(), mp, p_mass_irr, pow(), Scalar::raccord(), sqrt(), Scalar::std_spectral_base(), and Map::val_r().
double Black_hole::mass_kom | ( | ) | const [virtual] |
Komar mass.
Definition at line 285 of file blackhole_global.C.
References Tensor::annule_domain(), confo, Map::cosp, Map::cost, Scalar::dec_dzpuis(), Scalar::deriv(), Scalar::dsdr(), Map::get_bvect_cart(), Scalar::integrale(), Map_af::integrale_surface(), Map_af::integrale_surface_infini(), kerrschild, lapconf, mp, p_mass_kom, pow(), Map::r, Scalar::raccord(), Tensor::set(), Scalar::set_etat_qcq(), Map::sinp, Map::sint, Scalar::std_spectral_base(), taij_quad, and Map::val_r().
void Black_hole::operator= | ( | const Black_hole & | bh | ) |
Assignment to another Black_hole
.
Reimplemented in Hole_bhns.
Definition at line 232 of file blackhole.C.
References confo, del_deriv(), flat, kerrschild, lapconf, lapconf_bh, lapconf_rs, lapse, mass_bh, mp, shift, shift_bh, shift_rs, taij, taij_quad, taij_quad_rs, and taij_rs.
ostream & Black_hole::operator>> | ( | ostream & | ost | ) | const [protected, virtual] |
Operator >> (virtual function called by the operator <<).
Reimplemented in Hole_bhns.
Definition at line 283 of file blackhole.C.
References confo, Map::get_mg(), Mg3d::get_nt(), kerrschild, lapconf, lapse, mass_adm(), mass_bh, mass_irr(), mass_kom(), mp, shift, and Scalar::val_grid_point().
const Scalar Black_hole::r_coord | ( | bool | neumann, | |
bool | first | |||
) | const |
Expresses the areal radial coordinate by that in spatially isotropic coordinates.
neumann | [input] true for the neumann bc, false for the dirichlet one for the lapse | |
first | [input] true for the first type of bc, false for the second type |
Definition at line 62 of file blackhole_r_coord.C.
References Tensor::annule_domain(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), log(), mass_bh, mp, Map::r, Scalar::raccord(), Scalar::set_grid_point(), sqrt(), Scalar::std_spectral_base(), and Scalar::val_grid_point().
double Black_hole::rad_ah | ( | ) | const [virtual] |
Radius of the apparent horizon.
Definition at line 501 of file blackhole_global.C.
References mp, p_rad_ah, Map::r, Scalar::std_spectral_base(), and Scalar::val_grid_point().
double Black_hole::rah_iso | ( | bool | neumann, | |
bool | first | |||
) | const |
Computes a radius of apparent horizon (excised surface) in isotropic coordinates.
neumann | [input] true for the neumann bc, false for the dirichlet one for the lapse | |
first | [input] true for the first type of bc, false for the second type |
Definition at line 65 of file blackhole_rah_iso.C.
Tbl Black_hole::runge_kutta_phi_bh | ( | const Tbl & | xi_i, | |
const double & | phi_i, | |||
const int & | nrk | |||
) | const |
Compute a forth-order Runge-Kutta integration to the phi direction for the solution of the Killing vectors on the equator.
xi_i | [input] initial set of the Killing vectors at the starting point on the equator | |
phi_i | [input] initial phi coordinate at which the integration starts | |
nrk | [input] number of the Runge-Kutta integration between two successive collocation points |
Definition at line 63 of file blackhole_rk_phi.C.
References confo, Scalar::dsdt(), Map::get_mg(), Mg3d::get_np(), kerrschild, Scalar::lapang(), mp, rad_ah(), Tbl::set(), Tbl::set_etat_qcq(), Scalar::std_spectral_base(), and Scalar::val_point().
Tbl Black_hole::runge_kutta_theta_bh | ( | const Tbl & | xi_i, | |
const double & | theta_i, | |||
const double & | phi, | |||
const int & | nrk | |||
) | const |
Compute a forth-order Runge-Kutta integration to the theta direction for the solution of the Killing vectors.
xi_i | [input] initial set of the Killing vectors at the starting point on the equator | |
theta_i | [input] initial theta coordinate at which the integration starts | |
phi | [input] fixed phi coordinate during integration to the theta direction | |
nrk | [input] number of the Runge-Kutta integration between two successive collocation points |
Definition at line 63 of file blackhole_rk_theta.C.
References confo, Map::get_mg(), Mg3d::get_nt(), kerrschild, Scalar::lapang(), mp, rad_ah(), Tbl::set(), Tbl::set_etat_qcq(), Scalar::std_spectral_base(), Scalar::stdsdp(), and Scalar::val_point().
void Black_hole::sauve | ( | FILE * | fich | ) | const [virtual] |
Save in a file.
Reimplemented in Hole_bhns.
Definition at line 263 of file blackhole.C.
References confo, kerrschild, lapconf_rs, mass_bh, Tensor::sauve(), Scalar::sauve(), and shift_rs.
void Black_hole::set_der_0x0 | ( | ) | const [protected] |
Sets to 0x0
all the pointers on derived quantities.
Reimplemented in Hole_bhns.
Definition at line 214 of file blackhole.C.
References p_angu_mom_bh, p_mass_adm, p_mass_irr, p_mass_kom, p_rad_ah, and p_spin_am_bh.
double& Black_hole::set_mass_bh | ( | ) | [inline] |
Read/write of the gravitational mass of BH [{ m}].
Definition at line 203 of file blackhole.h.
References mass_bh.
Map& Black_hole::set_mp | ( | ) | [inline] |
double Black_hole::spin_am_bh | ( | bool | bclapconf_nd, | |
bool | bclapconf_fs, | |||
const Tbl & | xi_i, | |||
const double & | phi_i, | |||
const double & | theta_i, | |||
const int & | nrk_phi, | |||
const int & | nrk_theta | |||
) | const |
Spin angular momentum.
Definition at line 524 of file blackhole_global.C.
References Tensor::annule_domain(), confo, Map::cosp, Map::cost, Scalar::dec_dzpuis(), Map::get_bvect_cart(), Map::get_bvect_spher(), Map_af::integrale_surface(), kerrschild, killing_vect_bh(), mp, p_spin_am_bh, pow(), Scalar::raccord(), Tensor::set(), Map::sinp, Map::sint, Scalar::std_spectral_base(), taij, and Map::val_r().
void Black_hole::static_bh | ( | bool | neumann, | |
bool | first | |||
) |
Sets the metric quantities to a spherical, static blach-hole analytic solution.
neumann | [input] true for the neumann bc, false for the dirichlet one for the lapse | |
first | [input] true for the first type of bc, false for the second type |
Definition at line 56 of file blackhole_static.C.
References Tensor::annule_domain(), confo, Map::cosp, Map::cost, Map::get_bvect_cart(), kerrschild, lapconf, lapse, mass_bh, mp, pow(), Map::r, r_coord(), Scalar::raccord(), Vector::set(), Tensor::set(), shift, Map::sinp, Map::sint, sqrt(), Vector::std_spectral_base(), and Scalar::std_spectral_base().
ostream& operator<< | ( | ostream & | , | |
const Black_hole & | ||||
) | [friend] |
Display.
Scalar Black_hole::confo [protected] |
Conformal factor generated by the black hole.
Definition at line 114 of file blackhole.h.
Metric_flat Black_hole::flat [protected] |
Flat metric defined on the mapping (Spherical components with respect to the mapping of the black hole ).
Definition at line 139 of file blackhole.h.
bool Black_hole::kerrschild [protected] |
true
for a Kerr-Schild background, false
for a conformally flat background
Definition at line 81 of file blackhole.h.
Scalar Black_hole::lapconf [protected] |
A function (lapse function * conformal factor) lapconf generated by the black hole.
Definition at line 93 of file blackhole.h.
Scalar Black_hole::lapconf_bh [protected] |
Part of lapconf from the analytic background.
Definition at line 99 of file blackhole.h.
Scalar Black_hole::lapconf_rs [protected] |
Part of lapconf from the numerical computation.
Definition at line 96 of file blackhole.h.
Scalar Black_hole::lapse [protected] |
Lapse function generated by the black hole.
Definition at line 102 of file blackhole.h.
double Black_hole::mass_bh [protected] |
Gravitational mass of BH.
Definition at line 84 of file blackhole.h.
Map& Black_hole::mp [protected] |
Mapping associated with the black hole.
Definition at line 76 of file blackhole.h.
Tbl* Black_hole::p_angu_mom_bh [mutable, protected] |
Spin angular momentum.
Total angular momentum of the black hole
Definition at line 158 of file blackhole.h.
double* Black_hole::p_mass_adm [mutable, protected] |
Irreducible mass of the black hole.
Definition at line 149 of file blackhole.h.
double* Black_hole::p_mass_irr [mutable, protected] |
Conformal metric .
Definition at line 147 of file blackhole.h.
double* Black_hole::p_mass_kom [mutable, protected] |
ADM mass.
Definition at line 151 of file blackhole.h.
double* Black_hole::p_rad_ah [mutable, protected] |
Komar mass.
Definition at line 153 of file blackhole.h.
double* Black_hole::p_spin_am_bh [mutable, protected] |
Radius of the apparent horizon.
Definition at line 155 of file blackhole.h.
Vector Black_hole::shift [protected] |
Shift vector generated by the black hole.
Definition at line 105 of file blackhole.h.
Vector Black_hole::shift_bh [protected] |
Part of the shift vector from the analytic background.
Definition at line 111 of file blackhole.h.
Vector Black_hole::shift_rs [protected] |
Part of the shift vector from the numerical computation.
Definition at line 108 of file blackhole.h.
Sym_tensor Black_hole::taij [protected] |
Trace of the extrinsic curvature.
Extrinsic curvature tensor generated by
shift
, lapse
, and confo
.
Definition at line 123 of file blackhole.h.
Scalar Black_hole::taij_quad [protected] |
Part of the scalar generated by
.
Definition at line 131 of file blackhole.h.
Scalar Black_hole::taij_quad_rs [protected] |
Part of the scalar.
Definition at line 134 of file blackhole.h.
Sym_tensor Black_hole::taij_rs [protected] |
Part of the extrinsic curvature tensor.
Definition at line 126 of file blackhole.h.