LORENE
|
Class to compute single black hole spacetime excised slices. More...
#include <excised_slice.h>
Public Member Functions | |
Excised_slice (const Map &mp_i, int hor_type, int scheme_type) | |
Standard constructor. More... | |
Excised_slice (const Excised_slice &) | |
Copy constructor. More... | |
Excised_slice (const Map &mp_i, int hor_type, int scheme_type, FILE *fich) | |
Constructor from a file (see sauve(FILE* ) ). More... | |
virtual | ~Excised_slice () |
Destructor. More... | |
void | operator= (const Excised_slice &) |
Assignment to another Excised_slice . More... | |
void | compute_stat_metric (double precis, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i=true, double relax=0.2, int mer_max=2000, int mer_max2=200, bool isvoid=true) |
If hor_type=1, computes a quasi-stationary 3-slice from the chosen parameters. More... | |
void | secmembre_kerr (Sym_tensor &source_hh) |
If hor_type=1, computes the rhs of hyperbolic equation for conformal metric assuming stationarity; WARNING; up to now, we are only able to handle void spacetimes. More... | |
const Map & | get_mp () const |
Returns the mapping. More... | |
int | get_type_hor () const |
Returns the type of horizon that performs the excision. More... | |
int | get_field_set () const |
Returns the field set chosen for the data. 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 (or its TT part if scheme_type=2) More... | |
const Vector & | get_Xx () const |
Return the longitudinal part of Einstein Equations if the scheme is appropriate. More... | |
int | set_type_hor () |
Sets the horizon type. More... | |
Scalar & | set_lapse () |
Sets the lapse. More... | |
Scalar & | set_conf_fact () |
Sets the conformal factor. More... | |
Vector & | set_shift () |
Sets the shift vector. More... | |
Sym_tensor & | set_hij () |
Sets the deviation tensor. More... | |
Sym_tensor & | set_hatA () |
Sets the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2) More... | |
Vector & | set_Xx () |
Sets the longitudinal part of Einstein Equations if the scheme is apropriate. 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 slice. More... | |
int | type_hor |
Chosen horizon type. More... | |
int | field_set |
Chose field set type. 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: rescaled same way as Cordero et al. More... | |
Vector | Xx |
Longitudinal part of the extrinsic curvature. 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 single black hole spacetime excised slices.
It takes as arguments:
mp
Definition at line 41 of file excised_slice.h.
Lorene::Excised_slice::Excised_slice | ( | const Map & | mp_i, |
int | hor_type, | ||
int | scheme_type | ||
) |
Standard constructor.
mp_i | Mapping on which the black hole slice will be defined |
hor_type | Type of horizon used for excision:
|
scheme_type | Type of scheme use for metric field computation
|
Definition at line 58 of file excised_slice.C.
References hatA, hij, lapse, set_der_0x0(), Lorene::Tensor::set_etat_zero(), shift, and Xx.
Lorene::Excised_slice::Excised_slice | ( | const Excised_slice & | ih | ) |
Lorene::Excised_slice::Excised_slice | ( | const Map & | mp_i, |
int | hor_type, | ||
int | scheme_type, | ||
FILE * | fich | ||
) |
Constructor from a file (see sauve(FILE* )
).
mp_i | Mapping on which the black hole slice will be defined |
hor_type | Type of horizon used for excision:
|
scheme_type | Type of scheme use for metric field computation
|
fich | input file (must have been created by the function sauve ) |
Definition at line 101 of file excised_slice.C.
References set_der_0x0().
|
virtual |
double Lorene::Excised_slice::adm_mass | ( | ) |
Computation of the ADM mass of the BH spacetime.
Definition at line 288 of file excised_slice.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::Excised_slice::compute_stat_metric | ( | double | precis, |
double | Omega_i, | ||
bool | NorKappa_i, | ||
Scalar | NoK_i, | ||
bool | isCF_i = true , |
||
double | relax = 0.2 , |
||
int | mer_max = 2000 , |
||
int | mer_max2 = 200 , |
||
bool | isvoid = true |
||
) |
If hor_type=1, 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) |
Omega_i | Rotation parameter for the isolated horizon slice. The principal Killing vector associated to this quantity will be by default the $e_{}$ vector. |
NorKappa_i | boolean to tell whether boundary condition is put on the surface gravity (true) or the lapse (false) |
NoK_i | The corresponding Dirichlet value for this quantity. |
isCF | Tells whether the calculation is using the conformal flatness approximation (true) or not (false). |
relax | relaxation coefficient for all metric fields (0: no relaxation) |
mer_max | Maximum number of iterations allowed in the main loop |
mer_max | Maximum number of iterations allowed in the possible secondary loop for conformal geometry computation. |
isvoid | indicates whether or not the spacetime is vacuum. This is useless for now since matter terms are not handled in the class. |
Definition at line 38 of file Excised_slice/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(), type_hor, Lorene::Tensor::up_down(), Lorene::Scalar::val_grid_point(), and Lorene::Valeur::ylm().
|
protectedvirtual |
Deletes all the derived quantities.
Definition at line 135 of file excised_slice.C.
References p_adm_mass, p_hor, p_komar_angmom, p_virial_residue, and set_der_0x0().
void Lorene::Excised_slice::Einstein_errors | ( | ) |
Prints out errors in Einstein equations for the data obtained.
Definition at line 193 of file excised_slice.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 conformal factor.
Definition at line 210 of file excised_slice.h.
|
inline |
Returns the field set chosen for the data.
Definition at line 204 of file excised_slice.h.
References field_set.
|
inline |
Returns the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2)
Definition at line 219 of file excised_slice.h.
References hatA.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Returns the type of horizon that performs the excision.
Definition at line 200 of file excised_slice.h.
References type_hor.
|
inline |
Return the longitudinal part of Einstein Equations if the scheme is appropriate.
Definition at line 222 of file excised_slice.h.
Spheroid Lorene::Excised_slice::hor | ( | ) |
Spheroid at the excised boundary associated with the black hole MOTS on the slice.
Set by default at the position .
Definition at line 246 of file excised_slice.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::Excised_slice::komar_angmom | ( | ) |
Computation of the Komar angular momentum w.r.t.
assumed rotational symmetry
Definition at line 308 of file excised_slice.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::Excised_slice::operator= | ( | const Excised_slice & | ih | ) |
Assignment to another Excised_slice
.
Definition at line 159 of file excised_slice.C.
References del_deriv(), field_set, hatA, hij, lapse, mp, shift, type_hor, and Xx.
|
virtual |
Save in a file.
Definition at line 181 of file excised_slice.C.
References hatA, hij, lapse, Lorene::Tensor::sauve(), Lorene::Scalar::sauve(), Lorene::Tensor_sym::sauve(), shift, and Xx.
void Lorene::Excised_slice::secmembre_kerr | ( | Sym_tensor & | source_hh | ) |
If hor_type=1, computes the rhs of hyperbolic equation for conformal metric assuming stationarity; WARNING; up to now, we are only able to handle void spacetimes.
Definition at line 12 of file Excised_slice/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().
|
inline |
|
protectedvirtual |
Sets to 0x0
all the pointers on derived quantities.
Definition at line 144 of file excised_slice.C.
References p_adm_mass, p_hor, p_komar_angmom, and p_virial_residue.
|
inline |
Sets the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2)
Definition at line 243 of file excised_slice.h.
References del_deriv(), and hatA.
|
inline |
Sets the deviation tensor.
Definition at line 240 of file excised_slice.h.
References del_deriv(), and hij.
|
inline |
|
inline |
Sets the shift vector.
Definition at line 237 of file excised_slice.h.
References del_deriv(), and shift.
|
inline |
Sets the horizon type.
Definition at line 228 of file excised_slice.h.
References del_deriv(), and type_hor.
|
inline |
Sets the longitudinal part of Einstein Equations if the scheme is apropriate.
Definition at line 246 of file excised_slice.h.
References del_deriv(), and Xx.
double Lorene::Excised_slice::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 333 of file excised_slice.C.
References Lorene::Scalar::asymptot(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), lapse, and mp.
|
protected |
Chose field set type.
Definition at line 55 of file excised_slice.h.
|
protected |
Rescaled tracefree extrinsic curvature tensor: rescaled same way as Cordero et al.
( 2009). In the New FCF scheme, this member will represent the TT part of the extrinsic curvature only.
Definition at line 74 of file excised_slice.h.
|
protected |
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al.
(2004)).
Definition at line 69 of file excised_slice.h.
|
protected |
Lapse function.
Definition at line 58 of file excised_slice.h.
|
protected |
Mapping associated with the slice.
Definition at line 46 of file excised_slice.h.
|
mutableprotected |
Computation of the ADM mass of the BH spacetime.
Definition at line 91 of file excised_slice.h.
|
mutableprotected |
Computation of the spheroid associated with the black hole horizon.
Definition at line 89 of file excised_slice.h.
|
mutableprotected |
Computation of the Komar angular momentum w.r.t.
assumed rotational symmetry
Definition at line 96 of file excised_slice.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 101 of file excised_slice.h.
|
protected |
Shift vector.
Definition at line 64 of file excised_slice.h.
|
protected |
Chosen horizon type.
Definition at line 52 of file excised_slice.h.
|
protected |
Longitudinal part of the extrinsic curvature.
Set to zero if we are in the original FCF scheme.
Definition at line 79 of file excised_slice.h.