Isol_hor Class Reference
[Time evolution (under development)]

Spacelike time-slice of an Isolated Horizon in a 3+1 spacetime with conformal decomposition. More...

#include <isol_hor.h>

Inheritance diagram for Isol_hor:
Time_slice_conf Time_slice

List of all members.

Public Member Functions

 Isol_hor (Map_af &mpi, int depth_in=3)
 Standard constructor.
 Isol_hor (Map_af &mpi, const Scalar &lapse_in, const Scalar &psi_in, const Vector &shift_in, const Sym_tensor &aa_in, const Metric &gamt, const Sym_tensor &gamt_point, const Scalar &trK, const Scalar &trK_point, const Metric_flat &ff_in, int depth_in=3)
 Constructor from conformal decomposition.
 Isol_hor (const Isol_hor &)
 Copy constructor.
 Isol_hor (Map_af &mp, FILE *fich, bool partial_read, int depth_in=3)
 Constructor from a binary file.
virtual ~Isol_hor ()
 Destructor.
void operator= (const Isol_hor &)
 Assignment to another Isol_hor.
const Map_afget_mp () const
 Returns the mapping (readonly).
Map_afset_mp ()
 Read/write of the mapping.
double get_radius () const
 Returns the radius of the horizon.
void set_radius (double rad)
 Sets the radius of the horizon to rad .
double get_omega () const
 Returns the angular velocity.
void set_omega (double ome)
 Sets the angular velocity to ome .
double get_boost_x () const
 Returns the boost velocity in x-direction.
void set_boost_x (double bo)
 Sets the boost velocity in x-direction to bo .
double get_boost_z () const
 Returns the boost velocity in z-direction.
void set_boost_z (double bo)
 Sets the boost velocity in z-direction to bo .
virtual const Scalarn_auto () const
 Lapse function $ N_{auto} $ at the current time step jtime.
virtual const Scalarn_comp () const
 Lapse function $ N_{comp} $ at the current time step jtime.
virtual const Scalarpsi_auto () const
 Conformal factor $ \Psi_{auto} $ at the current time step jtime.
virtual const Scalarpsi_comp () const
 Conformal factor $ \Psi_{comp} $ at the current time step jtime.
virtual const Vectordnn () const
 Covariant derivative of the lapse function $ \overline\nabla_i N $ at the current time step jtime.
virtual const Vectordpsi () const
 Covariant derivative with respect to the flat metric of the conformal factor $ \overline\nabla_i \Psi $ at the current time step jtime.
virtual const Vectorbeta_auto () const
 Shift function $ \beta^i_{auto} $ at the current time step jtime.
virtual const Vectorbeta_comp () const
 Shift function $ \beta^i_{comp} $ at the current time step jtime.
virtual const Sym_tensoraa_auto () const
 Conformal representation $ A^{ij}_{auto} $ of the traceless part of the extrinsic curvature: Returns the value at the current time step jtime.
virtual const Sym_tensoraa_comp () const
 Conformal representation $ A^{ij}_{comp} $ of the traceless part of the extrinsic curvature: Returns the value at the current time step jtime.
virtual const Scalaraa_quad () const
 Conformal representation $ A^{ij}A_{ij} $.
virtual const Metrictgam () const
 Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $ Returns the value at the current time step (jtime ).
const Scalar get_decouple () const
 Returns the function used to construct tkij_auto from tkij_tot .
void n_comp (const Isol_hor &comp)
 Imports the part of N due to the companion hole comp .
void psi_comp (const Isol_hor &comp)
 Imports the part of $\Psi$ due to the companion hole comp .
void beta_comp (const Isol_hor &comp)
 Imports the part of $\beta^i$ due to the companion hole comp.
double viriel_seul () const
 Computes the viriel error, that is the difference between the ADM and the Komar masses, calculated by the asymptotic behaviours of respectively $\Psi$ and N .
void init_bhole ()
 Sets the values of the fields to :.
void init_met_trK ()
 Sets the 3-metric tilde to the flat metric and gamt_point, trK and trK_point to zero.
void init_bhole_seul ()
 Initiates for a single black hole.
void set_psi (const Scalar &psi_in)
 Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
void set_nn (const Scalar &nn_in)
 Sets the lapse.
void set_gamt (const Metric &gam_tilde)
 Sets the conformal metric to gam_tilde.
const Vector radial_vect_hor () const
 Vector radial normal.
const Vector tradial_vect_hor () const
 Vector radial normal tilde.
const Scalar b_tilde () const
 Radial component of the shift with respect to the conformal metric.
const Scalar darea_hor () const
 Element of area of the horizon.
double area_hor () const
 Area of the horizon.
double radius_hor () const
 Radius of the horizon.
double ang_mom_hor () const
 Angular momentum (modulo).
double mass_hor () const
 Mass computed at the horizon.
double kappa_hor () const
 Surface gravity.
double omega_hor () const
 Orbital velocity.
double ang_mom_adm () const
 ADM angular Momentum.
Scalar expansion () const
 Expansion of the outgoing null normal ($ \bf n + \bf s $).
void init_data (int bound_nn, double lim_nn, int bound_psi, int bound_beta, int solve_lapse, int solve_psi, int solve_shift, double precis=1.e-12, double relax_nn=0.5, double relax_psi=0.5, double relax_beta=0.5, int niter=100)
void init_data_loop (int bound_nn, double lim_nn, int bound_psi, int bound_beta, int solve_lapse, int solve_psi, int solve_shift, double precis=1.e-12, double precis_loop=1.e-12, double relax_nn=1., double relax_psi=1., double relax_beta=1., double relax_loop=1., int niter=100)
void init_data_spher (int bound_nn, double lim_nn, int bound_psi, int bound_beta, int solve_lapse, int solve_psi, int solve_shift, double precis=1.e-12, double relax=1., int niter=100)
void init_data_alt (int bound_nn, double lim_nn, int bound_psi, int bound_beta, int solve_lapse, int solve_psi, int solve_shift, double precis=1.e-12, double relax=1., int niter=100)
void init_data_CTS_gen (int bound_nn, double lim_nn, int bound_psi, int bound_beta, int solve_lapse, int solve_psi, int solve_shift, double precis=1.e-12, double relax_nn=1., double relax_psi=1., double relax_beta=1., int niter=100, double a=1., double zeta=4.)
const Scalar source_psi () const
 Source for $ \Psi $.
const Scalar source_nn () const
 Source for N.
const Vector source_beta () const
 Source for $ \beta $.
const Scalar source_b_tilde () const
 Source for b_tilde.
const Vector source_vector_b () const
 Source for vector_b.
const Valeur boundary_psi_Dir_evol () const
 Dirichlet boundary condition for $ \Psi $ (evolution).
const Valeur boundary_psi_Neu_evol () const
 Neumann boundary condition for $ \Psi $ (evolution).
const Valeur boundary_psi_Dir_spat () const
 Dirichlet boundary condition for $ \Psi $ (spatial).
const Valeur boundary_psi_Neu_spat () const
 Neumann boundary condition for $ \Psi $ (spatial).
const Valeur boundary_psi_app_hor () const
 Neumann boundary condition for $ \Psi $ (spatial).
const Valeur boundary_psi_Dir () const
 Dirichlet boundary condition for $ \Psi $ (spatial).
const Valeur boundary_nn_Dir_kk () const
 Dirichlet boundary condition for N using the extrinsic curvature.
const Valeur boundary_nn_Neu_kk (int nn=1) const
 Neumann boundary condition for N using the extrinsic curvature.
const Valeur boundary_nn_Neu_Cook () const
 Neumann boundary condition for N using Cook's boundary condition.
const Valeur boundary_nn_Dir_eff (double aa) const
 Dirichlet boundary condition for N (effectif) $ \partial_r N + a N = 0 $.
const Valeur boundary_nn_Dir_lapl (int mer=1) const
 Dirichlet boundary condition for N fixing the divergence of the connection form $ \omega $.
const Valeur boundary_nn_Neu_eff (double aa) const
 Neumann boundary condition on nn (effectif) $ \partial_r N + a N = 0 $.
const Valeur boundary_nn_Dir (double aa) const
 Dirichlet boundary condition $ N = a $.
const Valeur boundary_beta_r () const
 Component r of boundary value of $ \beta $.
const Valeur boundary_beta_theta () const
 Component theta of boundary value of $ \beta $.
const Valeur boundary_beta_phi (double om) const
 Component phi of boundary value of $ \beta $.
const Valeur boundary_beta_x (double om) const
 Component x of boundary value of $ \beta $.
