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. | |
| Excised_slice (const Excised_slice &) | |
| Copy constructor. | |
| Excised_slice (const Map &mp_i, int hor_type, int scheme_type, FILE *fich) | |
Constructor from a file (see sauve(FILE* )). | |
| virtual | ~Excised_slice () |
| Destructor. | |
| void | operator= (const Excised_slice &) |
Assignment to another Excised_slice. | |
| 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. | |
| 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. | |
| const Map & | get_mp () const |
| Returns the mapping. | |
| int | get_type_hor () const |
| Returns the type of horizon that performs the excision. | |
| int | get_field_set () const |
| Returns the field set chosen for the data. | |
| const Scalar & | get_lapse () const |
| Returns the lapse function N. | |
| const Scalar & | get_conf_fact () const |
| Returns the conformal factor. | |
| const Vector & | get_shift () const |
Returns the shift vector . | |
| const Sym_tensor & | get_hij () const |
Returns the deviation tensor . | |
| const Sym_tensor & | get_hatA () const |
Returns the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2). | |
| const Vector & | get_Xx () const |
| Return the longitudinal part of Einstein Equations if the scheme is appropriate. | |
| int | set_type_hor () |
| Sets the horizon type. | |
| Scalar & | set_lapse () |
| Sets the lapse. | |
| Scalar & | set_conf_fact () |
| Sets the conformal factor. | |
| Vector & | set_shift () |
| Sets the shift vector. | |
| Sym_tensor & | set_hij () |
| Sets the deviation tensor. | |
| Sym_tensor & | set_hatA () |
Sets the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2). | |
| Vector & | set_Xx () |
| Sets the longitudinal part of Einstein Equations if the scheme is apropriate. | |
| virtual void | sauve (FILE *) const |
| Save in a file. | |
| void | Einstein_errors () |
| Prints out errors in Einstein equations for the data obtained. | |
| Spheroid | hor () |
| Spheroid at the excised boundary associated with the black hole MOTS on the slice. | |
| double | adm_mass () |
| Computation of the ADM mass of the BH spacetime. | |
| double | komar_angmom () |
| Computation of the Komar angular momentum w.r.t. | |
| double | virial_residue () |
| Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar integral associated to the mass. | |
Protected Member Functions | |
| virtual void | del_deriv () const |
| Deletes all the derived quantities. | |
| virtual void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. | |
Protected Attributes | |
| const Map & | mp |
| Mapping associated with the slice. | |
| int | type_hor |
| Chosen horizon type. | |
| int | field_set |
| Chose field set type. | |
| Scalar | lapse |
| Lapse function. | |
| Scalar | conf_fact |
| Vector | shift |
| Shift vector. | |
| Sym_tensor | hij |
| Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al. | |
| Sym_tensor | hatA |
| Rescaled tracefree extrinsic curvature tensor: rescaled same way as Cordero et al. | |
| Vector | Xx |
| Longitudinal part of the extrinsic curvature. | |
| Spheroid * | p_hor |
| Computation of the spheroid associated with the black hole horizon. | |
| double * | p_adm_mass |
| Computation of the ADM mass of the BH spacetime. | |
| double * | p_komar_angmom |
| Computation of the Komar angular momentum w.r.t. | |
| 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. | |
Class to compute single black hole spacetime excised slices.
It takes as arguments:
mp Definition at line 40 of file excised_slice.h.
| 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(), Tensor::set_etat_zero(), shift, and Xx.
| Excised_slice::Excised_slice | ( | const Excised_slice & | ih | ) |
| 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().
| Excised_slice::~Excised_slice | ( | ) | [virtual] |
| double Excised_slice::adm_mass | ( | ) |
Computation of the ADM mass of the BH spacetime.
Definition at line 288 of file excised_slice.C.
References Metric_flat::con(), contract(), Metric::cov(), Tensor_sym::derive_con(), Scalar::dsdr(), Unites::ggrav, hij, mp, sqrt(), and Scalar::std_spectral_base().
| void 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 37 of file Excised_slice/isol_hole_compute_metric.C.
References Scalar::allocate_all(), Tensor::annule(), Tensor::annule_domain(), Scalar::annule_hard(), Valeur::c_cf, Valeur::coef(), Metric::con(), Metric_flat::con(), contract(), Metric::cov(), Scalar::dec_dzpuis(), Tensor::dec_dzpuis(), Scalar::derive_con(), Tensor::derive_con(), Tensor_sym::derive_con(), Tensor::derive_cov(), Scalar::derive_cov(), Tensor_sym::derive_cov(), Sym_tensor::derive_lie(), Metric::determinant(), Sym_tensor::divergence(), Vector::divergence(), Scalar::dsdr(), Vector::eta(), Scalar::get_etat(), Scalar::get_spectral_va(), hatA, hij, Tensor::inc_dzpuis(), Scalar::inc_dzpuis(), Scalar::laplacian(), lapse, log(), max(), maxabs(), mp, Vector::mu(), Scalar::mult_rsint(), secmembre_kerr(), Vector::set(), Tensor::set(), Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_grid_point(), Scalar::set_spectral_va(), shift, Scalar::sol_elliptic_boundary(), sqrt(), Vector::std_spectral_base(), Scalar::std_spectral_base(), Tensor::trace(), type_hor, Tensor::up(), Tensor::up_down(), Scalar::val_grid_point(), and Valeur::ylm().
| void Excised_slice::del_deriv | ( | ) | const [protected, virtual] |
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 Excised_slice::Einstein_errors | ( | ) |
Prints out errors in Einstein equations for the data obtained.
Definition at line 193 of file excised_slice.C.
References Metric_flat::con(), contract(), Scalar::dec_dzpuis(), Tensor::dec_dzpuis(), Scalar::derive_con(), Tensor::derive_cov(), Scalar::derive_cov(), Sym_tensor::derive_lie(), Sym_tensor::divergence(), hatA, hij, lapse, maxabs(), mp, pow(), Metric::ricci(), Metric::ricci_scal(), shift, Tensor::std_spectral_base(), Tensor::trace(), Tensor::up(), and Tensor::up_down().
| const Scalar& Excised_slice::get_conf_fact | ( | ) | const [inline] |
Returns the conformal factor.
Definition at line 209 of file excised_slice.h.
| int Excised_slice::get_field_set | ( | ) | const [inline] |
Returns the field set chosen for the data.
Definition at line 203 of file excised_slice.h.
References field_set.
| const Sym_tensor& Excised_slice::get_hatA | ( | ) | const [inline] |
Returns the rescaled tracefree extrinsic curvature
(or its TT part if scheme_type=2).
Definition at line 218 of file excised_slice.h.
References hatA.
| const Sym_tensor& Excised_slice::get_hij | ( | ) | const [inline] |
| const Scalar& Excised_slice::get_lapse | ( | ) | const [inline] |
| const Map& Excised_slice::get_mp | ( | ) | const [inline] |
| const Vector& Excised_slice::get_shift | ( | ) | const [inline] |
| int Excised_slice::get_type_hor | ( | ) | const [inline] |
Returns the type of horizon that performs the excision.
Definition at line 199 of file excised_slice.h.
References type_hor.
| const Vector& Excised_slice::get_Xx | ( | ) | const [inline] |
Return the longitudinal part of Einstein Equations if the scheme is appropriate.
Definition at line 221 of file excised_slice.h.
| Spheroid 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 Metric_flat::con(), hatA, hij, mp, pow(), Tensor::std_spectral_base(), Scalar::std_spectral_base(), Tensor::up_down(), and Scalar::val_grid_point().
| double 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 Metric_flat::con(), Metric::cov(), Tensor::dec_dzpuis(), Unites::ggrav, hatA, hij, mp, Scalar::mult_r_dzpuis(), Scalar::mult_sint(), pow(), sqrt(), Scalar::std_spectral_base(), Tensor::std_spectral_base(), and Tensor::up_down().
| void 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.
| void Excised_slice::sauve | ( | FILE * | fich | ) | const [virtual] |
Save in a file.
Definition at line 181 of file excised_slice.C.
References hatA, hij, lapse, Tensor_sym::sauve(), Tensor::sauve(), Scalar::sauve(), shift, and Xx.
| void 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 11 of file Excised_slice/secmembre_hij_stat.C.
References Tensor::annule_domain(), Scalar::annule_hard(), Metric_flat::con(), contract(), Metric::cov(), Tensor_sym::derive_con(), Scalar::derive_con(), Tensor::derive_cov(), Scalar::derive_cov(), Tensor_sym::derive_cov(), Scalar::derive_lie(), Sym_tensor::derive_lie(), Vector::divergence(), Map::get_bvect_spher(), Map::get_mg(), Tensor::get_mp(), hatA, hij, Tensor::inc_dzpuis(), Scalar::inc_dzpuis(), Tenseur::inc_dzpuis(), lapse, log(), Vector::ope_killing_conf(), Vector::set(), Tensor::set(), shift, sqrt(), Vector::std_spectral_base(), Tensor::std_spectral_base(), Scalar::std_spectral_base(), and Tensor::trace().
| Scalar& Excised_slice::set_conf_fact | ( | ) | [inline] |
| void Excised_slice::set_der_0x0 | ( | ) | const [protected, virtual] |
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.
| Sym_tensor& Excised_slice::set_hatA | ( | ) | [inline] |
Sets the rescaled tracefree extrinsic curvature
(or its TT part if scheme_type=2).
Definition at line 242 of file excised_slice.h.
References del_deriv(), and hatA.
| Sym_tensor& Excised_slice::set_hij | ( | ) | [inline] |
Sets the deviation tensor.
Definition at line 239 of file excised_slice.h.
References del_deriv(), and hij.
| Scalar& Excised_slice::set_lapse | ( | ) | [inline] |
| Vector& Excised_slice::set_shift | ( | ) | [inline] |
Sets the shift vector.
Definition at line 236 of file excised_slice.h.
References del_deriv(), and shift.
| int Excised_slice::set_type_hor | ( | ) | [inline] |
Sets the horizon type.
Definition at line 227 of file excised_slice.h.
References del_deriv(), and type_hor.
| Vector& Excised_slice::set_Xx | ( | ) | [inline] |
Sets the longitudinal part of Einstein Equations if the scheme is apropriate.
Definition at line 245 of file excised_slice.h.
References del_deriv(), and Xx.
| double 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 Scalar::asymptot(), Map::get_mg(), Mg3d::get_nzone(), lapse, and mp.
int Excised_slice::field_set [protected] |
Chose field set type.
Definition at line 54 of file excised_slice.h.
Sym_tensor Excised_slice::hatA [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 73 of file excised_slice.h.
Sym_tensor Excised_slice::hij [protected] |
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al.
(2004)).
Definition at line 68 of file excised_slice.h.
Scalar Excised_slice::lapse [protected] |
Lapse function.
Definition at line 57 of file excised_slice.h.
const Map& Excised_slice::mp [protected] |
Mapping associated with the slice.
Definition at line 45 of file excised_slice.h.
double* Excised_slice::p_adm_mass [mutable, protected] |
Computation of the ADM mass of the BH spacetime.
Definition at line 90 of file excised_slice.h.
Spheroid* Excised_slice::p_hor [mutable, protected] |
Computation of the spheroid associated with the black hole horizon.
Definition at line 88 of file excised_slice.h.
double* Excised_slice::p_komar_angmom [mutable, protected] |
Computation of the Komar angular momentum w.r.t.
assumed rotational symmetry
Definition at line 95 of file excised_slice.h.
double* Excised_slice::p_virial_residue [mutable, protected] |
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar integral associated to the mass.
Definition at line 100 of file excised_slice.h.
Vector Excised_slice::shift [protected] |
Shift vector.
Definition at line 63 of file excised_slice.h.
int Excised_slice::type_hor [protected] |
Chosen horizon type.
Definition at line 51 of file excised_slice.h.
Vector Excised_slice::Xx [protected] |
Longitudinal part of the extrinsic curvature.
Set to zero if we are in the original FCF scheme.
Definition at line 78 of file excised_slice.h.
1.6.1