246 virtual void sauve(FILE *)
const ;
253 virtual ostream&
operator>>(ostream& )
const ;
268 virtual double rad_ah()
const ;
271 double spin_am_bh(
bool bclapconf_nd,
bool bclapconf_fs,
272 const Tbl& xi_i,
const double& phi_i,
273 const double& theta_i,
const int& nrk_phi,
274 const int& nrk_theta)
const ;
336 double precis = 1.e-14,
337 double precis_shift = 1.e-8) ;
347 void static_bh(
bool neumann,
bool first) ;
357 double rah_iso(
bool neumann,
bool first)
const ;
381 const int& nrk)
const ;
397 const double& phi,
const int& nrk)
const ;
416 const double& theta_i,
const int& nrk_phi,
417 const int& nrk_theta)
const ;
420ostream& operator<<(ostream& ,
const Black_hole& ) ;
Base class for black holes.
Vector shift_bh
Part of the shift vector from the analytic background.
double rah_iso(bool neumann, bool first) const
Computes a radius of apparent horizon (excised surface) in isotropic coordinates.
virtual double mass_kom() const
Komar mass.
Vector shift_rs
Part of the shift vector from the numerical computation.
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 ve...
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 Valeur bc_confo() const
Boundary condition on the apparent horizon of the black hole for the conformal factor: 2-D Valeur.
double & set_mass_bh()
Read/write of the gravitational mass of BH [{\tt m_unit}].
Scalar taij_quad
Part of the scalar generated by .
const Scalar & get_lapse() const
Returns the lapse function generated by the black hole.
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 \fx\f direct...
Scalar taij_quad_rs
Part of the scalar.
Sym_tensor taij
Trace of the extrinsic curvature.
const Scalar r_coord(bool neumann, bool first) const
Expresses the areal radial coordinate by that in spatially isotropic coordinates.
const Vector & get_shift() const
Returns the shift vector generated by the black hole.
Map & set_mp()
Read/write of the mapping.
virtual double mass_irr() const
Irreducible mass of the black hole.
double * p_spin_am_bh
Radius of the apparent horizon.
Sym_tensor taij_rs
Part of the extrinsic curvature tensor.
Map & mp
Mapping associated with the black hole.
bool is_kerrschild() const
Returns true for a Kerr-Schild background, false for a Conformally flat one.
Scalar lapconf
A function (lapse function * conformal factor) lapconf generated by the black hole.
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.
Vector shift
Shift vector generated by the black hole.
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the black hol...
virtual double rad_ah() const
Radius of the apparent horizon.
const Scalar & get_lapconf_rs() const
Returns the part of lapconf from the numerical computation.
Scalar lapconf_rs
Part of lapconf from the numerical computation.
const Scalar & get_confo() const
Returns the conformal factor generated by the black hole.
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
void static_bh(bool neumann, bool first)
Sets the metric quantities to a spherical, static blach-hole analytic solution.
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 ...
friend ostream & operator<<(ostream &, const Black_hole &)
Display.
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.
const Map & get_mp() const
Returns the mapping.
Tbl * p_angu_mom_bh
Spin angular momentum.
virtual ~Black_hole()
Destructor.
const Scalar & get_lapconf() const
Returns lapconf generated by the black hole.
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.
virtual double mass_adm() const
ADM mass.
double * p_rad_ah
Komar mass.
virtual void sauve(FILE *) const
Save in a file.
const Valeur bc_shift_z() const
Boundary condition on the apparent horizon of the black hole for the shift vector of the \fz\f direct...
Scalar lapconf_bh
Part of lapconf from the analytic background.
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 \fy\f direct...
double * p_mass_kom
ADM mass.
void operator=(const Black_hole &)
Assignment to another Black_hole.
const Vector & get_shift_rs() const
Returns the part of the shift vector from the numerical computation.
const Tbl & angu_mom_bh() const
Total angular momentum.
Scalar lapse
Lapse function generated by the black hole.
bool kerrschild
true for a Kerr-Schild background, false for a conformally flat background
Scalar confo
Conformal factor generated by the black hole.
double mass_bh
Gravitational mass of BH.
double * p_mass_adm
Irreducible mass of the black hole.
void extr_curv_bh()
Computes taij , taij_quad from shift , lapse , confo .
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<)
double get_mass_bh() const
Returns the gravitational mass of BH [{\tt m_unit}].
double * p_mass_irr
Conformal metric .
virtual void del_deriv() const
Deletes all the derived quantities.
Base class for coordinate mappings.
Flat metric for tensor calculation.
Tensor field of valence 0 (or component of a tensorial field).
Class intended to describe valence-2 symmetric tensors.
Values and coefficients of a (real-value) function.
Tensor field of valence 1.