const Valeur boundary_beta_y (double om) const
 Component y of boundary value of $ \beta $.
const Valeur boundary_beta_z () const
 Component z of boundary value of $ \beta $.
const Valeur beta_boost_x () const
 Boundary value for a boost in x-direction.
const Valeur beta_boost_z () const
 Boundary value for a boost in z-direction.
const Vector vv_bound_cart (double om) const
 Vector $ V^i $ for boundary conditions in cartesian.
const Vector vv_bound_cart_bin (double om, int hole=0) const
 Vector $ V^i $ for boundary conditions in cartesian for binary systems.
const Valeur boundary_vv_x (double om) const
 Component x of boundary value of $ V^i $.
const Valeur boundary_vv_y (double om) const
 Component y of boundary value of $ V^i $.
const Valeur boundary_vv_z (double om) const
 Component z of boundary value of $ V^i $.
const Valeur boundary_vv_x_bin (double om, int hole=0) const
 Component x of boundary value of $ V^i $.
const Valeur boundary_vv_y_bin (double om, int hole=0) const
 Component y of boundary value of $ V^i $.
const Valeur boundary_vv_z_bin (double om, int hole=0) const
 Component z of boundary value of $ V^i $.
const Valeur boundary_b_tilde_Neu () const
 Neumann boundary condition for b_tilde.
const Valeur boundary_b_tilde_Dir () const
 Dirichlet boundary condition for b_tilde.
void update_aa ()
 Conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature: $ A^{ij} = \Psi^4 \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.
double regularisation (const Vector &shift_auto, const Vector &shift_comp, double ang_vel)
 Corrects shift_auto in such a way that the total $A^{ij}$ is equal to zero in the horizon, which should ensure the regularity of $K^{ij}$.
double regularise_one ()
 Corrects the shift in the innermost shell, so that it remains $ {\mathcal{C}}^2$ and that $A^{ij}$ equals zero on the horizon.
void met_kerr_perturb ()
 Initialisation of the metric tilde from equation (15) of Dain (2002).
void aa_kerr_ww (double mm, double aa)
double axi_break () const
 Breaking of the axial symmetry on the horizon.
void adapt_hor (double c_min, double c_max)
virtual void sauve (FILE *fich, bool partial_save) const
 Total or partial saves in a binary file.
virtual void set_psi_del_npsi (const Scalar &psi_in)
 Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
virtual void set_psi_del_n (const Scalar &psi_in)
 Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
virtual void set_npsi_del_psi (const Scalar &npsi_in)
 Sets the factor $ N\Psi $ at the current time step (jtime ) and deletes the value of $\Psi$.
virtual void set_npsi_del_n (const Scalar &npsi_in)
 Sets the factor $ N\Psi $ at the current time step (jtime ) and deletes the value of N.
virtual void set_hh (const Sym_tensor &hh_in)
 Sets the deviation $ h^{ij} $ of the conformal metric $ \tilde\gamma^{ij} $ from the flat metric $ f^{ij} $: $\tilde\gamma^{ij} = f^{ij} + h^{ij} $.
virtual void set_hata (const Sym_tensor &hata_in)
 Sets the conformal representation $ \hat{A}{ij} $ of the traceless part of the extrinsic curvature: $ \hat{A}^{ij} = \Psi^{10} \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.
virtual void set_hata_TT (const Sym_tensor_tt &hata_tt)
 Sets the TT part of $ \hat{A}^{ij} $ (see member hata_evol ).
virtual void set_hata_from_XAB (Param *par_bc=0x0, Param *par_mat=0x0)
 Sets the conformal representation $ \hat{A}{ij} $ of the traceless part of the extrinsic curvature from its potentials A, $ \tilde{B} $ and $ X^i $.
virtual const Scalarnn () const
 Lapse function N at the current time step (jtime ).
virtual const Sym_tensorgam_dd () const
 Induced metric (covariant components $ \gamma_{ij} $) at the current time step (jtime ).
virtual const Sym_tensorgam_uu () const
 Induced metric (contravariant components $ \gamma^{ij} $) at the current time step (jtime ).
virtual const Sym_tensork_dd () const
 Extrinsic curvature tensor (covariant components $ K_{ij} $) at the current time step (jtime ).
virtual const Sym_tensork_uu () const
 Extrinsic curvature tensor (contravariant components $ K^{ij} $) at the current time step (jtime ).
virtual const ScalarA_hata () const
 Returns the potential A of $ \hat{A}^{ij} $.
virtual const ScalarB_hata () const
 Returns the potential $\tilde{B}$ of $ \hat{A}^{ij} $.
virtual const Scalarpsi () const
 Conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
const Scalarpsi4 () const
 Factor $ \Psi^4 $ at the current time step (jtime ).
const Scalarln_psi () const
 Logarithm of $ \Psi $ at the current time step (jtime ).
virtual const Scalarnpsi () const
 Factor $ N\Psi $ at the current time step (jtime ).
virtual const Sym_tensorhh (Param *=0x0, Param *=0x0) const
 Deviation $ h^{ij} $ of the conformal metric $ \tilde\gamma^{ij} $ from the flat metric $ f^{ij} $: $\tilde\gamma^{ij} = f^{ij} + h^{ij} $.
virtual const Sym_tensorhata () const
 Conformal representation $ \hat{A}^{ij} $ of the traceless part of the extrinsic curvature: $ \hat{A}^{ij} = \Psi^{10} \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.
virtual Sym_tensor aa () const
 Conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature: $ A^{ij} = \Psi^4 \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.
virtual const Scalartrk () const
 Trace K of the extrinsic curvature at the current time step (jtime ).
virtual const Vectorhdirac () const
 Vector $ H^i = {\cal D}_j \tilde\gamma^{ij} $ which vanishes in Dirac gauge.
virtual const Vectorvec_X (int method_poisson=6) const
 Vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $.
void compute_X_from_momentum_constraint (const Vector &hat_S, const Sym_tensor_tt &hata_tt, int iter_max=200, double precis=1.e-12, double relax=0.8, int methode_poisson=6)
 Computes the vector $ X^i $ from the conformally-rescaled momentum $ \hat{S}^i = \Psi^6 S^i $, using the momentum constraint.
void set_scheme_order (int ord)
 Sets the order of the finite-differences scheme.
int get_scheme_order () const
 Gets the order of the finite-differences scheme.
int get_latest_j () const
 Gets the latest value of time step index.
const Evolution_std< double > & get_time () const
 Gets the time coordinate t at successive time steps.
virtual const Vectorbeta () const
 shift vector $ \beta^i $ at the current time step (jtime )
const Metricgam () const
 Induced metric $ \mathbf{\gamma} $ at the current time step (jtime ).
Tbl check_hamiltonian_constraint (const Scalar *energy_density=0x0, ostream &ost=cout, bool verb=true) const
 Checks the level at which the hamiltonian constraint is verified.
Tbl check_momentum_constraint (const Vector *momentum_density=0x0, ostream &ost=cout, bool verb=true) const
 Checks the level at which the momentum constraints are verified.
Tbl check_dynamical_equations (const Sym_tensor *strain_tensor=0x0, const Scalar *energy_density=0x0, ostream &ost=cout, bool verb=true) const
 Checks the level at which the dynamical equations are verified.
virtual double adm_mass () const
 Returns the ADM mass (geometrical units) at the current step.
void save (const char *rootname) const
 Saves in a binary file.

Protected Member Functions

virtual ostream & operator>> (ostream &) const
 Operator >> (virtual function called by the operator<<).
virtual void del_deriv () const
 Deletes all the derived quantities.
void set_der_0x0 () const
 Sets to 0x0 all the pointers on derived quantities.

Protected Attributes

Map_afmp
 Affine mapping.
int nz
 Number of zones.
double radius
 Radius of the horizon in LORENE's units.
double omega
 Angular velocity in LORENE's units.
double boost_x
 Boost velocity in x-direction.
double boost_z
 Boost velocity in z-direction.
double regul
 Intensity of the correction on the shift vector.
Evolution_std< Scalarn_auto_evol
 Values at successive time steps of the lapse function $ N_{auto} $.
Evolution_std< Scalarn_comp_evol
 Values at successive time steps of the lapse function $ N_{comp} $.
Evolution_std< Scalarpsi_auto_evol
 Values at successive time steps of the conformal factor $ \Psi_{auto} $.
Evolution_std< Scalarpsi_comp_evol
 Values at successive time steps of the lapse function $ \Psi_{comp} $.
Evolution_std< Vectordn_evol
 Values at successive time steps of the covariant derivative of the lapse with respect to the flat metric $ \overline\nabla_i N $.
