19 #include "excision_surf.h" 20 #include "excision_hor.h" 135 Isol_hole(
const Map& mp_i,
double Omega_i,
bool NorKappa_i,
Scalar NoK_i,
bool isCF_i =
false) ;
154 Isol_hole(
const Map& mp_i,
double Omega_i,
bool NorKappa_i,
Scalar NoK_i,
bool isCF_i, FILE* fich) ;
181 void compute_stat_metric(
double precis,
double relax,
int mer_max,
int mer_max2,
bool isvoid =
true) ;
206 else cout <<
"The boundary condition is imposed on the surface gravity!" << endl;
214 else cout <<
"The boundary condition is imposed on the lapse!" <<endl;
239 virtual void sauve(FILE* )
const ;
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.
const Map & get_mp() const
Returns the mapping.
double * p_adm_mass
Computation of the ADM mass of the BH spacetime.
Class to compute quasistationary single black hole spacetimes in vacuum.
double komar_angmom()
Computation of the Komar angular momentum w.r.t.
double get_Omega() const
Returns the rotation rate.
const Scalar & get_boundN() const
Returns the boundary value used for the lapse (if it is the one used)
Tensor field of valence 0 (or component of a tensorial field).
Base class for coordinate mappings.
virtual void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
void secmembre_kerr(Sym_tensor &source_hh)
Computes the rhs of hyperbolic equation for conformal metric assuming statioarity; WARNING; up to now...
Sym_tensor hij
Deviation tensor( non-flat part of the conformal 3-metric on the slice; see Bonazzola et al...
virtual void del_deriv() const
Deletes all the derived quantities.
Tensor field of valence 1.
double * p_komar_angmom
Computation of the Komar angular momentum w.r.t.
const Vector & get_shift() const
Returns the shift vector .
Spheroidal 2-surfaces embedded in a time-slice of the 3+1 formalism.
Scalar boundNoK
Indicates if the boundary value for the lapse or the surface gravity is used.
void Einstein_errors()
Prints out errors in Einstein equations for the data obtained.
const Scalar & get_Kappa() const
Returns the surface gravity value at the boundary (if it is the one used)
Vector shift
Shift vector.
const Scalar & get_lapse() const
Returns the lapse function N.
const Map & mp
Mapping associated with the star.
Spheroid * p_hor
Computation of the spheroid associated with the black hole horizon.
double * p_virial_residue
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
const Sym_tensor & get_hatA() const
Returns the rescaled tracefree extrinsic curvature .
Spheroid hor()
Spheroid at the excised boundary associated with the black hole MOTS on the slice.
bool NorKappa
Rotation rate of the horizon in the azimuthal direction.
Scalar lapse
Lapse function.
double adm_mass()
Computation of the ADM mass of the BH spacetime.
double virial_residue()
Computation of the Virial residual, as difference at infinity between the ADM mass and the Komar inte...
void operator=(const Isol_hole &)
Assignment to another Isol_hole.
virtual void sauve(FILE *) const
Save in a file.
const Sym_tensor & get_hij() const
Returns the deviation tensor .
virtual ~Isol_hole()
Destructor.
bool isCF
Stores the boundary value of the lapse or surface gravity.
Isol_hole(const Map &mp_i, double Omega_i, bool NorKappa_i, Scalar NoK_i, bool isCF_i=false)
Standard constructor.
Class intended to describe valence-2 symmetric tensors.
Sym_tensor hatA
Rescaled tracefree extrinsic curvature tensor: see Cordero et al.
const Scalar & get_conf_fact() const
Returns the conformal factor.