LORENE
|
Class to compute quasistationary single black hole spacetimes in vacuum. More...
#include <isol_hole.h>
Public Member Functions | |
Isol_hole (const Map &mp_i, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i=false) | |
Standard constructor. More... | |
Isol_hole (const Isol_hole &) | |
Copy constructor. More... | |
Isol_hole (const Map &mp_i, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i, FILE *fich) | |
Constructor from a file (see sauve(FILE* ) ). More... | |
virtual | ~Isol_hole () |
Destructor. More... | |
void | operator= (const Isol_hole &) |
Assignment to another Isol_hole . More... | |
void | compute_stat_metric (double precis, double relax, int mer_max, int mer_max2, bool isvoid=true) |
Computes a quasi-stationary 3-slice from the chosen parameters. More... | |
void | secmembre_kerr (Sym_tensor &source_hh) |
Computes the rhs of hyperbolic equation for conformal metric assuming statioarity; WARNING; up to now, we are only able to handle void spacetimes. More... | |
const Map & | get_mp () const |
Returns the mapping. More... | |
double | get_Omega () const |
Returns the rotation rate. More... | |
const Scalar & | get_boundN () const |
Returns the boundary value used for the lapse (if it is the one used) More... | |
const Scalar & | get_Kappa () const |
Returns the surface gravity value at the boundary (if it is the one used) More... | |
const Scalar & | get_lapse () const |
Returns the lapse function N. More... | |
const Scalar & | get_conf_fact () const |
Returns the conformal factor. More... | |
const Vector & | get_shift () const |
Returns the shift vector . More... | |
const Sym_tensor & | get_hij () const |
Returns the deviation tensor . More... | |
const Sym_tensor & | get_hatA () const |
Returns the rescaled tracefree extrinsic curvature . More... | |
virtual void | sauve (FILE *) const |
Save in a file. More... | |
void | Einstein_errors () |
Prints out errors in Einstein equations for the data obtained. More... | |
Spheroid | hor () |
Spheroid at the excised boundary associated with the black hole MOTS on the slice. More... | |
double | adm_mass () |
Computation of the ADM mass of the BH spacetime. More... | |
double | komar_angmom () |
Computation of the Komar angular momentum w.r.t. More... | |
double | virial_residue () |
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar integral associated to the mass. More... | |
Protected Member Functions | |
virtual void | del_deriv () const |
Deletes all the derived quantities. More... | |
virtual void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. More... | |
Protected Attributes | |
const Map & | mp |
Mapping associated with the star. More... | |
double | Omega |
bool | NorKappa |
Rotation rate of the horizon in the azimuthal direction. More... | |
Scalar | boundNoK |
Indicates if the boundary value for the lapse or the surface gravity is used. More... | |
bool | isCF |
Stores the boundary value of the lapse or surface gravity. More... | |
Scalar | lapse |
Lapse function. More... | |
Scalar | conf_fact |
Vector | shift |
Shift vector. More... | |
Sym_tensor | hij |
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al. More... | |
Sym_tensor | hatA |
Rescaled tracefree extrinsic curvature tensor: see Cordero et al. More... | |
Spheroid * | p_hor |
Computation of the spheroid associated with the black hole horizon. More... | |
double * | p_adm_mass |
Computation of the ADM mass of the BH spacetime. More... | |
double * | p_komar_angmom |
Computation of the Komar angular momentum w.r.t. More... | |
double * | p_virial_residue |
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar integral associated to the mass. More... | |
Class to compute quasistationary single black hole spacetimes in vacuum.
It takes as arguments:
mp
A boolean isCF, to determine whether or not one adopts the conformally flat approximation in the spacetime construction. Set to FALSE by default.
The main goal of this class is the computation of metric data from those parameters, as well as some global quantities related to the physical characteristics of the spacetime. Those metric data are given on a spacelike 3-slice, assuming a global timelike Killing field and using adapted coordinates.
Definition at line 50 of file isol_hole.h.
Lorene::Isol_hole::Isol_hole | ( | const Map & | mp_i, |
double | Omega_i, | ||
bool | NorKappa_i, | ||
Scalar | NoK_i, | ||
bool | isCF_i = false |
||
) |
Standard constructor.
mp_i | Mapping on which the black hole slice will be defined |
Omega_i | rotation rate of the horizon |
NorKappa_i | FALSE: use of boundary value for the lapse TRUE use of boundary value for the surface gravity (Warning! not implemented yet!) |
NoK_i | value either for the lapse or surface gravity. |
isCF_i | FALSE: Full GR 3+1 equations to be verified TRUE: IWM approximation used in determination of the spacetime geometry; |
Definition at line 51 of file isol_hole.C.
References boundNoK, Lorene::Tensor::get_mp(), hatA, hij, lapse, set_der_0x0(), Lorene::Tensor::set_etat_zero(), and shift.
Lorene::Isol_hole::Isol_hole | ( | const Isol_hole & | ih | ) |
Lorene::Isol_hole::Isol_hole | ( | const Map & | mp_i, |
double | Omega_i, | ||
bool | NorKappa_i, | ||
Scalar | NoK_i, | ||
bool | isCF_i, | ||
FILE * | fich | ||
) |
Constructor from a file (see sauve(FILE* )
).
mp_i | Mapping on which the star will be defined |
Omega_i | rotation rate of the horizon |
NorKappa_i | FALSE: use of boundary value for the lapse TRUE use of boundary value for the surface gravity (Warning! not implemented yet!) |
NoK_i | value either for the lapse or surface gravity. |
isCF_i | FALSE: Full GR 3+1 equations to be verified TRUE: IWM approximation used in determination of the spacetime geometry; |
fich | input file (must have been created by the function sauve ) |
Definition at line 101 of file isol_hole.C.
References set_der_0x0().
|
virtual |
double Lorene::Isol_hole::adm_mass | ( | ) |
Computation of the ADM mass of the BH spacetime.
Definition at line 289 of file isol_hole.C.
References Lorene::Metric_flat::con(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Tensor_sym::derive_con(), Lorene::Scalar::dsdr(), hij, mp, Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
void Lorene::Isol_hole::compute_stat_metric | ( | double | precis, |
double | relax, | ||
int | mer_max, | ||
int | mer_max2, | ||
bool | isvoid = true |
||
) |
Computes a quasi-stationary 3-slice from the chosen parameters.
precis | [input] threshold in the relative difference between the radial shift for two consecutive steps to stop the iterative procedure (default value: 1.e-11) |
Definition at line 50 of file Isol_hole/isol_hole_compute_metric.C.
References Lorene::Tensor::annule_domain(), Lorene::Scalar::annule_hard(), Lorene::Valeur::c_cf, Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Tensor_sym::derive_con(), Lorene::Scalar::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::Vector::divergence(), Lorene::Mg3d::get_angu_1dom(), Lorene::Scalar::get_etat(), Lorene::Mg3d::get_nzone(), hij, Lorene::Scalar::inc_dzpuis(), lapse, Lorene::log(), mp, Lorene::Vector::set(), Lorene::Tensor::set(), Lorene::Scalar::set_spectral_va(), Lorene::sqrt(), Lorene::Vector::std_spectral_base(), Lorene::Scalar::std_spectral_base(), Lorene::Tensor::up_down(), Lorene::Scalar::val_grid_point(), and Lorene::Valeur::ylm().
|
protectedvirtual |
Deletes all the derived quantities.
Definition at line 136 of file isol_hole.C.
References p_adm_mass, p_hor, p_komar_angmom, p_virial_residue, and set_der_0x0().
void Lorene::Isol_hole::Einstein_errors | ( | ) |
Prints out errors in Einstein equations for the data obtained.
Definition at line 194 of file isol_hole.C.
References Lorene::Metric_flat::con(), Lorene::contract(), Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::dec_dzpuis(), Lorene::Scalar::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Sym_tensor::derive_lie(), Lorene::Sym_tensor::divergence(), hatA, hij, lapse, Lorene::maxabs(), mp, Lorene::pow(), Lorene::Metric::ricci(), Lorene::Metric::ricci_scal(), shift, Lorene::Tensor::std_spectral_base(), Lorene::Tensor::trace(), Lorene::Tensor::up(), and Lorene::Tensor::up_down().
|
inline |
Returns the boundary value used for the lapse (if it is the one used)
Definition at line 202 of file isol_hole.h.
|
inline |
Returns the conformal factor.
Definition at line 222 of file isol_hole.h.
|
inline |
Returns the rescaled tracefree extrinsic curvature .
Definition at line 231 of file isol_hole.h.
References hatA.
|
inline |
|
inline |
Returns the surface gravity value at the boundary (if it is the one used)
Definition at line 210 of file isol_hole.h.
|
inline |
|
inline |
|
inline |
Returns the rotation rate.
Definition at line 199 of file isol_hole.h.
|
inline |
Spheroid Lorene::Isol_hole::hor | ( | ) |
Spheroid at the excised boundary associated with the black hole MOTS on the slice.
Set by default at the position .
Definition at line 247 of file isol_hole.C.
References Lorene::Metric_flat::con(), Lorene::Mg3d::get_angu_1dom(), hatA, hij, mp, Lorene::pow(), Lorene::Tensor::std_spectral_base(), Lorene::Scalar::std_spectral_base(), Lorene::Tensor::up_down(), and Lorene::Scalar::val_grid_point().
double Lorene::Isol_hole::komar_angmom | ( | ) |
Computation of the Komar angular momentum w.r.t.
assumed rotational symmetry
Definition at line 309 of file isol_hole.C.
References Lorene::Metric_flat::con(), Lorene::Metric::cov(), Lorene::Tensor::dec_dzpuis(), hatA, hij, mp, Lorene::Scalar::mult_r_dzpuis(), Lorene::Scalar::mult_sint(), Lorene::pow(), Lorene::sqrt(), Lorene::Tensor::std_spectral_base(), Lorene::Scalar::std_spectral_base(), and Lorene::Tensor::up_down().
void Lorene::Isol_hole::operator= | ( | const Isol_hole & | ih | ) |
|
virtual |
Save in a file.
Definition at line 183 of file isol_hole.C.
References hatA, hij, lapse, Lorene::Tensor::sauve(), Lorene::Scalar::sauve(), Lorene::Tensor_sym::sauve(), and shift.
void Lorene::Isol_hole::secmembre_kerr | ( | Sym_tensor & | source_hh | ) |
Computes the rhs of hyperbolic equation for conformal metric assuming statioarity; WARNING; up to now, we are only able to handle void spacetimes.
Definition at line 18 of file Isol_hole/secmembre_hij_stat.C.
References Lorene::Tensor::annule_domain(), Lorene::Scalar::annule_hard(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Scalar::derive_con(), Lorene::Tensor_sym::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::Sym_tensor::derive_lie(), Lorene::Scalar::derive_lie(), Lorene::Vector::divergence(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), hatA, hij, Lorene::Tenseur::inc_dzpuis(), Lorene::Tensor::inc_dzpuis(), Lorene::Scalar::inc_dzpuis(), lapse, Lorene::log(), Lorene::Vector::ope_killing_conf(), Lorene::Vector::set(), Lorene::Tensor::set(), shift, Lorene::sqrt(), Lorene::Vector::std_spectral_base(), Lorene::Tensor::std_spectral_base(), Lorene::Scalar::std_spectral_base(), and Lorene::Tensor::trace().
|
protectedvirtual |
Sets to 0x0
all the pointers on derived quantities.
Definition at line 145 of file isol_hole.C.
References p_adm_mass, p_hor, p_komar_angmom, and p_virial_residue.
double Lorene::Isol_hole::virial_residue | ( | ) |
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar integral associated to the mass.
Definition at line 334 of file isol_hole.C.
References Lorene::Scalar::asymptot(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), lapse, and mp.
|
protected |
Indicates if the boundary value for the lapse or the surface gravity is used.
Definition at line 65 of file isol_hole.h.
|
protected |
Rescaled tracefree extrinsic curvature tensor: see Cordero et al.
( * 2009).
Definition at line 92 of file isol_hole.h.
|
protected |
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al.
(2004)).
Definition at line 87 of file isol_hole.h.
|
protected |
Stores the boundary value of the lapse or surface gravity.
Indicates if the CF approximation is used.
Definition at line 69 of file isol_hole.h.
|
protected |
Lapse function.
Definition at line 76 of file isol_hole.h.
|
protected |
Mapping associated with the star.
Definition at line 55 of file isol_hole.h.
|
protected |
Rotation rate of the horizon in the azimuthal direction.
Definition at line 61 of file isol_hole.h.
|
mutableprotected |
Computation of the ADM mass of the BH spacetime.
Definition at line 103 of file isol_hole.h.
|
mutableprotected |
Computation of the spheroid associated with the black hole horizon.
Definition at line 101 of file isol_hole.h.
|
mutableprotected |
Computation of the Komar angular momentum w.r.t.
assumed rotational symmetry
Definition at line 108 of file isol_hole.h.
|
mutableprotected |
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar integral associated to the mass.
Definition at line 113 of file isol_hole.h.
|
protected |
Shift vector.
Definition at line 82 of file isol_hole.h.