Evolution_std< Vectordpsi_evol
 Values at successive time steps of the covariant derivative of the conformal factor $ \overline\nabla_i \Psi $.
Evolution_std< Vectorbeta_auto_evol
 Values at successive time steps of the shift function $ \beta^i_{auto} $.
Evolution_std< Vectorbeta_comp_evol
 Values at successive time steps of the shift function $ \beta^i_{comp} $.
Evolution_std< Sym_tensoraa_auto_evol
 Values at successive time steps of the components $ A^{ij}_{auto} $ of the conformal representation of the traceless part of the extrinsic curvature:.
Evolution_std< Sym_tensoraa_comp_evol
 Values at successive time steps of the components $ A^{ij}_{comp} $ of the conformal representation of the traceless part of the extrinsic curvature:.
Evolution_std< Sym_tensoraa_nn
 Values at successive time steps of the components $ A^{ij}*2N $.
Evolution_std< Scalaraa_quad_evol
 Values at successive time steps of the components $ A^{ij}A_{ij} $.
Metric met_gamt
 3 metric tilde
Sym_tensor gamt_point
 Time derivative of the 3-metric tilde.
Scalar trK
 Trace of the extrinsic curvature.
Scalar trK_point
 Time derivative of the trace of the extrinsic curvature.
Scalar decouple
 Function used to construct $ A^{ij}_{auto} $ from the total $A^{ij}$.
const Metric_flatff
 Pointer on the flat metric $ f_{ij} $ with respect to which the conformal decomposition is performed.
Evolution_std< Scalarpsi_evol
 Values at successive time steps of the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
Evolution_std< Scalarnpsi_evol
 Values at successive time steps of the factor $ N\Psi $.
Evolution_std< Sym_tensorhh_evol
 Values at successive time steps of the components $ h^{ij} $.
Evolution_std< Sym_tensorhata_evol
 Values at successive time steps of the components $ \hat{A}^{ij} $.
Evolution_std< ScalarA_hata_evol
 Potential A associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $.
Evolution_std< ScalarB_hata_evol
 Potential $ \tilde{B} $ associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $.
Metricp_tgamma
 Pointer on the conformal metric $ \tilde\gamma_{ij} $ at the current time step (jtime).
Scalarp_psi4
 Pointer on the factor $ \Psi^4 $ at the current time step (jtime).
Scalarp_ln_psi
 Pointer on the logarithm of $ \Psi $ at the current time step (jtime).
Vectorp_hdirac
 Pointer on the vector $ H^i = {\cal D}_j \tilde\gamma^{ij} $ (which vanishes in Dirac gauge), at the current time step (jtime).
Vectorp_vec_X
 Pointer on the vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $.
int depth
 Number of stored time slices.
int scheme_order
 Order of the finite-differences scheme for the computation of time derivatives.
int jtime
 Time step index of the latest slice.
Evolution_std< double > the_time
 Time label of each slice.
Evolution_std< Sym_tensorgam_dd_evol
 Values at successive time steps of the covariant components of the induced metric $ \gamma_{ij} $.
Evolution_std< Sym_tensorgam_uu_evol
 Values at successive time steps of the contravariant components of the induced metric $ \gamma^{ij} $.
Evolution_std< Sym_tensork_dd_evol
 Values at successive time steps of the covariant components of the extrinsic curvature tensor $ K_{ij} $.
Evolution_std< Sym_tensork_uu_evol
 Values at successive time steps of the contravariant components of the extrinsic curvature tensor $ K^{ij} $.
Evolution_std< Scalarn_evol
 Values at successive time steps of the lapse function N.
Evolution_std< Vectorbeta_evol
 Values at successive time steps of the shift vector $ \beta^i $.
Evolution_std< Scalartrk_evol
 Values at successive time steps of the trace K of the extrinsic curvature.
Evolution_full< Tbladm_mass_evol
 ADM mass at each time step, since the creation of the slice.
Metricp_gamma
 Pointer on the induced metric at the current time step (jtime).

Friends

class Bin_hor
ostream & operator<< (ostream &, const Time_slice &)
 Display.

Detailed Description

Spacelike time-slice of an Isolated Horizon in a 3+1 spacetime with conformal decomposition.

No gauge choice imposed. ()

Definition at line 252 of file isol_hor.h.


Constructor & Destructor Documentation

Isol_hor::Isol_hor ( Map_af mpi,
int  depth_in = 3 
)

Standard constructor.

Parameters:
mpi affine mapping
depth_in number of stored time slices; this parameter is used to set the scheme_order member with scheme_order = depth_in - 1. scheme_order can be changed afterwards by the method set_scheme_order(int).

Definition at line 172 of file isol_hor.C.

Isol_hor::Isol_hor ( Map_af mpi,
const Scalar lapse_in,
const Scalar psi_in,
const Vector shift_in,
const Sym_tensor aa_in,
const Metric gamt,
const Sym_tensor gamt_point,
const Scalar trK,
const Scalar trK_point,
const Metric_flat ff_in,
int  depth_in = 3 
)

Constructor from conformal decomposition.

Parameters:
mpi affine mapping
lapse_in lapse function N
psi_in conformal factor $\Psi$ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $
shift_in shift vector
aa_in conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature: $ A^{ij} = \Psi^4 \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $
gamt 3-metric tilde
gamt_point time derivative of the 3-metric tilde
trK trace K of the extrinsic curvature
trK_point time derivative of the trace K of the extrinsic curvature
ff_in reference flat metric with respect to which the conformal decomposition is performed
depth_in number of stored time slices; this parameter is used to set the scheme_order member with scheme_order = depth_in - 1. scheme_order can be changed afterwards by the method set_scheme_order(int).

Definition at line 189 of file isol_hor.C.

References Metric_flat::con(), Metric::con(), Time_slice_conf::ff, Time_slice_conf::hh_evol, Time_slice::jtime, met_gamt, Time_slice::the_time, trK, Time_slice::trk_evol, and Evolution_std< TyT >::update().

Isol_hor::Isol_hor ( const Isol_hor isolhor_in  ) 

Copy constructor.

Definition at line 218 of file isol_hor.C.

Isol_hor::Isol_hor ( Map_af mp,
FILE *  fich,
bool  partial_read,
int  depth_in = 3 
)

Constructor from a binary file.

Parameters:
mpi affine mapping
fich file containing the saved isol_hor
partial_read indicates whether the full object must be read in file or whether the final construction is devoted to a constructor of a derived class
depth_in number of stored time slices; this parameter is used to set the scheme_order member with scheme_order = depth_in - 1. scheme_order can be changed afterwards by the method set_scheme_order(int).

Definition at line 249 of file isol_hor.C.

References beta_auto_evol, boost_x, boost_z, Metric_flat::con(), Time_slice::depth, Time_slice_conf::ff, fread_be(), gamt_point, Map::get_bvect_spher(), Map::get_mg(), Time_slice_conf::hh_evol, Time_slice::jtime, met_gamt, mp, n_auto_evol, omega, psi_auto_evol, Time_slice_conf::psi_evol, Time_slice::the_time, trK, Time_slice::trk_evol, trK_point, and Evolution_std< TyT >::update().

Isol_hor::~Isol_hor (  )  [virtual]

Destructor.

Definition at line 332 of file isol_hor.C.


Member Function Documentation

const Scalar & Time_slice_conf::A_hata (  )  const [virtual, inherited]

Returns the potential A of $ \hat{A}^{ij} $.

See the documentation of Sym_tensor for details. Returns the value at the current time step (jtime ).

Definition at line 660 of file time_slice_conf.C.

References Time_slice_conf::A_hata_evol, Time_slice_conf::hata_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::the_time, and Evolution_std< TyT >::update().

Sym_tensor Time_slice_conf::aa (  )  const [virtual, inherited]

Conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature: $ A^{ij} = \Psi^4 \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.

Returns the value at the current time step (jtime ).

Definition at line 761 of file time_slice_conf.C.

References Time_slice_conf::hata(), Time_slice_conf::psi(), and Time_slice_conf::psi4().

const Sym_tensor & Isol_hor::aa_auto (  )  const [virtual]

Conformal representation $ A^{ij}_{auto} $ of the traceless part of the extrinsic curvature: Returns the value at the current time step jtime.

Definition at line 499 of file isol_hor.C.

References aa_auto_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

const Sym_tensor & Isol_hor::aa_comp (  )  const [virtual]

Conformal representation $ A^{ij}_{comp} $ of the traceless part of the extrinsic curvature: Returns the value at the current time step jtime.

Definition at line 505 of file isol_hor.C.

References aa_comp_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

const Scalar & Isol_hor::aa_quad (  )  const [virtual]

Conformal representation $ A^{ij}A_{ij} $.

Returns the value at the current time step jtime.

Definition at line 880 of file isol_hor.C.

References Time_slice_conf::aa(), aa_quad_evol, contract(), Evolution< TyT >::is_known(), Time_slice::jtime, met_gamt, Time_slice_conf::psi4(), Time_slice::the_time, and Evolution_std< TyT >::update().

double Time_slice::adm_mass (  )  const [virtual, inherited]
double Isol_hor::ang_mom_adm (  )  const
double Isol_hor::ang_mom_hor (  )  const
double Isol_hor::area_hor (  )  const

Area of the horizon.

Definition at line 153 of file phys_param.C.

References darea_hor(), Map_af::integrale_surface(), mp, Scalar::raccord(), and radius.

double Isol_hor::axi_break (  )  const
const Scalar & Time_slice_conf::B_hata (  )  const [virtual, inherited]

Returns the potential $\tilde{B}$ of $ \hat{A}^{ij} $.

See the documentation of Sym_tensor_tt for details. Returns the value at the current time step (jtime ).

Definition at line 674 of file time_slice_conf.C.

References Time_slice_conf::A_hata_evol, Time_slice_conf::B_hata_evol, Time_slice_conf::hata_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::the_time, and Evolution_std< TyT >::update().

const Scalar Isol_hor::b_tilde (  )  const

Radial component of the shift with respect to the conformal metric.

Definition at line 132 of file phys_param.C.

References Time_slice::beta(), contract(), Tensor::down(), met_gamt, and Metric::radial_vect().

const Vector & Time_slice::beta (  )  const [virtual, inherited]

shift vector $ \beta^i $ at the current time step (jtime )

Definition at line 83 of file time_slice_access.C.

References Time_slice::beta_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

const Vector & Isol_hor::beta_auto (  )  const [virtual]

Shift function $ \beta^i_{auto} $ at the current time step jtime.

Definition at line 487 of file isol_hor.C.

References beta_auto_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

const Valeur Isol_hor::beta_boost_x (  )  const

Boundary value for a boost in x-direction.

Definition at line 1140 of file bound_hor.C.

References boost_x, Mg3d::get_angu(), Map::get_mg(), mp, and Mg3d::std_base_vect_cart().

const Valeur Isol_hor::beta_boost_z (  )  const

Boundary value for a boost in z-direction.

Definition at line 1151 of file bound_hor.C.

References boost_z, Mg3d::get_angu(), Map::get_mg(), mp, and Mg3d::std_base_vect_cart().

void Isol_hor::beta_comp ( const Isol_hor comp  ) 

Imports the part of $\beta^i$ due to the companion hole comp.

The total $\beta^i$ is then calculated.

Definition at line 708 of file isol_hor.C.

References beta_auto(), beta_comp(), beta_comp_evol, Time_slice::beta_evol, Map::get_bvect_cart(), Map::get_bvect_spher(), Time_slice::jtime, mp, Time_slice::the_time, and Evolution_std< TyT >::update().

const Vector & Isol_hor::beta_comp (  )  const [virtual]

Shift function $ \beta^i_{comp} $ at the current time step jtime.

Definition at line 493 of file isol_hor.C.

References beta_comp_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

const Valeur Isol_hor::boundary_b_tilde_Dir (  )  const
const Valeur Isol_hor::boundary_b_tilde_Neu (  )  const
const Valeur Isol_hor::boundary_beta_phi ( double  om  )  const
const Valeur Isol_hor::boundary_beta_r (  )  const
const Valeur Isol_hor::boundary_beta_theta (  )  const
const Valeur Isol_hor::boundary_beta_x ( double  om  )  const
const Valeur Isol_hor::boundary_beta_y ( double  om  )  const
const Valeur Isol_hor::boundary_beta_z (  )  const
const Valeur Isol_hor::boundary_nn_Dir ( double  aa  )  const

Dirichlet boundary condition $ N = a $.

Definition at line 643 of file bound_hor.C.

References Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), mp, Time_slice_conf::nn(), and Scalar::val_grid_point().

const Valeur Isol_hor::boundary_nn_Dir_eff ( double  aa  )  const
const Valeur Isol_hor::boundary_nn_Dir_kk (  )  const
const Valeur Isol_hor::boundary_nn_Dir_lapl ( int  mer = 1  )  const
const Valeur Isol_hor::boundary_nn_Neu_Cook (  )  const
const Valeur Isol_hor::boundary_nn_Neu_eff ( double  aa  )  const

Neumann boundary condition on nn (effectif) $ \partial_r N + a N = 0 $.

Definition at line 618 of file bound_hor.C.

References Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), mp, Time_slice_conf::nn(), and Scalar::val_grid_point().

const Valeur Isol_hor::boundary_nn_Neu_kk ( int  nn = 1  )  const
const Valeur Isol_hor::boundary_psi_app_hor (  )  const

Neumann boundary condition for $ \Psi $ (spatial).

Definition at line 293 of file bound_hor.C.

References Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), mp, Time_slice_conf::psi(), radius, and Scalar::val_grid_point().

const Valeur Isol_hor::boundary_psi_Dir (  )  const

Dirichlet boundary condition for $ \Psi $ (spatial).

Definition at line 320 of file bound_hor.C.

References Mg3d::get_angu(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), mp, Time_slice_conf::psi(), Scalar::std_spectral_base(), and Scalar::val_grid_point().

const Valeur Isol_hor::boundary_psi_Dir_evol (  )  const
const Valeur Isol_hor::boundary_psi_Dir_spat (  )  const
const Valeur Isol_hor::boundary_psi_Neu_evol (  )  const
const Valeur Isol_hor::boundary_psi_Neu_spat (  )  const
const Valeur Isol_hor::boundary_vv_x ( double  om  )  const
const Valeur Isol_hor::boundary_vv_x_bin ( double  om,
int  hole = 0 
) const
const Valeur Isol_hor::boundary_vv_y ( double  om  )  const
const Valeur Isol_hor::boundary_vv_y_bin ( double  om,
int  hole = 0 
) const
const Valeur Isol_hor::boundary_vv_z ( double  om  )  const
const Valeur Isol_hor::boundary_vv_z_bin ( double  om,
int  hole = 0 
) const
Tbl Time_slice::check_dynamical_equations ( const Sym_tensor strain_tensor = 0x0,
const Scalar energy_density = 0x0,
ostream &  ost = cout,
bool  verb = true 
) const [inherited]

Checks the level at which the dynamical equations are verified.

\[ \frac{\partial K_{ij}}{\partial t} - \pounds_{\mbox{\boldmath{$\beta $}}} K_{ij} = - D_i D_j N + N \left[ R_{ij} - 2 K_{ik} K^k_{\ j} + K K_{ij} + 4\pi \left( (S-E)\gamma_{ij} - 2 S_{ij} \right) \right] \]

Parameters:
strain_tensor : a pointer on the strain_tensor $ S_{ij} $ measured by the Eulerian observer of 4-velocity $\mbox{\boldmath{$n $}}$ ; if this is the null pointer, it is assumed that $ S_{ij} $ = 0 (vacuum).
energy_density : a pointer on the energy density E (see check_hamiltonian_constraint)
ost : output stream for a formatted output of the result
Returns:
Tbl 3D of size the number of domains times 3 times 3 (corresponding to the rank-2 tensor, with the symmetry in the components) containing the absolute ( if $ J_i $ = 0 ) or the relative (in presence of matter) error in max version.

Definition at line 135 of file tslice_check_einstein.C.

