11 #ifndef __EXCISEDSLICE_H_ 12 #define __EXCISEDSLICE_H_ 183 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) ;
224 else if (
field_set==1) cout <<
"Error: the scheme used is the original FCF; this variable is irrelevant" << endl;
225 else cout <<
"error in the scheme definition; please check the class consistency" << endl;}
252 virtual void sauve(FILE* )
const ;
Vector shift
Shift vector.
Sym_tensor hatA
Rescaled tracefree extrinsic curvature tensor: rescaled same way as Cordero et al.
void Einstein_errors()
Prints out errors in Einstein equations for the data obtained.
const Scalar & get_conf_fact() const
Returns the conformal factor.
Vector & set_shift()
Sets the shift vector.
Tensor field of valence 0 (or component of a tensorial field).
Base class for coordinate mappings.
const Map & get_mp() const
Returns the mapping.
virtual void sauve(FILE *) const
Save in a file.
Class to compute single black hole spacetime excised slices.
Tensor field of valence 1.
int field_set
Chose field set type.
void secmembre_kerr(Sym_tensor &source_hh)
If hor_type=1, computes the rhs of hyperbolic equation for conformal metric assuming stationarity; WA...
double * p_adm_mass
Computation of the ADM mass of the BH spacetime.
void operator=(const Excised_slice &)
Assignment to another Excised_slice.
virtual void del_deriv() const
Deletes all the derived quantities.
const Vector & get_shift() const
Returns the shift vector .
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
double * p_virial_residue
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
double komar_angmom()
Computation of the Komar angular momentum w.r.t.
Vector & set_Xx()
Sets the longitudinal part of Einstein Equations if the scheme is apropriate.
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.
Vector Xx
Longitudinal part of the extrinsic curvature.
Sym_tensor hij
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al...
Spheroid * p_hor
Computation of the spheroid associated with the black hole horizon.
const Vector & get_Xx() const
Return the longitudinal part of Einstein Equations if the scheme is appropriate.
const Map & mp
Mapping associated with the slice.
Spheroid hor()
Spheroid at the excised boundary associated with the black hole MOTS on the slice.
int type_hor
Chosen horizon type.
Scalar lapse
Lapse function.
Sym_tensor & set_hij()
Sets the deviation tensor.
int get_type_hor() const
Returns the type of horizon that performs the excision.
Sym_tensor & set_hatA()
Sets the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2) ...
double virial_residue()
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
Excised_slice(const Map &mp_i, int hor_type, int scheme_type)
Standard constructor.
double adm_mass()
Computation of the ADM mass of the BH spacetime.
Scalar & set_conf_fact()
Sets the conformal factor.
int get_field_set() const
Returns the field set chosen for the data.
const Scalar & get_lapse() const
Returns the lapse function N.
virtual ~Excised_slice()
Destructor.
Scalar & set_lapse()
Sets the lapse.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
int set_type_hor()
Sets the horizon type.
const Sym_tensor & get_hatA() const
Returns the rescaled tracefree extrinsic curvature (or its TT part if scheme_type=2) ...
double * p_komar_angmom
Computation of the Komar angular momentum w.r.t.
Class intended to describe valence-2 symmetric tensors.
const Sym_tensor & get_hij() const
Returns the deviation tensor .