LORENE
Lorene::Black_hole Class Reference

Base class for black holes. More...

#include <blackhole.h>

Inheritance diagram for Lorene::Black_hole:
Lorene::Hole_bhns

Public Member Functions

 Black_hole (Map &mp_i, bool Kerr_schild, double massbh)
 Standard constructor. More...
 
 Black_hole (const Black_hole &)
 Copy constructor. More...
 
 Black_hole (Map &mp_i, FILE *fich)
 Constructor from a file (see sauve(FILE*) ) More...
 
virtual ~Black_hole ()
 Destructor. More...
 
void operator= (const Black_hole &)
 Assignment to another Black_hole. More...
 
Mapset_mp ()
 Read/write of the mapping. More...
 
double & set_mass_bh ()
 Read/write of the gravitational mass of BH [{ m_unit}]. More...
 
const Mapget_mp () const
 Returns the mapping. More...
 
bool is_kerrschild () const
 Returns true for a Kerr-Schild background, false for a Conformally flat one. More...
 
double get_mass_bh () const
 Returns the gravitational mass of BH [{ m_unit}]. More...
 
const Scalarget_lapconf () const
 Returns lapconf generated by the black hole. More...
 
const Scalarget_lapconf_rs () const
 Returns the part of lapconf from the numerical computation. More...
 
const Scalarget_lapse () const
 Returns the lapse function generated by the black hole. More...
 
const Vectorget_shift () const
 Returns the shift vector generated by the black hole. More...
 
const Vectorget_shift_rs () const
 Returns the part of the shift vector from the numerical computation. More...
 
const Scalarget_confo () const
 Returns the conformal factor generated by the black hole. More...
 
virtual void sauve (FILE *) const
 Save in a file. More...
 
virtual double mass_irr () const
 Irreducible mass of the black hole. More...
 
virtual double mass_adm () const
 ADM mass. More...
 
virtual double mass_kom () const
 Komar mass. More...
 
virtual double rad_ah () const
 Radius of the apparent horizon. More...
 
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. More...
 
const Tblangu_mom_bh () const
 Total angular momentum. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
const Valeur bc_confo () const
 Boundary condition on the apparent horizon of the black hole for the conformal factor: 2-D Valeur. More...
 
void extr_curv_bh ()
 Computes taij , taij_quad from shift , lapse , confo . More...
 
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. More...
 
void static_bh (bool neumann, bool first)
 Sets the metric quantities to a spherical, static blach-hole analytic solution. More...
 
double rah_iso (bool neumann, bool first) const
 Computes a radius of apparent horizon (excised surface) in isotropic coordinates. More...
 
const Scalar r_coord (bool neumann, bool first) const
 Expresses the areal radial coordinate by that in spatially isotropic coordinates. More...
 
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. More...
 
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. More...
 
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. More...
 

Protected Member Functions

virtual void del_deriv () const
 Deletes all the derived quantities. More...
 
void set_der_0x0 () const
 Sets to 0x0 all the pointers on derived quantities. More...
 
virtual ostream & operator>> (ostream &) const
 Operator >> (virtual function called by the operator <<) More...
 

Protected Attributes

Mapmp
 Mapping associated with the black hole. More...
 
bool kerrschild
 true for a Kerr-Schild background, false for a conformally flat background More...
 
double mass_bh
 Gravitational mass of BH. More...
 
Scalar lapconf
 A function (lapse function * conformal factor) lapconf generated by the black hole. More...
 
Scalar lapconf_rs
 Part of lapconf from the numerical computation. More...
 
Scalar lapconf_bh
 Part of lapconf from the analytic background. More...
 
Scalar lapse
 Lapse function generated by the black hole. More...
 
Vector shift
 Shift vector generated by the black hole. More...
 
Vector shift_rs
 Part of the shift vector from the numerical computation. More...
 
Vector shift_bh
 Part of the shift vector from the analytic background. More...
 
Scalar confo
 Conformal factor generated by the black hole. More...
 
Sym_tensor taij
 Trace of the extrinsic curvature. More...
 