References Tensor::annule_domain(), Time_slice::beta(), Metric::con(), contract(), Tensor::derive_con(), Scalar::derive_con(), Sym_tensor::derive_lie(), Time_slice::gam(), Tbl::get_dim(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Time_slice::jtime, Time_slice::k_dd(), Time_slice::k_dd_evol, Time_slice::k_uu(), maxabs(), Time_slice::nn(), Tensor_sym::position(), Metric::ricci(), Time_slice::scheme_order, Itbl::set(), Evolution< TyT >::time_derive(), Tensor::trace(), Time_slice::trk(), and Tensor::up_down().

Tbl Time_slice::check_hamiltonian_constraint ( const Scalar energy_density = 0x0,
ostream &  ost = cout,
bool  verb = true 
) const [inherited]

Checks the level at which the hamiltonian constraint is verified.

\[ R + K^2 - K_{ij}K^{ij} = 16\pi E \]

Parameters:
energy_density : a pointer on the energy density E measured by the Eulerian observer of 4-velocity $\mbox{\boldmath{$n $}}$ ; if this is the null pointer, it is assumed that E = 0 (vacuum).
ost : output stream for a formatted output of the result
Returns:
Tbl of size the number of domains containing the absolute ( if E = 0 ) or the relative (in presence of matter) error in max version.

Definition at line 75 of file tslice_check_einstein.C.

References contract(), Scalar::dec_dzpuis(), Time_slice::gam(), Time_slice::k_dd(), Time_slice::k_uu(), maxabs(), Metric::ricci_scal(), and Time_slice::trk().

Tbl Time_slice::check_momentum_constraint ( const Vector momentum_density = 0x0,
ostream &  ost = cout,
bool  verb = true 
) const [inherited]

Checks the level at which the momentum constraints are verified.

\[ D_j K_i^{\ j} - D_i K = 8 \pi J_i \]

Parameters:
momentum_density : a pointer on the momentum density $ J_i $ measured by the Eulerian observer of 4-velocity $\mbox{\boldmath{$n $}}$ ; if this is the null pointer, it is assumed that$ J_i $ = 0 (vacuum).
ost : output stream for a formatted output of the result
Returns:
Tbl 2D of size the number of domains times 3 (components)containing the absolute ( if $ J_i $ = 0 ) or the relative (in presence of matter) error in max version.

Definition at line 105 of file tslice_check_einstein.C.

References Scalar::derive_con(), Sym_tensor::divergence(), Time_slice::gam(), Tensor::get_index_type(), Time_slice::k_uu(), maxabs(), and Time_slice::trk().

void Time_slice_conf::compute_X_from_momentum_constraint ( const Vector hat_S,
const Sym_tensor_tt hata_tt,
int  iter_max = 200,
double  precis = 1.e-12,
double  relax = 0.8,
int  methode_poisson = 6 
) [inherited]

Computes the vector $ X^i $ from the conformally-rescaled momentum $ \hat{S}^i = \Psi^6 S^i $, using the momentum constraint.

Definition at line 833 of file time_slice_conf.C.

References abs(), contract(), Tensor::get_index_type(), max(), and Vector::poisson().

const Scalar Isol_hor::darea_hor (  )  const

Element of area of the horizon.

Definition at line 142 of file phys_param.C.

References Time_slice_conf::gam_dd(), sqrt(), and Scalar::std_spectral_base().

void Time_slice_conf::del_deriv (  )  const [protected, virtual, inherited]
const Vector & Isol_hor::dnn (  )  const [virtual]

Covariant derivative of the lapse function $ \overline\nabla_i N $ at the current time step jtime.

Definition at line 475 of file isol_hor.C.

References dn_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

const Vector & Isol_hor::dpsi (  )  const [virtual]

Covariant derivative with respect to the flat metric of the conformal factor $ \overline\nabla_i \Psi $ at the current time step jtime.

Definition at line 481 of file isol_hor.C.

References dpsi_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

Scalar Isol_hor::expansion (  )  const

Expansion of the outgoing null normal ($ \bf n + \bf s $).

Definition at line 259 of file phys_param.C.

References contract(), Time_slice::gam(), Time_slice_conf::k_dd(), and Time_slice_conf::trk().

const Metric & Time_slice::gam (  )  const [inherited]

Induced metric $ \mathbf{\gamma} $ at the current time step (jtime ).

Definition at line 91 of file time_slice_access.C.

References Time_slice::gam_dd(), and Time_slice::p_gamma.

const Sym_tensor & Time_slice_conf::gam_dd (  )  const [virtual, inherited]

Induced metric (covariant components $ \gamma_{ij} $) at the current time step (jtime ).

Reimplemented from Time_slice.

Definition at line 604 of file time_slice_conf.C.

References Time_slice::gam_dd_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice_conf::psi4(), Time_slice_conf::tgam(), Time_slice::the_time, and Evolution_std< TyT >::update().

const Sym_tensor & Time_slice_conf::gam_uu (  )  const [virtual, inherited]

Induced metric (contravariant components $ \gamma^{ij} $) at the current time step (jtime ).

Reimplemented from Time_slice.

Definition at line 615 of file time_slice_conf.C.

References Time_slice::gam_uu_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice_conf::psi4(), Time_slice_conf::tgam(), Time_slice::the_time, and Evolution_std< TyT >::update().

double Isol_hor::get_boost_x (  )  const [inline]

Returns the boost velocity in x-direction.

Definition at line 443 of file isol_hor.h.

References boost_x.

double Isol_hor::get_boost_z (  )  const [inline]

Returns the boost velocity in z-direction.

Definition at line 452 of file isol_hor.h.

References boost_z.

const Scalar Isol_hor::get_decouple (  )  const [inline]

Returns the function used to construct tkij_auto from tkij_tot .

Definition at line 517 of file isol_hor.h.

References decouple.

int Time_slice::get_latest_j (  )  const [inline, inherited]

Gets the latest value of time step index.

Definition at line 339 of file time_slice.h.

References Time_slice::jtime.

const Map_af& Isol_hor::get_mp (  )  const [inline]

Returns the mapping (readonly).

Definition at line 416 of file isol_hor.h.

References mp.

double Isol_hor::get_omega (  )  const [inline]

Returns the angular velocity.

Definition at line 434 of file isol_hor.h.

References omega.

double Isol_hor::get_radius (  )  const [inline]

Returns the radius of the horizon.

Definition at line 424 of file isol_hor.h.

References radius.

int Time_slice::get_scheme_order (  )  const [inline, inherited]

Gets the order of the finite-differences scheme.

Definition at line 336 of file time_slice.h.

References Time_slice::scheme_order.

const Evolution_std<double>& Time_slice::get_time (  )  const [inline, inherited]

Gets the time coordinate t at successive time steps.

Definition at line 342 of file time_slice.h.

References Time_slice::the_time.

const Sym_tensor & Time_slice_conf::hata (  )  const [virtual, inherited]
const Vector & Time_slice_conf::hdirac (  )  const [virtual, inherited]

Vector $ H^i = {\cal D}_j \tilde\gamma^{ij} $ which vanishes in Dirac gauge.

Reimplemented in Tslice_dirac_max.

Definition at line 811 of file time_slice_conf.C.

References Time_slice_conf::ff, Time_slice_conf::hh(), and Time_slice_conf::p_hdirac.

const Sym_tensor & Time_slice_conf::hh ( Param = 0x0,
Param = 0x0 
) const [virtual, inherited]

Deviation $ h^{ij} $ of the conformal metric $ \tilde\gamma^{ij} $ from the flat metric $ f^{ij} $: $\tilde\gamma^{ij} = f^{ij} + h^{ij} $.

Returns the value at the current time step (jtime ).

Reimplemented in Tslice_dirac_max.

Definition at line 754 of file time_slice_conf.C.

References Time_slice_conf::hh_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

void Isol_hor::init_bhole (  ) 
void Isol_hor::init_bhole_seul (  ) 
void Isol_hor::init_met_trK (  ) 

Sets the 3-metric tilde to the flat metric and gamt_point, trK and trK_point to zero.

Definition at line 781 of file isol_hor.C.

References Map::flat_met_spher(), gamt_point, met_gamt, mp, Scalar::set_etat_zero(), Tensor::set_etat_zero(), trK, and trK_point.

const Sym_tensor & Time_slice_conf::k_dd (  )  const [virtual, inherited]

Extrinsic curvature tensor (covariant components $ K_{ij} $) at the current time step (jtime ).

Reimplemented from Time_slice.

Definition at line 626 of file time_slice_conf.C.

References Time_slice::gam(), Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::k_dd_evol, Time_slice_conf::k_uu(), Time_slice::the_time, and Evolution_std< TyT >::update().

const Sym_tensor & Time_slice_conf::k_uu (  )  const [virtual, inherited]

Extrinsic curvature tensor (contravariant components $ K^{ij} $) at the current time step (jtime ).

Reimplemented from Time_slice.

Definition at line 639 of file time_slice_conf.C.

References Time_slice::gam(), Time_slice_conf::hata(), Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::k_uu_evol, Time_slice_conf::psi(), Time_slice_conf::psi4(), Time_slice::the_time, Time_slice_conf::trk(), and Evolution_std< TyT >::update().

double Isol_hor::kappa_hor (  )  const

Surface gravity.

Definition at line 214 of file phys_param.C.

References ang_mom_hor(), pow(), radius_hor(), and sqrt().

const Scalar & Time_slice_conf::ln_psi (  )  const [inherited]

Logarithm of $ \Psi $ at the current time step (jtime ).

Definition at line 715 of file time_slice_conf.C.

References log(), Time_slice_conf::p_ln_psi, Time_slice_conf::psi(), and Scalar::std_spectral_base().

double Isol_hor::mass_hor (  )  const

Mass computed at the horizon.

Definition at line 202 of file phys_param.C.

References ang_mom_hor(), pow(), radius_hor(), and sqrt().

void Isol_hor::met_kerr_perturb (  ) 

Initialisation of the metric tilde from equation (15) of Dain (2002).

The determinant of this conformal metric is not one.

Definition at line 891 of file isol_hor.C.

References Metric_flat::con(), Metric::con(), Metric::cov(), Metric::determinant(), Time_slice_conf::ff, Time_slice::gam(), Time_slice_conf::hh_evol, Time_slice::jtime, met_gamt, norme(), pow(), set_psi(), Scalar::std_spectral_base(), Time_slice::the_time, and Evolution_std< TyT >::update().

const Scalar & Isol_hor::n_auto (  )  const [virtual]

Lapse function $ N_{auto} $ at the current time step jtime.

Definition at line 451 of file isol_hor.C.

References Evolution< TyT >::is_known(), Time_slice::jtime, and n_auto_evol.

void Isol_hor::n_comp ( const Isol_hor comp  ) 

Imports the part of N due to the companion hole comp .

The total N is then calculated.

It also imports the covariant derivative of N and construct the total $\nabla_i N$.

Definition at line 573 of file isol_hor.C.

References Scalar::derive_cov(), dn_evol, Time_slice_conf::ff, Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::import(), Time_slice::jtime, mp, n_auto(), n_comp_evol, Time_slice::n_evol, Scalar::std_spectral_base(), Time_slice::the_time, and Evolution_std< TyT >::update().

const Scalar & Isol_hor::n_comp (  )  const [virtual]

Lapse function $ N_{comp} $ at the current time step jtime.

Definition at line 457 of file isol_hor.C.

References Evolution< TyT >::is_known(), Time_slice::jtime, and n_comp_evol.

const Scalar & Time_slice_conf::nn (  )  const [virtual, inherited]

Lapse function N at the current time step (jtime ).

Reimplemented from Time_slice.

Definition at line 587 of file time_slice_conf.C.

References Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::n_evol, Time_slice_conf::npsi_evol, Time_slice_conf::psi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

const Scalar & Time_slice_conf::npsi (  )  const [virtual, inherited]
double Isol_hor::omega_hor (  )  const

Orbital velocity.

Definition at line 229 of file phys_param.C.

References ang_mom_hor(), pow(), radius_hor(), and sqrt().

void Isol_hor::operator= ( const Isol_hor isolhor_in  ) 
ostream & Isol_hor::operator>> ( ostream &  flux  )  const [protected, virtual]

Operator >> (virtual function called by the operator<<).

Reimplemented from Time_slice.

Definition at line 374 of file isol_hor.C.

References Time_slice::adm_mass(), ang_mom_adm(), ang_mom_hor(), area_hor(), boost_x, boost_z, mass_hor(), omega_hor(), and radius.

const Scalar & Time_slice_conf::psi (  )  const [virtual, inherited]

Conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.

$ \Psi $ is defined by

\[ \Psi := \left( \frac{\det\gamma_{ij}}{\det f_{ij}} \right) ^{1/12} \]

Returns the value at the current time step (jtime ).

Definition at line 689 of file time_slice_conf.C.

References Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::n_evol, Time_slice_conf::npsi_evol, Time_slice_conf::psi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

const Scalar & Time_slice_conf::psi4 (  )  const [inherited]

Factor $ \Psi^4 $ at the current time step (jtime ).

Definition at line 703 of file time_slice_conf.C.

References Time_slice_conf::p_psi4, pow(), Time_slice_conf::psi(), and Scalar::std_spectral_base().

const Scalar & Isol_hor::psi_auto (  )  const [virtual]

Conformal factor $ \Psi_{auto} $ at the current time step jtime.

Definition at line 463 of file isol_hor.C.

References Evolution< TyT >::is_known(), Time_slice::jtime, and psi_auto_evol.

void Isol_hor::psi_comp ( const Isol_hor comp  ) 

Imports the part of $\Psi$ due to the companion hole comp .

The total $\Psi$ is then calculated.

It also imports the covariant derivative of $\Psi$ and construct the total $\nabla_i \Psi$.

Definition at line 665 of file isol_hor.C.

References Scalar::derive_cov(), dpsi_evol, Time_slice_conf::ff, Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::import(), Time_slice::jtime, mp, psi_auto(), psi_comp_evol, Time_slice_conf::psi_evol, Scalar::std_spectral_base(), Time_slice::the_time, and Evolution_std< TyT >::update().

const Scalar & Isol_hor::psi_comp (  )  const [virtual]

Conformal factor $ \Psi_{comp} $ at the current time step jtime.

Definition at line 469 of file isol_hor.C.

References Evolution< TyT >::is_known(), Time_slice::jtime, and psi_comp_evol.

const Vector Isol_hor::radial_vect_hor (  )  const
double Isol_hor::radius_hor (  )  const

Radius of the horizon.

Definition at line 163 of file phys_param.C.

References area_hor(), and pow().

double Isol_hor::regularisation ( const Vector shift_auto,
const Vector shift_comp,
double  ang_vel 
)

Corrects shift_auto in such a way that the total $A^{ij}$ is equal to zero in the horizon, which should ensure the regularity of $K^{ij}$.

WARNING : this should only be used for a black hole in a binary system Bin_hor.

Parameters:
comp [input]: the part of $\beta^i$ generated by the companion hole.

Definition at line 84 of file regularisation.C.

References Cmp::annule(), Tensor::annule_domain(), beta_auto_evol, Valeur::c, Vector::change_triad(), Valeur::coef_i(), diffrelmax(), Map::get_bvect_cart(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Map::get_rot_phi(), Scalar::get_spectral_va(), Tensor::get_triad(), Time_slice::jtime, norme(), nz, pow(), Vector::set(), Valeur::set_etat_c_qcq(), Scalar::set_spectral_va(), Vector::std_spectral_base(), Time_slice::the_time, and Evolution_std< TyT >::update().

double Isol_hor::regularise_one (  ) 

Corrects the shift in the innermost shell, so that it remains $ {\mathcal{C}}^2$ and that $A^{ij}$ equals zero on the horizon.

return the relative difference between the shift before and after the regularisation.

WARNING this should only be used for an isolated black hole.

Definition at line 179 of file regularisation.C.

References Cmp::annule(), Tensor::annule_domain(), Time_slice::beta(), Time_slice::beta_evol, boost_x, boost_z, Vector::change_triad(), Valeur::coef_i(), diffrelmax(), Map::get_bvect_cart(), Map::get_bvect_spher(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Time_slice::jtime, mp, norme(), nz, omega, pow(), Map::r, Vector::set(), Valeur::set_etat_c_qcq(), Scalar::set_spectral_va(), Vector::std_spectral_base(), Time_slice::the_time, Evolution_std< TyT >::update(), Map_af::val_r(), Map::x, and Map::y.

void Isol_hor::sauve ( FILE *  fich,
bool  partial_save 
) const [virtual]

Total or partial saves in a binary file.

Parameters:
fich binary file
partial_save indicates whether the whole object must be saved.

Reimplemented from Time_slice.

Definition at line 397 of file isol_hor.C.

References beta_auto_evol, boost_x, boost_z, Metric::con(), Time_slice::depth, fwrite_be(), gamt_point, Evolution< TyT >::is_known(), Time_slice::jtime, met_gamt, n_auto_evol, omega, psi_auto_evol, Time_slice_conf::psi_evol, Scalar::sauve(), Tensor_sym::sauve(), trK, and trK_point.

void Time_slice::save ( const char *  rootname  )  const [inherited]

Saves in a binary file.

The saved data is sufficient to restore the whole time slice via the constructor from file.

Parameters:
rootname root for the file name; the current time step index will be appended to it.

Definition at line 457 of file time_slice.C.

References Time_slice::beta(), Time_slice::depth, fwrite_be(), Map::get_mg(), Tensor::get_mp(), Tensor::get_triad(), Time_slice::jtime, Time_slice::nn(), Time_slice::sauve(), Base_vect::sauve(), Map::sauve(), and Mg3d::sauve().

void Isol_hor::set_boost_x ( double  bo  )  [inline]

Sets the boost velocity in x-direction to bo .

Definition at line 447 of file isol_hor.h.

References boost_x.

void Isol_hor::set_boost_z ( double  bo  )  [inline]

Sets the boost velocity in z-direction to bo .

Definition at line 456 of file isol_hor.h.

References boost_z.

void Time_slice_conf::set_der_0x0 (  )  const [protected, inherited]

Sets to 0x0 all the pointers on derived quantities.

Reimplemented from Time_slice.

Definition at line 357 of file time_slice_conf.C.

References Time_slice_conf::p_hdirac, Time_slice_conf::p_ln_psi, Time_slice_conf::p_psi4, Time_slice_conf::p_tgamma, and Time_slice_conf::p_vec_X.

void Isol_hor::set_gamt ( const Metric gam_tilde  ) 
void Time_slice_conf::set_hata ( const Sym_tensor hata_in  )  [virtual, inherited]

Sets the conformal representation $ \hat{A}{ij} $ of the traceless part of the extrinsic curvature: $ \hat{A}^{ij} = \Psi^{10} \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.

Sets the value at the current time step (jtime ), and updates the potentials A_hata_evol, B_hata_evol and p_vec_X accordingly.

Definition at line 530 of file time_slice_conf.C.

References Time_slice_conf::A_hata_evol, Time_slice_conf::B_hata_evol, Evolution< TyT >::downdate(), Time_slice_conf::hata_evol, Time_slice::jtime, Time_slice::k_dd_evol, Time_slice::k_uu_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

void Time_slice_conf::set_hata_from_XAB ( Param par_bc = 0x0,
Param par_mat = 0x0 
) [virtual, inherited]

Sets the conformal representation $ \hat{A}{ij} $ of the traceless part of the extrinsic curvature from its potentials A, $ \tilde{B} $ and $ X^i $.

These potentials must be up-to-date. It sets the value at the current time step (jtime ).

Definition at line 562 of file time_slice_conf.C.

References Time_slice_conf::A_hata_evol, Time_slice_conf::B_hata_evol, Evolution< TyT >::downdate(), Time_slice_conf::ff, Tensor::get_mp(), Time_slice_conf::hata_evol, Map::inc_dzpuis(), Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::k_dd_evol, Time_slice::k_uu_evol, Vector::ope_killing_conf(), Time_slice_conf::p_vec_X, Time_slice::the_time, and Evolution_std< TyT >::update().

void Time_slice_conf::set_hata_TT ( const Sym_tensor_tt hata_tt  )  [virtual, inherited]
void Time_slice_conf::set_hh ( const Sym_tensor hh_in  )  [virtual, inherited]
Map_af& Isol_hor::set_mp (  )  [inline]

Read/write of the mapping.

Definition at line 419 of file isol_hor.h.

References mp.

void Isol_hor::set_nn ( const Scalar nn_in  ) 
void Time_slice_conf::set_npsi_del_n ( const Scalar npsi_in  )  [virtual, inherited]

Sets the factor $ N\Psi $ at the current time step (jtime ) and deletes the value of N.

Definition at line 474 of file time_slice_conf.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::jtime, Time_slice::n_evol, Time_slice_conf::npsi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

void Time_slice_conf::set_npsi_del_psi ( const Scalar npsi_in  )  [virtual, inherited]
void Isol_hor::set_omega ( double  ome  )  [inline]

Sets the angular velocity to ome .

Definition at line 438 of file isol_hor.h.

References omega.

void Isol_hor::set_psi ( const Scalar psi_in  ) 

Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.

$ \Psi $ is defined by

\[ \Psi := \left( \frac{\det\gamma_{ij}}{\det f_{ij}} \right) ^{1/12} \]

Sets the value at the current time step (jtime ) and delete all quantities which depend on $ \Psi $.

Definition at line 511 of file isol_hor.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::gam_dd_evol, Time_slice::gam_uu_evol, Time_slice::jtime, Time_slice::k_dd_evol, Time_slice::k_uu_evol, Time_slice::p_gamma, Time_slice_conf::p_ln_psi, Time_slice_conf::p_psi4, Time_slice_conf::psi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

void Time_slice_conf::set_psi_del_n ( const Scalar psi_in  )  [virtual, inherited]

Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.

$ \Psi $ is defined by

\[ \Psi := \left( \frac{\det\gamma_{ij}}{\det f_{ij}} \right) ^{1/12} \]

Sets the value at the current time step (jtime ) and deletes the value of N.

Definition at line 424 of file time_slice_conf.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::gam_dd_evol, Time_slice::gam_uu_evol, Time_slice::jtime, Time_slice::n_evol, Time_slice::p_gamma, Time_slice_conf::p_ln_psi, Time_slice_conf::p_psi4, Time_slice_conf::psi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

void Time_slice_conf::set_psi_del_npsi ( const Scalar psi_in  )  [virtual, inherited]

Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.

$ \Psi $ is defined by

\[ \Psi := \left( \frac{\det\gamma_{ij}}{\det f_{ij}} \right) ^{1/12} \]

Sets the value at the current time step (jtime ) and deletes the value of $N\Psi$.

Definition at line 400 of file time_slice_conf.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::gam_dd_evol, Time_slice::gam_uu_evol, Time_slice::jtime, Time_slice_conf::npsi_evol, Time_slice::p_gamma, Time_slice_conf::p_ln_psi, Time_slice_conf::p_psi4, Time_slice_conf::psi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

void Isol_hor::set_radius ( double  rad  )  [inline]

Sets the radius of the horizon to rad .

Definition at line 429 of file isol_hor.h.

References radius.

void Time_slice::set_scheme_order ( int  ord  )  [inline, inherited]

Sets the order of the finite-differences scheme.

Definition at line 327 of file time_slice.h.

References Time_slice::scheme_order.

const Scalar Isol_hor::source_b_tilde (  )  const

Source for b_tilde.

const Vector Isol_hor::source_beta (  )  const
const Scalar Isol_hor::source_nn (  )  const
const Scalar Isol_hor::source_psi (  )  const
const Vector Isol_hor::source_vector_b (  )  const

Source for vector_b.

virtual const Metric& Isol_hor::tgam (  )  const [inline, virtual]

Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $ Returns the value at the current time step (jtime ).

Reimplemented from Time_slice_conf.

Definition at line 512 of file isol_hor.h.

References met_gamt.

const Vector Isol_hor::tradial_vect_hor (  )  const
const Scalar & Time_slice_conf::trk (  )  const [virtual, inherited]
void Isol_hor::update_aa (  ) 
const Vector & Time_slice_conf::vec_X ( int  method_poisson = 6  )  const [virtual, inherited]

Vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $.

(see the documentation of hata_evol)

Definition at line 821 of file time_slice_conf.C.

References Time_slice_conf::ff, Time_slice_conf::hata_evol, Evolution< TyT >::is_known(), Time_slice::jtime, and Time_slice_conf::p_vec_X.

double Isol_hor::viriel_seul (  )  const

Computes the viriel error, that is the difference between the ADM and the Komar masses, calculated by the asymptotic behaviours of respectively $\Psi$ and N .

WARNING this should only be used for an isolated black hole.

const Vector Isol_hor::vv_bound_cart ( double  om  )  const
const Vector Isol_hor::vv_bound_cart_bin ( double  om,
int  hole = 0 
) const

Vector $ V^i $ for boundary conditions in cartesian for binary systems.

Definition at line 1375 of file bound_hor.C.

References Map::get_bvect_cart(), and mp.


Friends And Related Function Documentation

ostream& operator<< ( ostream &  ,
const Time_slice  
) [friend, inherited]

Display.


Member Data Documentation

Evolution_std<Scalar> Time_slice_conf::A_hata_evol [mutable, protected, inherited]

Potential A associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $.

(see documentation of Sym_tensor::p_aaa).

Definition at line 543 of file time_slice.h.

Values at successive time steps of the components $ A^{ij}_{auto} $ of the conformal representation of the traceless part of the extrinsic curvature:.

Definition at line 308 of file isol_hor.h.

Values at successive time steps of the components $ A^{ij}_{comp} $ of the conformal representation of the traceless part of the extrinsic curvature:.

Definition at line 314 of file isol_hor.h.

Evolution_std<Sym_tensor> Isol_hor::aa_nn [mutable, protected]

Values at successive time steps of the components $ A^{ij}*2N $.

Definition at line 318 of file isol_hor.h.

Values at successive time steps of the components $ A^{ij}A_{ij} $.

Definition at line 321 of file isol_hor.h.

Evolution_full<Tbl> Time_slice::adm_mass_evol [mutable, protected, inherited]

ADM mass at each time step, since the creation of the slice.

At a given time step j, adm_mass_evol[j] is a 1-D Tbl of size the number nz of domains, containing the "ADM mass" evaluated at the outer boundary of each domain. The true ADM mass is thus the last value, i.e. adm_mass_evol[j](nz-1).

Definition at line 229 of file time_slice.h.

Evolution_std<Scalar> Time_slice_conf::B_hata_evol [mutable, protected, inherited]

Potential $ \tilde{B} $ associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $.

(see documentation of Sym_tensor::p_tilde_b).

Definition at line 548 of file time_slice.h.

Values at successive time steps of the shift function $ \beta^i_{auto} $.

Definition at line 299 of file isol_hor.h.

Values at successive time steps of the shift function $ \beta^i_{comp} $.

Definition at line 302 of file isol_hor.h.

Evolution_std<Vector> Time_slice::beta_evol [mutable, protected, inherited]

Values at successive time steps of the shift vector $ \beta^i $.

Definition at line 215 of file time_slice.h.

double Isol_hor::boost_x [protected]

Boost velocity in x-direction.

Definition at line 270 of file isol_hor.h.

double Isol_hor::boost_z [protected]

Boost velocity in z-direction.

Definition at line 273 of file isol_hor.h.

Function used to construct $ A^{ij}_{auto} $ from the total $A^{ij}$.

Only used for a binary system.

Mainly this Scalar is 1 around the hole and 0 around the companion and the sum of decouple for the hole and his companion is 1 everywhere.

Definition at line 343 of file isol_hor.h.

int Time_slice::depth [protected, inherited]

Number of stored time slices.

Definition at line 175 of file time_slice.h.

Evolution_std<Vector> Isol_hor::dn_evol [mutable, protected]

Values at successive time steps of the covariant derivative of the lapse with respect to the flat metric $ \overline\nabla_i N $.

Definition at line 292 of file isol_hor.h.

Evolution_std<Vector> Isol_hor::dpsi_evol [mutable, protected]

Values at successive time steps of the covariant derivative of the conformal factor $ \overline\nabla_i \Psi $.

Definition at line 296 of file isol_hor.h.

const Metric_flat& Time_slice_conf::ff [protected, inherited]

Pointer on the flat metric $ f_{ij} $ with respect to which the conformal decomposition is performed.

Definition at line 503 of file time_slice.h.

Evolution_std<Sym_tensor> Time_slice::gam_dd_evol [mutable, protected, inherited]

Values at successive time steps of the covariant components of the induced metric $ \gamma_{ij} $.

Definition at line 194 of file time_slice.h.

Evolution_std<Sym_tensor> Time_slice::gam_uu_evol [mutable, protected, inherited]

Values at successive time steps of the contravariant components of the induced metric $ \gamma^{ij} $.

Definition at line 199 of file time_slice.h.

Time derivative of the 3-metric tilde.

Definition at line 327 of file isol_hor.h.

Evolution_std<Sym_tensor> Time_slice_conf::hata_evol [mutable, protected, inherited]

Values at successive time steps of the components $ \hat{A}^{ij} $.

It is the conformal representation of the traceless part of the extrinsic curvature: $ \hat{A}^{ij} = \Psi^{10} \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $. One can uniquely (up to the boundary conditions) define the decomposition: $ \hat{A}^{ij} = {\cal D}^i X^j + {\cal D}^j X^i - \frac{2}{3} {\cal D}_k X^k f^{ij} + \hat{A}^{ij}_{TT} $, where $ X^i $ represents the longitudinal part and $ \hat{A}^{ij}_{TT} $ is the transverse-traceless part.

Definition at line 538 of file time_slice.h.

Evolution_std<Sym_tensor> Time_slice_conf::hh_evol [mutable, protected, inherited]

Values at successive time steps of the components $ h^{ij} $.

It is the deviation of the conformal metric $ \tilde\gamma^{ij} $ from the flat metric $ f^{ij} $: $\tilde\gamma^{ij} = f^{ij} + h^{ij} $.

Definition at line 526 of file time_slice.h.

int Time_slice::jtime [protected, inherited]

Time step index of the latest slice.

Definition at line 186 of file time_slice.h.

Evolution_std<Sym_tensor> Time_slice::k_dd_evol [mutable, protected, inherited]

Values at successive time steps of the covariant components of the extrinsic curvature tensor $ K_{ij} $.

Definition at line 204 of file time_slice.h.

Evolution_std<Sym_tensor> Time_slice::k_uu_evol [mutable, protected, inherited]

Values at successive time steps of the contravariant components of the extrinsic curvature tensor $ K^{ij} $.

Definition at line 209 of file time_slice.h.

3 metric tilde

Definition at line 324 of file isol_hor.h.

Map_af& Isol_hor::mp [protected]

Affine mapping.

Definition at line 258 of file isol_hor.h.

Values at successive time steps of the lapse function $ N_{auto} $.

Definition at line 279 of file isol_hor.h.

Values at successive time steps of the lapse function $ N_{comp} $.

Definition at line 282 of file isol_hor.h.

Evolution_std<Scalar> Time_slice::n_evol [mutable, protected, inherited]

Values at successive time steps of the lapse function N.

Definition at line 212 of file time_slice.h.

Evolution_std<Scalar> Time_slice_conf::npsi_evol [mutable, protected, inherited]

Values at successive time steps of the factor $ N\Psi $.

Definition at line 518 of file time_slice.h.

int Isol_hor::nz [protected]

Number of zones.

Definition at line 261 of file isol_hor.h.

double Isol_hor::omega [protected]

Angular velocity in LORENE's units.

Definition at line 267 of file isol_hor.h.

Metric* Time_slice::p_gamma [mutable, protected, inherited]

Pointer on the induced metric at the current time step (jtime).

Definition at line 235 of file time_slice.h.

Vector* Time_slice_conf::p_hdirac [mutable, protected, inherited]

Pointer on the vector $ H^i = {\cal D}_j \tilde\gamma^{ij} $ (which vanishes in Dirac gauge), at the current time step (jtime).

Definition at line 567 of file time_slice.h.

Scalar* Time_slice_conf::p_ln_psi [mutable, protected, inherited]

Pointer on the logarithm of $ \Psi $ at the current time step (jtime).

Definition at line 562 of file time_slice.h.

Scalar* Time_slice_conf::p_psi4 [mutable, protected, inherited]

Pointer on the factor $ \Psi^4 $ at the current time step (jtime).

Definition at line 559 of file time_slice.h.

Metric* Time_slice_conf::p_tgamma [mutable, protected, inherited]

Pointer on the conformal metric $ \tilde\gamma_{ij} $ at the current time step (jtime).

Definition at line 556 of file time_slice.h.

Vector* Time_slice_conf::p_vec_X [mutable, protected, inherited]

Pointer on the vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $.

(see the documentation of hata_evol)

Definition at line 573 of file time_slice.h.

Values at successive time steps of the conformal factor $ \Psi_{auto} $.

Definition at line 285 of file isol_hor.h.

Values at successive time steps of the lapse function $ \Psi_{comp} $.

Definition at line 288 of file isol_hor.h.

Evolution_std<Scalar> Time_slice_conf::psi_evol [mutable, protected, inherited]

Values at successive time steps of the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.

$ \Psi $ is defined by

\[ \Psi := \left( \frac{\det\gamma_{ij}}{\det f_{ij}} \right) ^{1/12} \]

Definition at line 513 of file time_slice.h.

double Isol_hor::radius [protected]

Radius of the horizon in LORENE's units.

Definition at line 264 of file isol_hor.h.

double Isol_hor::regul [protected]

Intensity of the correction on the shift vector.

Definition at line 276 of file isol_hor.h.

int Time_slice::scheme_order [protected, inherited]

Order of the finite-differences scheme for the computation of time derivatives.

This order is not constant and can be adjusted via set_scheme_order() .

Definition at line 183 of file time_slice.h.

Evolution_std<double> Time_slice::the_time [protected, inherited]

Time label of each slice.

Definition at line 189 of file time_slice.h.

Scalar Isol_hor::trK [protected]

Trace of the extrinsic curvature.

Definition at line 330 of file isol_hor.h.

Evolution_std<Scalar> Time_slice::trk_evol [mutable, protected, inherited]

Values at successive time steps of the trace K of the extrinsic curvature.

Definition at line 220 of file time_slice.h.

Time derivative of the trace of the extrinsic curvature.

Definition at line 333 of file isol_hor.h.


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

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1