LORENE
Lorene::Isol_hole Class Reference

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 Mapget_mp () const
 Returns the mapping. More...
 
double get_Omega () const
 Returns the rotation rate. More...
 
const Scalarget_boundN () const
 Returns the boundary value used for the lapse (if it is the one used) More...
 
const Scalarget_Kappa () const
 Returns the surface gravity value at the boundary (if it is the one used) More...
 
const Scalarget_lapse () const
 Returns the lapse function N. More...
 
const Scalarget_conf_fact () const
 Returns the conformal factor. More...
 
const Vectorget_shift () const
 Returns the shift vector $\beta^i$. More...
 
const Sym_tensorget_hij () const
 Returns the deviation tensor $ h^{ij} $. More...
 
const Sym_tensorget_hatA () const
 Returns the rescaled tracefree extrinsic curvature $\hat{A}^{ij}$. 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 Mapmp
 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...
 
Spheroidp_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...
 

Detailed Description

Class to compute quasistationary single black hole spacetimes in vacuum.

It takes as arguments:

  • A Mapping mp
  • A rotation rate $ \Omega$ of the horizon in the $ \varphi $ direction, assumed to be the direction of the rotational symmetry.
  • A Scalar NoK, storing the boundary value to be verified either for the lapse or the surface gravity at the excised boundary (which is a MOTS in the philosophy of this class). Whether it corresponds to the lapse or the surface gravity is controlled by the boolean NorKappa. (Warning: setting a boundary on Kappa is not implemented yet).
  • 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.

Constructor & Destructor Documentation

◆ Isol_hole() [1/3]

Lorene::Isol_hole::Isol_hole ( const Map mp_i,
double  Omega_i,
bool  NorKappa_i,
Scalar  NoK_i,
bool  isCF_i = false 
)

Standard constructor.

Parameters
mp_iMapping on which the black hole slice will be defined
Omega_irotation rate of the horizon
NorKappa_iFALSE: use of boundary value for the lapse TRUE use of boundary value for the surface gravity (Warning! not implemented yet!)
NoK_ivalue either for the lapse or surface gravity.
isCF_iFALSE: 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.

◆ Isol_hole() [2/3]

Lorene::Isol_hole::Isol_hole ( const Isol_hole ih)

Copy constructor.

Definition at line 83 of file isol_hole.C.

References set_der_0x0().

◆ Isol_hole() [3/3]

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* )).

Parameters
mp_iMapping on which the star will be defined
Omega_irotation rate of the horizon
NorKappa_iFALSE: use of boundary value for the lapse TRUE use of boundary value for the surface gravity (Warning! not implemented yet!)
NoK_ivalue either for the lapse or surface gravity.
isCF_iFALSE: Full GR 3+1 equations to be verified TRUE: IWM approximation used in determination of the spacetime geometry;
fichinput file (must have been created by the function sauve)

Definition at line 101 of file isol_hole.C.

References set_der_0x0().

◆ ~Isol_hole()

Lorene::Isol_hole::~Isol_hole ( )
virtual

Destructor.

Definition at line 125 of file isol_hole.C.

References del_deriv().

Member Function Documentation

◆ adm_mass()

double Lorene::Isol_hole::adm_mass ( )

◆ compute_stat_metric()

void Lorene::Isol_hole::compute_stat_metric ( double  precis,
double  relax,
int  mer_max,
int  mer_max2,
bool  isvoid = true 
)

◆ del_deriv()

void Lorene::Isol_hole::del_deriv ( ) const
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().

◆ Einstein_errors()

◆ get_boundN()

const Scalar& Lorene::Isol_hole::get_boundN ( ) const
inline

Returns the boundary value used for the lapse (if it is the one used)

Definition at line 202 of file isol_hole.h.

References boundNoK, and NorKappa.

◆ get_conf_fact()

const Scalar& Lorene::Isol_hole::get_conf_fact ( ) const
inline

Returns the conformal factor.

Definition at line 222 of file isol_hole.h.

◆ get_hatA()

const Sym_tensor& Lorene::Isol_hole::get_hatA ( ) const
inline