Sym_tensor taij_rs
 Part of the extrinsic curvature tensor. More...
 
Scalar taij_quad
 Part of the scalar $\eta_{ik} \eta_{jl} A^{ij} A^{kl}$ generated by $A_{ij}$. More...
 
Scalar taij_quad_rs
 Part of the scalar. More...
 
Metric_flat flat
 Flat metric defined on the mapping (Spherical components with respect to the mapping of the black hole ). More...
 
double * p_mass_irr
 Conformal metric $\tilde \gamma_{ij}$. More...
 
double * p_mass_adm
 Irreducible mass of the black hole. More...
 
double * p_mass_kom
 ADM mass. More...
 
double * p_rad_ah
 Komar mass. More...
 
double * p_spin_am_bh
 Radius of the apparent horizon. More...
 
Tblp_angu_mom_bh
 Spin angular momentum. More...
 

Friends

ostream & operator<< (ostream &, const Black_hole &)
 Display. More...
 

Detailed Description

Base class for black holes.

()

According to the 3+1 formalism, the spacetime metric is written

\[ ds^2 = - \alpha^2 dt^2 + \gamma_{ij} ( dx^i + \beta^i dt ) ( dx^j + \beta^j dt ) \]

where $\gamma_{ij}$ is the spatial metric.

Definition at line 74 of file blackhole.h.

Constructor & Destructor Documentation

◆ Black_hole() [1/3]

Lorene::Black_hole::Black_hole ( Map mp_i,
bool  Kerr_schild,
double  massbh 
)

Standard constructor.

Parameters
mp_iMapping on which the black hole will be defined

Definition at line 74 of file blackhole.C.

References confo, lapconf, lapconf_bh, lapconf_rs, lapse, set_der_0x0(), Lorene::Tensor::set_etat_zero(), shift, shift_bh, shift_rs, Lorene::Scalar::std_spectral_base(), taij, taij_quad, taij_quad_rs, and taij_rs.

◆ Black_hole() [2/3]

Lorene::Black_hole::Black_hole ( const Black_hole bh)

Copy constructor.

Definition at line 121 of file blackhole.C.

References set_der_0x0().

◆ Black_hole() [3/3]

Lorene::Black_hole::Black_hole ( Map mp_i,
FILE *  fich 
)

Constructor from a file (see sauve(FILE*) )

Parameters
mp_iMapping on which the black hole will be defined
fichinput file (must have been created by the function sauve )

Definition at line 145 of file blackhole.C.

References confo, kerrschild, lapconf, lapconf_bh, lapconf_rs, lapse, mass_bh, set_der_0x0(), Lorene::Tensor::set_etat_zero(), shift, shift_bh, shift_rs, Lorene::Vector::std_spectral_base(), Lorene::Scalar::std_spectral_base(), taij, taij_quad, taij_quad_rs, and taij_rs.

◆ ~Black_hole()

Lorene::Black_hole::~Black_hole ( )
virtual

Destructor.

Definition at line 197 of file blackhole.C.

References del_deriv().

Member Function Documentation

◆ angu_mom_bh()

◆ bc_confo()

◆ bc_lapconf()

const Valeur Lorene::Black_hole::bc_lapconf ( bool  neumann,
bool  first 
) const

◆ bc_shift_x()

const Valeur Lorene::Black_hole::bc_shift_x ( double  omega_r) const

◆ bc_shift_y()

const Valeur Lorene::Black_hole::bc_shift_y ( double  omega_r) const

◆ bc_shift_z()

const Valeur Lorene::Black_hole::bc_shift_z ( ) const

◆ del_deriv()

void Lorene::Black_hole::del_deriv ( ) const
protectedvirtual

Deletes all the derived quantities.

Reimplemented in Lorene::Hole_bhns.

Definition at line 208 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().

◆ equilibrium_spher()

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

Parameters
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 73 of file blackhole_eq_spher.C.

References Lorene::Tensor::annule_domain(), bc_confo(), bc_lapconf(), confo, Lorene::contract(), Lorene::Map::cosp, Lorene::Map::cost, Lorene::Scalar::deriv(), extr_curv_bh(), Lorene::Mg3d::get_angu(), Lorene::Map::get_bvect_cart(), Lorene::Scalar::get_dzpuis(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::inc_dzpuis(), kerrschild, lapconf, lapconf_bh, lapconf_rs, lapse, mass_bh, mp, Lorene::Scalar::poisson_dirichlet(), Lorene::Scalar::poisson_neumann(), Lorene::pow(), Lorene::Map::r, r_coord(), Lorene::Scalar::raccord(), Lorene::Vector::set(), Lorene::Tensor::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Tensor::set_etat_qcq(), Lorene::Tensor::set_etat_zero(), shift, shift_bh, shift_rs, Lorene::Map::sinp, Lorene::Map::sint, Lorene::sqrt(), Lorene::Vector::std_spectral_base(), Lorene::Scalar::std_spectral_base(), taij, and taij_quad.

◆ extr_curv_bh()

◆ get_confo()

const Scalar& Lorene::Black_hole::get_confo ( ) const
inline

Returns the conformal factor generated by the black hole.

Definition at line 241 of file blackhole.h.

References confo.

◆ get_lapconf()

const Scalar& Lorene::Black_hole::get_lapconf ( ) const
inline

Returns lapconf generated by the black hole.

Definition at line 224 of file blackhole.h.

References lapconf.

◆ get_lapconf_rs()

const Scalar& Lorene::Black_hole::get_lapconf_rs ( ) const
inline

Returns the part of lapconf from the numerical computation.

Definition at line 227 of file blackhole.h.

References lapconf_rs.

◆ get_lapse()

const Scalar& Lorene::Black_hole::get_lapse ( ) const
inline

Returns the lapse function generated by the black hole.

Definition at line 230 of file blackhole.h.

References lapse.

◆ get_mass_bh()

double Lorene::Black_hole::get_mass_bh ( ) const
inline

Returns the gravitational mass of BH [{ m_unit}].

Definition at line 221 of file blackhole.h.

References mass_bh.

◆ get_mp()

const Map& Lorene::Black_hole::get_mp ( ) const
inline

Returns the mapping.

Definition at line 213 of file blackhole.h.

References mp.

◆ get_shift()

const Vector& Lorene::Black_hole::get_shift ( ) const
inline

Returns the shift vector generated by the black hole.

Definition at line 233 of file blackhole.h.

References shift.

◆ get_shift_rs()

const Vector& Lorene::Black_hole::get_shift_rs ( ) const
inline

Returns the part of the shift vector from the numerical computation.

Definition at line 238 of file blackhole.h.

References shift_rs.

◆ is_kerrschild()

bool Lorene::Black_hole::is_kerrschild ( ) const
inline

Returns true for a Kerr-Schild background, false for a Conformally flat one.

Definition at line 218 of file blackhole.h.

References kerrschild.

◆ killing_vect_bh()

Vector Lorene::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.

Parameters
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 69 of file blackhole_killing.C.

References confo, Lorene::Map::get_bvect_spher(), Lorene::Tbl::get_dim(), Lorene::Map::get_mg(), Lorene::Tbl::get_ndim(), Lorene::Mg3d::get_nr(), kerrschild, mp, Lorene::pow(), Lorene::Map::r, rad_ah(), runge_kutta_phi_bh(), runge_kutta_theta_bh(), Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_grid_point(), Lorene::Map::sint, Lorene::Scalar::std_spectral_base(), Lorene::Scalar::val_grid_point(), and Lorene::Scalar::val_point().

◆ mass_adm()

◆ mass_irr()

double Lorene::Black_hole::mass_irr ( ) const
virtual

◆ mass_kom()

◆ operator=()

void Lorene::Black_hole::operator= ( const Black_hole bh)

Assignment to another Black_hole.

Definition at line 239 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.

◆ operator>>()

ostream & Lorene::Black_hole::operator>> ( ostream &  ost) const
protectedvirtual

Operator >> (virtual function called by the operator <<)

Reimplemented in Lorene::Hole_bhns.

Definition at line 290 of file blackhole.C.

References confo, Lorene::Map::get_mg(), Lorene::Mg3d::get_nt(), kerrschild, lapconf, lapse, mass_adm(), mass_bh, mass_irr(), mass_kom(), mp, shift, and Lorene::Scalar::val_grid_point().

◆ r_coord()

const Scalar Lorene::Black_hole::r_coord ( bool  neumann,
bool  first 
) const

Expresses the areal radial coordinate by that in spatially isotropic coordinates.

Parameters
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 69 of file blackhole_r_coord.C.

References Lorene::Tensor::annule_domain(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::log(), mass_bh, mp, Lorene::Map::r, Lorene::Scalar::raccord(), Lorene::Scalar::set_grid_point(), Lorene::sqrt(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().

◆ rad_ah()

double Lorene::Black_hole::rad_ah ( ) const
virtual

Radius of the apparent horizon.

Definition at line 508 of file blackhole_global.C.

References mp, p_rad_ah, Lorene::Map::r, Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().

◆ rah_iso()

double Lorene::Black_hole::rah_iso ( bool  neumann,
bool  first 
) const

Computes a radius of apparent horizon (excised surface) in isotropic coordinates.

Parameters
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 72 of file blackhole_rah_iso.C.

References Lorene::exp(), and Lorene::sqrt().

◆ runge_kutta_phi_bh()

Tbl Lorene::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.

Parameters
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 70 of file blackhole_rk_phi.C.

References confo, Lorene::Scalar::dsdt(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), kerrschild, Lorene::Scalar::lapang(), mp, rad_ah(), Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_point().

◆ runge_kutta_theta_bh()

Tbl Lorene::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.

Parameters
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 70 of file blackhole_rk_theta.C.

References confo, Lorene::Map::get_mg(), Lorene::Mg3d::get_nt(), kerrschild, Lorene::Scalar::lapang(), mp, rad_ah(), Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::std_spectral_base(), Lorene::Scalar::stdsdp(), and Lorene::Scalar::val_point().

◆ sauve()

void Lorene::Black_hole::sauve ( FILE *  fich) const
virtual

Save in a file.

Reimplemented in Lorene::Hole_bhns.

Definition at line 270 of file blackhole.C.

References confo, kerrschild, lapconf_rs, mass_bh, Lorene::Tensor::sauve(), Lorene::Scalar::sauve(), and shift_rs.

◆ set_der_0x0()

void Lorene::Black_hole::set_der_0x0 ( ) const
protected

Sets to 0x0 all the pointers on derived quantities.

Definition at line 221 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.

◆ set_mass_bh()

double& Lorene::Black_hole::set_mass_bh ( )
inline

Read/write of the gravitational mass of BH [{ m_unit}].

Definition at line 207 of file blackhole.h.

References mass_bh.

◆ set_mp()

Map& Lorene::Black_hole::set_mp ( )
inline

Read/write of the mapping.

Definition at line 204 of file blackhole.h.

References mp.

◆ spin_am_bh()

double Lorene::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

◆ static_bh()

void Lorene::Black_hole::static_bh ( bool  neumann,
bool  first 
)

Sets the metric quantities to a spherical, static blach-hole analytic solution.

Parameters
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 63 of file blackhole_static.C.

References Lorene::Tensor::annule_domain(), confo, Lorene::Map::cosp, Lorene::Map::cost, Lorene::Map::get_bvect_cart(), kerrschild, lapconf, lapse, mass_bh, mp, Lorene::pow(), Lorene::Map::r, r_coord(), Lorene::Scalar::raccord(), Lorene::Vector::set(), Lorene::Tensor::set(), Lorene::Tensor::set_etat_qcq(), shift, Lorene::Map::sinp, Lorene::Map::sint, Lorene::sqrt(), Lorene::Vector::std_spectral_base(), and Lorene::Scalar::std_spectral_base().

Friends And Related Function Documentation

◆ operator<<

ostream& operator<< ( ostream &  ost,
const Black_hole bh 
)
friend

Display.

Definition at line 283 of file blackhole.C.

Member Data Documentation

◆ confo

Scalar Lorene::Black_hole::confo
protected

Conformal factor generated by the black hole.

Definition at line 118 of file blackhole.h.

◆ flat

Metric_flat Lorene::Black_hole::flat
protected

Flat metric defined on the mapping (Spherical components with respect to the mapping of the black hole ).

Definition at line 143 of file blackhole.h.

◆ kerrschild

bool Lorene::Black_hole::kerrschild
protected

true for a Kerr-Schild background, false for a conformally flat background

Definition at line 85 of file blackhole.h.

◆ lapconf

Scalar Lorene::Black_hole::lapconf
protected

A function (lapse function * conformal factor) lapconf generated by the black hole.

Definition at line 97 of file blackhole.h.

◆ lapconf_bh

Scalar Lorene::Black_hole::lapconf_bh
protected

Part of lapconf from the analytic background.

Definition at line 103 of file blackhole.h.

◆ lapconf_rs

Scalar Lorene::Black_hole::lapconf_rs
protected

Part of lapconf from the numerical computation.

Definition at line 100 of file blackhole.h.

◆ lapse

Scalar Lorene::Black_hole::lapse
protected

Lapse function generated by the black hole.

Definition at line 106 of file blackhole.h.

◆ mass_bh

double Lorene::Black_hole::mass_bh
protected

Gravitational mass of BH.

Definition at line 88 of file blackhole.h.

◆ mp

Map& Lorene::Black_hole::mp
protected

Mapping associated with the black hole.

Definition at line 80 of file blackhole.h.

◆ p_angu_mom_bh

Tbl* Lorene::Black_hole::p_angu_mom_bh
mutableprotected

Spin angular momentum.

Total angular momentum of the black hole

Definition at line 162 of file blackhole.h.

◆ p_mass_adm

double* Lorene::Black_hole::p_mass_adm
mutableprotected

Irreducible mass of the black hole.

Definition at line 153 of file blackhole.h.

◆ p_mass_irr

double* Lorene::Black_hole::p_mass_irr
mutableprotected

Conformal metric $\tilde \gamma_{ij}$.

Definition at line 151 of file blackhole.h.

◆ p_mass_kom

double* Lorene::Black_hole::p_mass_kom
mutableprotected

ADM mass.

Definition at line 155 of file blackhole.h.

◆ p_rad_ah

double* Lorene::Black_hole::p_rad_ah
mutableprotected

Komar mass.

Definition at line 157 of file blackhole.h.

◆ p_spin_am_bh

double* Lorene::Black_hole::p_spin_am_bh
mutableprotected

Radius of the apparent horizon.

Definition at line 159 of file blackhole.h.

◆ shift

Vector Lorene::Black_hole::shift
protected

Shift vector generated by the black hole.

Definition at line 109 of file blackhole.h.

◆ shift_bh

Vector Lorene::Black_hole::shift_bh
protected

Part of the shift vector from the analytic background.

Definition at line 115 of file blackhole.h.

◆ shift_rs

Vector Lorene::Black_hole::shift_rs
protected

Part of the shift vector from the numerical computation.

Definition at line 112 of file blackhole.h.

◆ taij

Sym_tensor Lorene::Black_hole::taij
protected

Trace of the extrinsic curvature.

Extrinsic curvature tensor $\ A^{ij}$ generated by shift , lapse , and confo .

Definition at line 127 of file blackhole.h.

◆ taij_quad

Scalar Lorene::Black_hole::taij_quad
protected

Part of the scalar $\eta_{ik} \eta_{jl} A^{ij} A^{kl}$ generated by $A_{ij}$.

Definition at line 135 of file blackhole.h.

◆ taij_quad_rs

Scalar Lorene::Black_hole::taij_quad_rs
protected

Part of the scalar.

Definition at line 138 of file blackhole.h.

◆ taij_rs

Sym_tensor Lorene::Black_hole::taij_rs
protected

Part of the extrinsic curvature tensor.

Definition at line 130 of file blackhole.h.


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