Returns the rescaled tracefree extrinsic curvature $\hat{A}^{ij}$.

Definition at line 231 of file isol_hole.h.

References hatA.

◆ get_hij()

const Sym_tensor& Lorene::Isol_hole::get_hij ( ) const
inline

Returns the deviation tensor $ h^{ij} $.

Definition at line 228 of file isol_hole.h.

References hij.

◆ get_Kappa()

const Scalar& Lorene::Isol_hole::get_Kappa ( ) const
inline

Returns the surface gravity value at the boundary (if it is the one used)

Definition at line 210 of file isol_hole.h.

References boundNoK, and NorKappa.

◆ get_lapse()

const Scalar& Lorene::Isol_hole::get_lapse ( ) const
inline

Returns the lapse function N.

Definition at line 219 of file isol_hole.h.

References lapse.

◆ get_mp()

const Map& Lorene::Isol_hole::get_mp ( ) const
inline

Returns the mapping.

Definition at line 196 of file isol_hole.h.

References mp.

◆ get_Omega()

double Lorene::Isol_hole::get_Omega ( ) const
inline

Returns the rotation rate.

Definition at line 199 of file isol_hole.h.

◆ get_shift()

const Vector& Lorene::Isol_hole::get_shift ( ) const
inline

Returns the shift vector $\beta^i$.

Definition at line 225 of file isol_hole.h.

References shift.

◆ hor()

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 $ r=1 $.

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().

◆ komar_angmom()

double Lorene::Isol_hole::komar_angmom ( )

◆ operator=()

void Lorene::Isol_hole::operator= ( const Isol_hole ih)

Assignment to another Isol_hole.

Definition at line 160 of file isol_hole.C.

References boundNoK, del_deriv(), hatA, hij, isCF, lapse, mp, NorKappa, and shift.

◆ sauve()

void Lorene::Isol_hole::sauve ( FILE *  fich) const
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.

◆ secmembre_kerr()

◆ set_der_0x0()

void Lorene::Isol_hole::set_der_0x0 ( ) const
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.

◆ 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.

Member Data Documentation

◆ boundNoK

Scalar Lorene::Isol_hole::boundNoK
protected

Indicates if the boundary value for the lapse or the surface gravity is used.

Definition at line 65 of file isol_hole.h.

◆ hatA

Sym_tensor Lorene::Isol_hole::hatA
protected

Rescaled tracefree extrinsic curvature tensor: see Cordero et al.

( * 2009).

Definition at line 92 of file isol_hole.h.

◆ hij

Sym_tensor Lorene::Isol_hole::hij
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.

◆ isCF

bool Lorene::Isol_hole::isCF
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.

◆ lapse

Scalar Lorene::Isol_hole::lapse
protected

Lapse function.

Definition at line 76 of file isol_hole.h.

◆ mp

const Map& Lorene::Isol_hole::mp
protected

Mapping associated with the star.

Definition at line 55 of file isol_hole.h.

◆ NorKappa

bool Lorene::Isol_hole::NorKappa
protected

Rotation rate of the horizon in the azimuthal direction.

Definition at line 61 of file isol_hole.h.

◆ p_adm_mass

double* Lorene::Isol_hole::p_adm_mass
mutableprotected

Computation of the ADM mass of the BH spacetime.

Definition at line 103 of file isol_hole.h.

◆ p_hor

Spheroid* Lorene::Isol_hole::p_hor
mutableprotected

Computation of the spheroid associated with the black hole horizon.

Definition at line 101 of file isol_hole.h.

◆ p_komar_angmom

double* Lorene::Isol_hole::p_komar_angmom
mutableprotected

Computation of the Komar angular momentum w.r.t.

assumed rotational symmetry

Definition at line 108 of file isol_hole.h.

◆ p_virial_residue

double* Lorene::Isol_hole::p_virial_residue
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.

◆ shift

Vector Lorene::Isol_hole::shift
protected

Shift vector.

Definition at line 82 of file isol_hole.h.


The documentation for this class was generated from the following files: