LORENE
|
Spacelike time-slice of an Isolated Horizon in a 3+1 spacetime with conformal decomposition. More...
#include <isol_hor.h>
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_af & | get_mp () const |
Returns the mapping (readonly). | |
Map_af & | set_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 Scalar & | n_auto () const |
Lapse function ![]() jtime . | |
virtual const Scalar & | n_comp () const |
Lapse function ![]() jtime . | |
virtual const Scalar & | psi_auto () const |
Conformal factor ![]() jtime . | |
virtual const Scalar & | psi_comp () const |
Conformal factor ![]() jtime . | |
virtual const Vector & | dnn () const |
Covariant derivative of the lapse function ![]() jtime | |
virtual const Vector & | dpsi () const |
Covariant derivative with respect to the flat metric of the conformal factor ![]() jtime | |
virtual const Vector & | beta_auto () const |
Shift function ![]() jtime . | |
virtual const Vector & | beta_comp () const |
Shift function ![]() jtime . | |
virtual const Sym_tensor & | aa_auto () const |
Conformal representation ![]() jtime . | |
virtual const Sym_tensor & | aa_comp () const |
Conformal representation ![]() jtime . | |
virtual const Scalar & | aa_quad () const |
Conformal representation ![]() | |
virtual const Metric & | tgam () const |
Conformal metric ![]() 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 ![]() comp . | |
void | beta_comp (const Isol_hor &comp) |
Imports the part of ![]() 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 ![]() | |
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 ![]() ![]() ![]() | |
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 ( ![]() | |
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 ![]() | |
const Scalar | source_nn () const |
Source for N . | |
const Vector | source_beta () const |
Source for ![]() | |
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 ![]() | |
const Valeur | boundary_psi_Neu_evol () const |
Neumann boundary condition for ![]() | |
const Valeur | boundary_psi_Dir_spat () const |
Dirichlet boundary condition for ![]() | |
const Valeur | boundary_psi_Neu_spat () const |
Neumann boundary condition for ![]() | |
const Valeur | boundary_psi_app_hor () const |
Neumann boundary condition for ![]() | |
const Valeur | boundary_psi_Dir () const |
Dirichlet boundary condition for ![]() | |
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) ![]() | |
const Valeur | boundary_nn_Dir_lapl (int mer=1) const |
Dirichlet boundary condition for N fixing the divergence of the connection form ![]() | |
const Valeur | boundary_nn_Neu_eff (double aa) const |
Neumann boundary condition on nn (effectif) ![]() | |
const Valeur | boundary_nn_Dir (double aa) const |
Dirichlet boundary condition ![]() | |
const Valeur | boundary_beta_r () const |
Component r of boundary value of ![]() | |
const Valeur | boundary_beta_theta () const |
Component theta of boundary value of ![]() | |
const Valeur | boundary_beta_phi (double om) const |
Component phi of boundary value of ![]() | |
const Valeur | boundary_beta_x (double om) const |
Component x of boundary value of ![]() | |
const Valeur | boundary_beta_y (double om) const |
Component y of boundary value of ![]() | |
const Valeur | boundary_beta_z () const |
Component z of boundary value of ![]() | |
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 ![]() | |
const Vector | vv_bound_cart_bin (double om, int hole=0) const |
Vector ![]() | |
const Valeur | boundary_vv_x (double om) const |
Component x of boundary value of ![]() | |
const Valeur | boundary_vv_y (double om) const |
Component y of boundary value of ![]() | |
const Valeur | boundary_vv_z (double om) const |
Component z of boundary value of ![]() | |
const Valeur | boundary_vv_x_bin (double om, int hole=0) const |
Component x of boundary value of ![]() | |
const Valeur | boundary_vv_y_bin (double om, int hole=0) const |
Component y of boundary value of ![]() | |
const Valeur | boundary_vv_z_bin (double om, int hole=0) const |
Component z of boundary value of ![]() | |
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 ![]() ![]() | |
double | regularisation (const Vector &shift_auto, const Vector &shift_comp, double ang_vel) |
Corrects shift_auto in such a way that the total ![]() ![]() | |
double | regularise_one () |
Corrects the shift in the innermost shell, so that it remains ![]() ![]() | |
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 ![]() ![]() ![]() | |
virtual void | set_psi_del_n (const Scalar &psi_in) |
Sets the conformal factor ![]() ![]() ![]() | |
virtual void | set_npsi_del_psi (const Scalar &npsi_in) |
Sets the factor ![]() jtime ) and deletes the value of ![]() | |
virtual void | set_npsi_del_n (const Scalar &npsi_in) |
Sets the factor ![]() jtime ) and deletes the value of N. | |
virtual void | set_hh (const Sym_tensor &hh_in) |
Sets the deviation ![]() ![]() ![]() ![]() | |
virtual void | set_hata (const Sym_tensor &hata_in) |
Sets the conformal representation ![]() ![]() | |
virtual void | set_hata_TT (const Sym_tensor_tt &hata_tt) |
Sets the TT part of ![]() hata_evol ). | |
virtual void | set_hata_from_XAB (Param *par_bc=0x0, Param *par_mat=0x0) |
Sets the conformal representation ![]() ![]() ![]() | |
virtual const Scalar & | nn () const |
Lapse function N at the current time step (jtime ) | |
virtual const Sym_tensor & | gam_dd () const |
Induced metric (covariant components ![]() jtime ) | |
virtual const Sym_tensor & | gam_uu () const |
Induced metric (contravariant components ![]() jtime ) | |
virtual const Sym_tensor & | k_dd () const |
Extrinsic curvature tensor (covariant components ![]() jtime ) | |
virtual const Sym_tensor & | k_uu () const |
Extrinsic curvature tensor (contravariant components ![]() jtime ) | |
virtual const Scalar & | A_hata () const |
Returns the potential A of ![]() | |
virtual const Scalar & | B_hata () const |
Returns the potential ![]() ![]() | |
virtual const Scalar & | psi () const |
Conformal factor ![]() ![]() ![]() | |
const Scalar & | psi4 () const |
Factor ![]() jtime ). | |
const Scalar & | ln_psi () const |
Logarithm of ![]() jtime ). | |
virtual const Scalar & | npsi () const |
Factor ![]() jtime ). | |
virtual const Sym_tensor & | hh (Param *=0x0, Param *=0x0) const |
Deviation ![]() ![]() ![]() ![]() | |
virtual const Sym_tensor & | hata () const |
Conformal representation ![]() ![]() | |
virtual Sym_tensor | aa () const |
Conformal representation ![]() ![]() | |
virtual const Scalar & | trk () const |
Trace K of the extrinsic curvature at the current time step (jtime ) | |
virtual const Vector & | hdirac () const |
Vector ![]() | |
virtual const Vector & | vec_X (int method_poisson=6) const |
Vector ![]() ![]() | |
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 ![]() ![]() | |
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 Vector & | beta () const |
shift vector ![]() jtime ) | |
const Metric & | gam () const |
Induced metric ![]() 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_af & | mp |
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< Scalar > | n_auto_evol |
Values at successive time steps of the lapse function ![]() | |
Evolution_std< Scalar > | n_comp_evol |
Values at successive time steps of the lapse function ![]() | |
Evolution_std< Scalar > | psi_auto_evol |
Values at successive time steps of the conformal factor ![]() | |
Evolution_std< Scalar > | psi_comp_evol |
Values at successive time steps of the lapse function ![]() | |
Evolution_std< Vector > | dn_evol |
Values at successive time steps of the covariant derivative of the lapse with respect to the flat metric ![]() | |
Evolution_std< Vector > | dpsi_evol |
Values at successive time steps of the covariant derivative of the conformal factor ![]() | |
Evolution_std< Vector > | beta_auto_evol |
Values at successive time steps of the shift function ![]() | |
Evolution_std< Vector > | beta_comp_evol |
Values at successive time steps of the shift function ![]() | |
Evolution_std< Sym_tensor > | aa_auto_evol |
Values at successive time steps of the components ![]() | |
Evolution_std< Sym_tensor > | aa_comp_evol |
Values at successive time steps of the components ![]() | |
Evolution_std< Sym_tensor > | aa_nn |
Values at successive time steps of the components ![]() | |
Evolution_std< Scalar > | aa_quad_evol |
Values at successive time steps of the components ![]() | |
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 ![]() ![]() | |
const Metric_flat & | ff |
Pointer on the flat metric ![]() | |
Evolution_std< Scalar > | psi_evol |
Values at successive time steps of the conformal factor ![]() ![]() ![]() | |
Evolution_std< Scalar > | npsi_evol |
Values at successive time steps of the factor ![]() | |
Evolution_std< Sym_tensor > | hh_evol |
Values at successive time steps of the components ![]() | |
Evolution_std< Sym_tensor > | hata_evol |
Values at successive time steps of the components ![]() | |
Evolution_std< Scalar > | A_hata_evol |
Potential A associated with the symmetric tensor ![]() | |
Evolution_std< Scalar > | B_hata_evol |
Potential ![]() ![]() | |
Metric * | p_tgamma |
Pointer on the conformal metric ![]() jtime ) | |
Scalar * | p_psi4 |
Pointer on the factor ![]() jtime ) | |
Scalar * | p_ln_psi |
Pointer on the logarithm of ![]() jtime ) | |
Vector * | p_hdirac |
Pointer on the vector ![]() jtime ). | |
Vector * | p_vec_X |
Pointer on the vector ![]() ![]() | |
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_tensor > | gam_dd_evol |
Values at successive time steps of the covariant components of the induced metric ![]() | |
Evolution_std< Sym_tensor > | gam_uu_evol |
Values at successive time steps of the contravariant components of the induced metric ![]() | |
Evolution_std< Sym_tensor > | k_dd_evol |
Values at successive time steps of the covariant components of the extrinsic curvature tensor ![]() | |
Evolution_std< Sym_tensor > | k_uu_evol |
Values at successive time steps of the contravariant components of the extrinsic curvature tensor ![]() | |
Evolution_std< Scalar > | n_evol |
Values at successive time steps of the lapse function N. | |
Evolution_std< Vector > | beta_evol |
Values at successive time steps of the shift vector ![]() | |
Evolution_std< Scalar > | trk_evol |
Values at successive time steps of the trace K of the extrinsic curvature. | |
Evolution_full< Tbl > | adm_mass_evol |
ADM mass at each time step, since the creation of the slice. | |
Metric * | p_gamma |
Pointer on the induced metric at the current time step (jtime ) | |
Friends | |
class | Bin_hor |
Spacelike time-slice of an Isolated Horizon in a 3+1 spacetime with conformal decomposition.
No gauge choice imposed. ()
Definition at line 254 of file isol_hor.h.
Standard constructor.
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 176 of file isol_hor.C.
Lorene::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.
mpi | affine mapping |
lapse_in | lapse function N |
psi_in | conformal factor ![]() ![]() ![]() |
shift_in | shift vector |
aa_in | conformal representation ![]() ![]() |
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 193 of file isol_hor.C.
References Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, met_gamt, Lorene::Time_slice::the_time, trK, and Lorene::Time_slice::trk_evol.
Copy constructor.
Definition at line 222 of file isol_hor.C.
Constructor from a binary file.
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 253 of file isol_hor.C.
References beta_auto_evol, boost_x, boost_z, Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::Time_slice::depth, Lorene::Time_slice_conf::ff, Lorene::fread_be(), gamt_point, Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, met_gamt, mp, n_auto_evol, omega, psi_auto_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::the_time, trK, Lorene::Time_slice::trk_evol, and trK_point.
|
virtual |
Destructor.
Definition at line 336 of file isol_hor.C.
Returns the potential A of
See the documentation of Sym_tensor
for details. Returns the value at the current time step (jtime
).
Definition at line 664 of file time_slice_conf.C.
References Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::hata_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
virtualinherited |
Conformal representation
Returns the value at the current time step (jtime
).
Definition at line 765 of file time_slice_conf.C.
References Lorene::Time_slice_conf::hata(), Lorene::Time_slice_conf::psi(), and Lorene::Time_slice_conf::psi4().
|
virtual |
Conformal representation jtime
.
Definition at line 503 of file isol_hor.C.
References aa_auto_evol, and Lorene::Time_slice::jtime.
|
virtual |
Conformal representation jtime
.
Definition at line 509 of file isol_hor.C.
References aa_comp_evol, and Lorene::Time_slice::jtime.
Definition at line 918 of file isol_hor.C.
Conformal representation
Returns the value at the current time step jtime
.
Definition at line 884 of file isol_hor.C.
References Lorene::Time_slice_conf::aa(), aa_quad_evol, Lorene::contract(), Lorene::Time_slice::jtime, met_gamt, Lorene::Time_slice_conf::psi4(), and Lorene::Time_slice::the_time.
Definition at line 1125 of file isol_hor.C.
|
virtualinherited |
Returns the ADM mass (geometrical units) at the current step.
Moreover this method updates adm_mass_evol
if necessary.
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 74 of file tslice_adm_mass.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Vector::flux(), Lorene::Time_slice::gam_dd(), Lorene::Map::get_mg(), Lorene::Time_slice::jtime, Lorene::Time_slice::the_time, Lorene::Evolution_std< TyT >::update(), and Lorene::Map::val_r().
double Lorene::Isol_hor::ang_mom_adm | ( | ) | const |
ADM angular Momentum
Definition at line 248 of file phys_param.C.
References Lorene::Time_slice_conf::gam_dd(), Lorene::Map_af::integrale_surface_infini(), Lorene::Time_slice_conf::k_dd(), mp, Lorene::Scalar::mult_rsint(), and Lorene::Time_slice_conf::trk().
double Lorene::Isol_hor::ang_mom_hor | ( | ) | const |
Angular momentum (modulo)
Definition at line 178 of file phys_param.C.
References Lorene::contract(), darea_hor(), Lorene::Time_slice_conf::ff, Lorene::Metric::get_mp(), Lorene::Metric_flat::get_triad(), Lorene::Map_af::integrale_surface(), Lorene::Time_slice_conf::k_dd(), mp, Lorene::Scalar::mult_rsint(), radial_vect_hor(), radius, Lorene::Vector::set(), and Lorene::Scalar::std_spectral_base().
double Lorene::Isol_hor::area_hor | ( | ) | const |
Area of the horizon.
Definition at line 157 of file phys_param.C.
References darea_hor(), Lorene::Map_af::integrale_surface(), mp, Lorene::Scalar::raccord(), and radius.
double Lorene::Isol_hor::axi_break | ( | ) | const |
Breaking of the axial symmetry on the horizon.
Definition at line 1087 of file isol_hor.C.
References Lorene::contract(), darea_hor(), Lorene::Sym_tensor::derive_lie(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam(), Lorene::Time_slice_conf::gam_uu(), Lorene::Metric::get_mp(), Lorene::Metric_flat::get_triad(), Lorene::Map_af::integrale_surface(), mp, Lorene::Scalar::mult_rsint(), radius_hor(), Lorene::Vector::set(), Lorene::Scalar::std_spectral_base(), and Lorene::Tensor::up_down().
Returns the potential
See the documentation of Sym_tensor_tt
for details. Returns the value at the current time step (jtime
).
Definition at line 678 of file time_slice_conf.C.
References Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Time_slice_conf::hata_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
Radial component of the shift with respect to the conformal metric.
Definition at line 136 of file phys_param.C.
References Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Tensor::down(), met_gamt, and Lorene::Metric::radial_vect().
shift vector jtime
)
Definition at line 87 of file time_slice_access.C.
References Lorene::Time_slice::beta_evol, Lorene::Evolution< TyT >::is_known(), and Lorene::Time_slice::jtime.
Shift function jtime
.
Definition at line 491 of file isol_hor.C.
References beta_auto_evol, and Lorene::Time_slice::jtime.
Boundary value for a boost in x-direction.
Definition at line 1144 of file bound_hor.C.
References boost_x, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), mp, Lorene::Valeur::set_base(), and Lorene::Mg3d::std_base_vect_cart().
Boundary value for a boost in z-direction.
Definition at line 1155 of file bound_hor.C.
References boost_z, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), mp, Lorene::Valeur::set_base(), and Lorene::Mg3d::std_base_vect_cart().
Shift function jtime
.
Definition at line 497 of file isol_hor.C.
References beta_comp_evol, and Lorene::Time_slice::jtime.
Imports the part of comp
.
The total
Definition at line 712 of file isol_hor.C.
References beta_auto(), beta_comp(), beta_comp_evol, Lorene::Time_slice::beta_evol, Lorene::Vector::change_triad(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Tensor::get_triad(), Lorene::Scalar::import(), Lorene::Time_slice::jtime, mp, Lorene::Vector::set(), Lorene::Vector::std_spectral_base(), and Lorene::Time_slice::the_time.
Dirichlet boundary condition for b_tilde
.
Definition at line 1212 of file bound_hor.C.
References b_tilde(), Lorene::contract(), Lorene::Tensor::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::inc_dzpuis(), met_gamt, mp, Lorene::Metric::radial_vect(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().
Neumann boundary condition for b_tilde
.
Definition at line 1169 of file bound_hor.C.
References b_tilde(), Lorene::contract(), Lorene::Scalar::derive_cov(), Lorene::Tensor::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::inc_dzpuis(), met_gamt, mp, Lorene::Metric::radial_vect(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().
Component phi of boundary value of
Definition at line 978 of file bound_hor.C.
References Lorene::Valeur::base, Lorene::Mg3d::get_angu(), Lorene::Base_val::get_base_p(), Lorene::Base_val::get_base_r(), Lorene::Base_val::get_base_t(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::get_spectral_va(), mp, Lorene::Scalar::mult_rsint(), Lorene::Time_slice_conf::nn(), radial_vect_hor(), Lorene::Valeur::set(), Lorene::Valeur::set_base_p(), Lorene::Valeur::set_base_r(), Lorene::Valeur::set_base_t(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().
Component r of boundary value of
Definition at line 897 of file bound_hor.C.
References Lorene::Valeur::base, Lorene::Mg3d::get_angu(), Lorene::Valeur::get_base(), Lorene::Base_val::get_base_p(), Lorene::Base_val::get_base_r(), Lorene::Base_val::get_base_t(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::get_spectral_va(), mp, Lorene::Time_slice_conf::nn(), Lorene::norme(), radial_vect_hor(), Lorene::Valeur::set(), Lorene::Valeur::set_base_p(), Lorene::Valeur::set_base_r(), Lorene::Valeur::set_base_t(), and Lorene::Scalar::val_grid_point().
Component theta of boundary value of
Definition at line 938 of file bound_hor.C.
References Lorene::Valeur::base, Lorene::Mg3d::get_angu(), Lorene::Base_val::get_base_p(), Lorene::Base_val::get_base_r(), Lorene::Base_val::get_base_t(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::nn(), radial_vect_hor(), Lorene::Valeur::set(), Lorene::Valeur::set_base_p(), Lorene::Valeur::set_base_r(), Lorene::Valeur::set_base_t(), and Lorene::Scalar::val_grid_point().
Component x of boundary value of
Definition at line 1021 of file bound_hor.C.
References Lorene::Vector::change_triad(), Lorene::Map::cosp, Lorene::Map::cost, Lorene::Time_slice::gam(), Lorene::Mg3d::get_angu(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Map::get_rot_phi(), mp, Lorene::Time_slice_conf::nn(), Lorene::Metric::radial_vect(), Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mtbl::set_etat_qcq(), Lorene::Map::sinp, Lorene::Map::sint, Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and Lorene::Map::ya.
Component y of boundary value of
Definition at line 1071 of file bound_hor.C.
References Lorene::Vector::change_triad(), Lorene::Map::cosp, Lorene::Map::cost, Lorene::Time_slice::gam(), Lorene::Mg3d::get_angu(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Map::get_rot_phi(), mp, Lorene::Time_slice_conf::nn(), Lorene::Metric::radial_vect(), Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mtbl::set_etat_qcq(), Lorene::Map::sinp, Lorene::Map::sint, Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and Lorene::Map::xa.
Component z of boundary value of
Definition at line 1121 of file bound_hor.C.
References Lorene::Vector::change_triad(), Lorene::Time_slice::gam(), Lorene::Mg3d::get_angu(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::nn(), Lorene::Metric::radial_vect(), Lorene::Valeur::set(), Lorene::Valeur::set_base(), and Lorene::Mg3d::std_base_vect_cart().
Dirichlet boundary condition
Definition at line 647 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::nn(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), and Lorene::Scalar::val_grid_point().
Dirichlet boundary condition for N
(effectif)
Definition at line 586 of file bound_hor.C.
References Lorene::Scalar::dec_dzpuis(), Lorene::Scalar::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::nn(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), and Lorene::Scalar::val_grid_point().
Dirichlet boundary condition for N
using the extrinsic curvature.
Definition at line 371 of file bound_hor.C.
References Lorene::contract(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice_conf::k_dd(), mp, Lorene::Time_slice_conf::nn(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().
Dirichlet boundary condition for N
fixing the divergence of the connection form
Definition at line 693 of file bound_hor.C.
References Lorene::Metric_flat::con(), Lorene::contract(), Lorene::Metric_flat::cov(), Lorene::Scalar::dec_dzpuis(), Lorene::Tensor::derive_con(), Lorene::Scalar::derive_cov(), Lorene::Tensor::derive_cov(), Lorene::Vector::divergence(), Lorene::Tensor::down(), Lorene::exp(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam(), Lorene::Time_slice_conf::gam_uu(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::inc_dzpuis(), Lorene::Map_af::integrale_surface(), Lorene::Time_slice_conf::k_dd(), Lorene::Scalar::lapang(), Lorene::log(), met_gamt, mp, Lorene::Time_slice_conf::nn(), Lorene::Scalar::poisson_angu(), Lorene::pow(), Lorene::Time_slice_conf::psi(), Lorene::Map::r, Lorene::Metric::radial_vect(), radius, Lorene::Valeur::set(), Lorene::Scalar::set_dzpuis(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), Lorene::Time_slice_conf::trk(), Lorene::Tensor::up_down(), and Lorene::Scalar::val_grid_point().
Neumann boundary condition for N
using Cook's boundary condition.
Definition at line 489 of file bound_hor.C.
References Lorene::contract(), Lorene::Scalar::derive_cov(), Lorene::Vector::divergence(), Lorene::Tensor::down(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam(), Lorene::Time_slice_conf::gam_uu(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice_conf::k_dd(), Lorene::Scalar::lapang(), Lorene::log(), met_gamt, mp, Lorene::Time_slice_conf::nn(), Lorene::pow(), Lorene::Time_slice_conf::psi(), Lorene::Map::r, Lorene::Metric::radial_vect(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), Lorene::Time_slice_conf::trk(), and Lorene::Scalar::val_grid_point().
Neumann boundary condition on nn (effectif)
Definition at line 622 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::nn(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), and Lorene::Scalar::val_grid_point().
Neumann boundary condition for N
using the extrinsic curvature.
Definition at line 420 of file bound_hor.C.
References aa_auto(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Scalar::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam(), Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice_conf::k_dd(), kappa_hor(), met_gamt, mp, Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi(), Lorene::Metric::radial_vect(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), Lorene::Time_slice_conf::trk(), and Lorene::Scalar::val_grid_point().
Neumann boundary condition for
Definition at line 297 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::psi(), radius, Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), and Lorene::Scalar::val_grid_point().
Dirichlet boundary condition for
Definition at line 324 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::psi(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), and Lorene::Scalar::val_grid_point().
Dirichlet boundary condition for
Definition at line 174 of file bound_hor.C.
References Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Time_slice_conf::ff, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Time_slice_conf::trk(), and Lorene::Scalar::val_grid_point().
Dirichlet boundary condition for
Definition at line 231 of file bound_hor.C.
References Lorene::contract(), Lorene::Vector::divergence(), Lorene::Time_slice_conf::ff, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice_conf::k_dd(), mp, Lorene::Time_slice_conf::psi(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Scalar::std_spectral_base(), tradial_vect_hor(), Lorene::Time_slice_conf::trk(), and Lorene::Scalar::val_grid_point().
Neumann boundary condition for
Definition at line 203 of file bound_hor.C.
References Lorene::Time_slice::beta(), Lorene::Scalar::derive_cov(), Lorene::Vector::divergence(), Lorene::Time_slice_conf::ff, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), Lorene::Time_slice_conf::trk(), and Lorene::Scalar::val_grid_point().
Neumann boundary condition for
Definition at line 267 of file bound_hor.C.
References Lorene::contract(), Lorene::Scalar::derive_cov(), Lorene::Vector::divergence(), Lorene::Time_slice_conf::ff, Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Time_slice_conf::k_dd(), mp, Lorene::Time_slice_conf::psi(), Lorene::Valeur::set(), Lorene::Valeur::std_base_scal(), tradial_vect_hor(), Lorene::Time_slice_conf::trk(), and Lorene::Scalar::val_grid_point().
Component x of boundary value of
Definition at line 1493 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and vv_bound_cart().
Component x of boundary value of
Definition at line 1572 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and vv_bound_cart_bin().
Component y of boundary value of
Definition at line 1521 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and vv_bound_cart().
Component y of boundary value of
Definition at line 1599 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and vv_bound_cart_bin().
Component z of boundary value of
Definition at line 1547 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and vv_bound_cart().
Component z of boundary value of
Definition at line 1625 of file bound_hor.C.
References Lorene::Mg3d::get_angu(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), mp, Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Scalar::val_grid_point(), and vv_bound_cart_bin().
|
inherited |
Checks the level at which the dynamical equations are verified.
strain_tensor | : a pointer on the strain_tensor ![]() ![]() ![]() |
energy_density | : a pointer on the energy density E (see check_hamiltonian_constraint ) |
ost | : output stream for a formatted output of the result |
Definition at line 139 of file tslice_check_einstein.C.
References Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Time_slice::gam(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd(), Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu(), Lorene::maxabs(), Lorene::Time_slice::nn(), Lorene::Evolution< TyT >::position(), Lorene::Time_slice::scheme_order, Lorene::Evolution< TyT >::time_derive(), and Lorene::Time_slice::trk().
|
inherited |
Checks the level at which the hamiltonian constraint is verified.
energy_density | : a pointer on the energy density E measured by the Eulerian observer of 4-velocity ![]() |
ost | : output stream for a formatted output of the result |
Definition at line 79 of file tslice_check_einstein.C.
References Lorene::contract(), Lorene::Time_slice::gam(), Lorene::Time_slice::k_dd(), Lorene::Time_slice::k_uu(), Lorene::maxabs(), and Lorene::Time_slice::trk().
|
inherited |
Checks the level at which the momentum constraints are verified.
momentum_density | : a pointer on the momentum density ![]() ![]() ![]() |
ost | : output stream for a formatted output of the result |
Definition at line 109 of file tslice_check_einstein.C.
References Lorene::Time_slice::gam(), Lorene::Time_slice::k_uu(), Lorene::maxabs(), and Lorene::Time_slice::trk().
|
inherited |
Computes the vector
Definition at line 836 of file time_slice_conf.C.
References Lorene::abs(), Lorene::contract(), Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice_conf::hh_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::max(), Lorene::Vector::ope_killing_conf(), Lorene::Time_slice_conf::p_vec_X, Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi4(), Lorene::Time_slice_conf::tgam(), Lorene::Time_slice::the_time, Lorene::Time_slice_conf::trk(), and Lorene::Evolution_std< TyT >::update().
Element of area of the horizon.
Definition at line 146 of file phys_param.C.
References Lorene::Time_slice_conf::gam_dd(), Lorene::sqrt(), and Lorene::Scalar::std_spectral_base().
|
protectedvirtualinherited |
Deletes all the derived quantities.
Reimplemented from Lorene::Time_slice.
Definition at line 347 of file time_slice_conf.C.
References Lorene::Time_slice::del_deriv(), Lorene::Time_slice_conf::p_hdirac, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::p_tgamma, Lorene::Time_slice_conf::p_vec_X, and Lorene::Time_slice_conf::set_der_0x0().
Covariant derivative of the lapse function jtime
Definition at line 479 of file isol_hor.C.
References dn_evol, and Lorene::Time_slice::jtime.
Covariant derivative with respect to the flat metric of the conformal factor jtime
Definition at line 485 of file isol_hor.C.
References dpsi_evol, and Lorene::Time_slice::jtime.
Scalar Lorene::Isol_hor::expansion | ( | ) | const |
Expansion of the outgoing null normal (
Definition at line 263 of file phys_param.C.
References Lorene::contract(), Lorene::Time_slice::gam(), Lorene::Time_slice_conf::k_dd(), and Lorene::Time_slice_conf::trk().
Induced metric jtime
)
Definition at line 95 of file time_slice_access.C.
References Lorene::Time_slice::gam_dd(), and Lorene::Time_slice::p_gamma.
|
virtualinherited |
Induced metric (covariant components jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 608 of file time_slice_conf.C.
References Lorene::Time_slice::gam_dd_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice_conf::psi4(), Lorene::Time_slice_conf::tgam(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
virtualinherited |
Induced metric (contravariant components jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 619 of file time_slice_conf.C.
References Lorene::Time_slice::gam_uu_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice_conf::psi4(), Lorene::Time_slice_conf::tgam(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
inline |
Returns the boost velocity in x-direction.
Definition at line 445 of file isol_hor.h.
References boost_x.
|
inline |
Returns the boost velocity in z-direction.
Definition at line 454 of file isol_hor.h.
References boost_z.
Returns the function used to construct tkij_auto
from tkij_tot
.
Definition at line 519 of file isol_hor.h.
References decouple.
|
inlineinherited |
Gets the latest value of time step index.
Definition at line 343 of file time_slice.h.
References Lorene::Time_slice::jtime.
|
inline |
|
inline |
|
inlineinherited |
Gets the order of the finite-differences scheme.
Definition at line 340 of file time_slice.h.
References Lorene::Time_slice::scheme_order.
|
inlineinherited |
Gets the time coordinate t at successive time steps.
Definition at line 346 of file time_slice.h.
References Lorene::Time_slice::the_time.
|
virtualinherited |
Conformal representation
Returns the value at the current time step (jtime
).
Definition at line 772 of file time_slice_conf.C.
References Lorene::Time_slice::beta(), Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice_conf::hh(), Lorene::Time_slice_conf::hh_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi4(), Lorene::Time_slice::scheme_order, Lorene::Time_slice::the_time, Lorene::Evolution< TyT >::time_derive(), and Lorene::Evolution_std< TyT >::update().
Vector
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 815 of file time_slice_conf.C.
References Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::hh(), and Lorene::Time_slice_conf::p_hdirac.
|
virtualinherited |
Deviation
Returns the value at the current time step (jtime
).
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 758 of file time_slice_conf.C.
References Lorene::Time_slice_conf::hh_evol, Lorene::Evolution< TyT >::is_known(), and Lorene::Time_slice::jtime.
void Lorene::Isol_hor::init_bhole | ( | ) |
Sets the values of the fields to :
n_auto
n_comp
psi_auto
psi_comp
a being the radius of the hole, the other fields being set to zero.
Definition at line 733 of file isol_hor.C.
References Lorene::Scalar::annule(), beta_auto_evol, beta_comp_evol, Lorene::Time_slice::beta_evol, dn_evol, dpsi_evol, Lorene::Time_slice_conf::ff, Lorene::Map::get_bvect_spher(), Lorene::Time_slice::jtime, mp, n_auto(), n_auto_evol, n_comp(), n_comp_evol, Lorene::Time_slice::n_evol, Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi(), psi_auto(), psi_auto_evol, psi_comp(), psi_comp_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Map::r, Lorene::Scalar::raccord(), radius, Lorene::Vector::set(), Lorene::Scalar::set_dzpuis(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::std_spectral_base(), Lorene::Vector::std_spectral_base(), and Lorene::Time_slice::the_time.
void Lorene::Isol_hor::init_bhole_seul | ( | ) |
Initiates for a single black hole.
WARNING It supposes that the boost is zero and should only be used for an isolated black hole..
Definition at line 797 of file isol_hor.C.
References Lorene::Scalar::annule(), beta_auto_evol, beta_comp_evol, Lorene::Time_slice::beta_evol, dn_evol, dpsi_evol, Lorene::Time_slice_conf::ff, Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Time_slice::jtime, mp, n_auto_evol, n_comp_evol, Lorene::Time_slice::n_evol, Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi(), psi_auto_evol, psi_comp_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Map::r, Lorene::Scalar::raccord(), radius, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::set_outer_boundary(), Lorene::Scalar::std_spectral_base(), and Lorene::Time_slice::the_time.
void Lorene::Isol_hor::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 |
||
) |
Definition at line 153 of file init_data.C.
void Lorene::Isol_hor::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 |
||
) |
Definition at line 1241 of file init_data.C.
void Lorene::Isol_hor::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 |
||
) |
Definition at line 922 of file init_data.C.
void Lorene::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 785 of file isol_hor.C.
References Lorene::Map::flat_met_spher(), gamt_point, met_gamt, mp, Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), trK, and trK_point.
|
virtualinherited |
Extrinsic curvature tensor (covariant components jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 630 of file time_slice_conf.C.
References Lorene::Time_slice::gam(), Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice_conf::k_uu(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
virtualinherited |
Extrinsic curvature tensor (contravariant components jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 643 of file time_slice_conf.C.
References Lorene::Time_slice::gam(), Lorene::Time_slice_conf::hata(), Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_uu_evol, Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi4(), Lorene::Time_slice::the_time, Lorene::Time_slice_conf::trk(), and Lorene::Evolution_std< TyT >::update().
double Lorene::Isol_hor::kappa_hor | ( | ) | const |
Surface gravity
Definition at line 218 of file phys_param.C.
References ang_mom_hor(), Lorene::pow(), radius_hor(), and Lorene::sqrt().
Logarithm of jtime
).
Definition at line 719 of file time_slice_conf.C.
References Lorene::log(), Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::psi(), and Lorene::Scalar::std_spectral_base().
double Lorene::Isol_hor::mass_hor | ( | ) | const |
Mass computed at the horizon
Definition at line 206 of file phys_param.C.
References ang_mom_hor(), Lorene::pow(), radius_hor(), and Lorene::sqrt().
void Lorene::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 895 of file isol_hor.C.
References Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::Metric::cov(), Lorene::Metric::determinant(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam(), Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, met_gamt, Lorene::norme(), Lorene::pow(), set_psi(), Lorene::Scalar::std_spectral_base(), and Lorene::Time_slice::the_time.
Lapse function jtime
.
Definition at line 455 of file isol_hor.C.
References Lorene::Time_slice::jtime, and n_auto_evol.
Lapse function jtime
.
Definition at line 461 of file isol_hor.C.
References Lorene::Time_slice::jtime, and n_comp_evol.
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
Definition at line 577 of file isol_hor.C.
References Lorene::Vector::change_triad(), Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::derive_cov(), dn_evol, Lorene::Time_slice_conf::ff, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Tensor::get_mp(), Lorene::Tensor::get_triad(), Lorene::Scalar::import(), Lorene::Tensor::inc_dzpuis(), Lorene::Time_slice::jtime, mp, n_auto(), n_comp_evol, Lorene::Time_slice::n_evol, Lorene::Scalar::raccord(), Lorene::Vector::set(), Lorene::Valeur::set_base(), Lorene::Tensor::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::std_spectral_base(), and Lorene::Time_slice::the_time.
Lapse function N at the current time step (jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 591 of file time_slice_conf.C.
References Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
Factor jtime
).
Definition at line 732 of file time_slice_conf.C.
References Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
double Lorene::Isol_hor::omega_hor | ( | ) | const |
Orbital velocity
Definition at line 233 of file phys_param.C.
References ang_mom_hor(), Lorene::pow(), radius_hor(), and Lorene::sqrt().
Assignment to another Isol_hor.
Definition at line 343 of file isol_hor.C.
References aa_auto_evol, aa_comp_evol, aa_nn, aa_quad_evol, beta_auto_evol, beta_comp_evol, boost_x, boost_z, decouple, dn_evol, dpsi_evol, gamt_point, met_gamt, mp, n_auto_evol, n_comp_evol, nz, omega, Lorene::Time_slice_conf::operator=(), psi_auto_evol, psi_comp_evol, radius, regul, trK, and trK_point.
Operator >> (virtual function called by the operator<<).
Reimplemented from Lorene::Time_slice.
Definition at line 378 of file isol_hor.C.
References Lorene::Time_slice::adm_mass(), ang_mom_adm(), ang_mom_hor(), area_hor(), boost_x, boost_z, mass_hor(), omega_hor(), Lorene::Time_slice::operator>>(), and radius.
Conformal factor
Returns the value at the current time step (jtime
).
Definition at line 693 of file time_slice_conf.C.
References Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
Factor jtime
).
Definition at line 707 of file time_slice_conf.C.
References Lorene::Time_slice_conf::p_psi4, Lorene::pow(), Lorene::Time_slice_conf::psi(), and Lorene::Scalar::std_spectral_base().
Conformal factor jtime
.
Definition at line 467 of file isol_hor.C.
References Lorene::Time_slice::jtime, and psi_auto_evol.
Conformal factor jtime
.
Definition at line 473 of file isol_hor.C.
References Lorene::Time_slice::jtime, and psi_comp_evol.
Imports the part of comp
.
The total
It also imports the covariant derivative of
Definition at line 669 of file isol_hor.C.
References Lorene::Vector::change_triad(), Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::derive_cov(), dpsi_evol, Lorene::Time_slice_conf::ff, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Tensor::get_mp(), Lorene::Tensor::get_triad(), Lorene::Scalar::import(), Lorene::Tensor::inc_dzpuis(), Lorene::Time_slice::jtime, mp, psi_auto(), psi_comp_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Scalar::raccord(), Lorene::Vector::set(), Lorene::Valeur::set_base(), Lorene::Tensor::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::std_spectral_base(), and Lorene::Time_slice::the_time.
Vector radial normal.
Definition at line 95 of file phys_param.C.
References Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::gam_uu(), Lorene::Metric::get_mp(), Lorene::Metric_flat::get_triad(), Lorene::Vector::set(), Lorene::sqrt(), and Lorene::Vector::std_spectral_base().
double Lorene::Isol_hor::radius_hor | ( | ) | const |
Radius of the horizon.
Definition at line 167 of file phys_param.C.
References area_hor(), and Lorene::pow().
double Lorene::Isol_hor::regularisation | ( | const Vector & | shift_auto, |
const Vector & | shift_comp, | ||
double | ang_vel | ||
) |
Corrects shift_auto
in such a way that the total
WARNING : this should only be used for a black hole in a binary system Bin_hor
.
comp | [input]: the part of ![]() |
Definition at line 88 of file regularisation.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::annule_hard(), beta_auto_evol, Lorene::Vector::change_triad(), Lorene::Valeur::coef_i(), Lorene::diffrelmax(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Map::get_rot_phi(), Lorene::Tensor::get_triad(), Lorene::Scalar::import(), Lorene::Time_slice::jtime, Lorene::norme(), nz, Lorene::pow(), Lorene::Map::r, Lorene::Valeur::set(), Lorene::Vector::set(), Lorene::Valeur::set_base(), Lorene::Valeur::set_etat_c_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Vector::std_spectral_base(), Lorene::Time_slice::the_time, Lorene::Evolution_std< TyT >::update(), Lorene::Map::val_r(), Lorene::Map::xa, and Lorene::Map::ya.
double Lorene::Isol_hor::regularise_one | ( | ) |
Corrects the shift in the innermost shell, so that it remains
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 183 of file regularisation.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::annule_hard(), Lorene::Valeur::base, Lorene::Time_slice::beta(), Lorene::Time_slice::beta_evol, boost_x, boost_z, Lorene::Vector::change_triad(), Lorene::Valeur::coef_i(), Lorene::diffrelmax(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Time_slice::jtime, mp, Lorene::norme(), nz, omega, Lorene::pow(), Lorene::Map::r, Lorene::Valeur::set(), Lorene::Vector::set(), Lorene::Valeur::set_etat_c_qcq(), Lorene::Tensor::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Vector::std_spectral_base(), Lorene::Time_slice::the_time, Lorene::Map_af::val_r(), Lorene::Map::x, and Lorene::Map::y.
Total or partial saves in a binary file.
fich | binary file |
partial_save | indicates whether the whole object must be saved. |
Reimplemented from Lorene::Time_slice.
Definition at line 401 of file isol_hor.C.
References beta_auto_evol, boost_x, boost_z, Lorene::Metric::con(), Lorene::Time_slice::depth, Lorene::fwrite_be(), gamt_point, Lorene::Time_slice::jtime, met_gamt, n_auto_evol, omega, psi_auto_evol, Lorene::Time_slice_conf::psi_evol, Lorene::Scalar::sauve(), Lorene::Tensor_sym::sauve(), Lorene::Time_slice::sauve(), trK, and trK_point.
Saves in a binary file.
The saved data is sufficient to restore the whole time slice via the constructor from file.
rootname | root for the file name; the current time step index will be appended to it. |
Definition at line 461 of file time_slice.C.
References Lorene::Time_slice::beta(), Lorene::Time_slice::depth, Lorene::fwrite_be(), Lorene::Map::get_mg(), Lorene::Time_slice::jtime, Lorene::Time_slice::nn(), Lorene::Base_vect::sauve(), Lorene::Map::sauve(), and Lorene::Time_slice::sauve().
Sets the boost velocity in x-direction to bo
.
Definition at line 449 of file isol_hor.h.
References boost_x.
Sets the boost velocity in z-direction to bo
.
Definition at line 458 of file isol_hor.h.
References boost_z.
|
protectedinherited |
Sets to 0x0
all the pointers on derived quantities.
Definition at line 361 of file time_slice_conf.C.
References Lorene::Time_slice_conf::p_hdirac, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::p_tgamma, and Lorene::Time_slice_conf::p_vec_X.
Sets the conformal metric to gam_tilde.
Definition at line 550 of file isol_hor.C.
References Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, met_gamt, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_tgamma, and Lorene::Time_slice::the_time.
|
virtualinherited |
Sets the conformal representation
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 534 of file time_slice_conf.C.
References Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
virtualinherited |
Sets the conformal representation
These potentials must be up-to-date. It sets the value at the current time step (jtime
).
Definition at line 566 of file time_slice_conf.C.
References Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice_conf::ff, Lorene::Map::get_bvect_spher(), Lorene::Tensor::get_mp(), Lorene::Time_slice_conf::hata_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Vector::ope_killing_conf(), Lorene::Time_slice_conf::p_vec_X, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
virtualinherited |
Sets the TT part of hata_evol
).
Sets the value at current time-step (jtime
) and updates the potentials A and
Definition at line 546 of file time_slice_conf.C.
References Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
virtualinherited |
Sets the deviation
jtime
).
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 489 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, Lorene::max(), Lorene::maxabs(), Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_hdirac, Lorene::Time_slice_conf::p_tgamma, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
|
inline |
Sets the lapse.
Definition at line 540 of file isol_hor.C.
References aa_quad_evol, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Time_slice::n_evol, and Lorene::Time_slice::the_time.
Sets the factor jtime
) and deletes the value of N.
Definition at line 478 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
Sets the factor jtime
) and deletes the value of
Definition at line 453 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
Sets the conformal factor
Sets the value at the current time step (jtime
) and delete all quantities which depend on
Definition at line 515 of file isol_hor.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::psi_evol, and Lorene::Time_slice::the_time.
Sets the conformal factor
Sets the value at the current time step (jtime
) and deletes the value of N.
Definition at line 428 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
Sets the conformal factor
Sets the value at the current time step (jtime
) and deletes the value of
Definition at line 404 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Evolution< TyT >::downdate(), Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().
Sets the radius of the horizon to rad
.
Definition at line 431 of file isol_hor.h.
References radius.
Sets the order of the finite-differences scheme.
Definition at line 331 of file time_slice.h.
References Lorene::Time_slice::scheme_order.
Source for
Definition at line 190 of file sources_hor.C.
References Lorene::Time_slice_conf::aa(), Lorene::Tensor::annule_domain(), Lorene::Time_slice::beta(), Lorene::Metric::connect(), Lorene::contract(), Lorene::Scalar::derive_con(), Lorene::Scalar::derive_cov(), Lorene::Vector::derive_lie(), Lorene::Sym_tensor::divergence(), Lorene::Vector::divergence(), Lorene::Time_slice_conf::ff, gamt_point, Lorene::Connection::get_delta(), Lorene::Metric_flat::get_triad(), Lorene::Time_slice_conf::hdirac(), Lorene::Time_slice_conf::hh(), Lorene::Tensor::inc_dzpuis(), Lorene::Time_slice_conf::ln_psi(), met_gamt, mp, Lorene::Time_slice_conf::nn(), and trK.
Source for N
.
Definition at line 145 of file sources_hor.C.
References aa_quad(), Lorene::Tensor::annule_domain(), Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Scalar::derive_con(), Lorene::Scalar::derive_cov(), Lorene::Tensor::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Metric_flat::get_triad(), Lorene::Time_slice_conf::hdirac(), Lorene::Time_slice_conf::hh(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice_conf::ln_psi(), met_gamt, mp, Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi4(), trK, and trK_point.
Source for
Definition at line 108 of file sources_hor.C.
References aa_quad(), Lorene::Tensor::annule_domain(), Lorene::contract(), Lorene::Scalar::derive_cov(), Lorene::Tensor::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Metric_flat::get_triad(), Lorene::Time_slice_conf::hdirac(), Lorene::Time_slice_conf::hh(), Lorene::Scalar::inc_dzpuis(), met_gamt, mp, Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi4(), Lorene::Metric::ricci_scal(), and trK.
Conformal metric jtime
).
Reimplemented from Lorene::Time_slice_conf.
Definition at line 514 of file isol_hor.h.
References met_gamt.
Vector radial normal tilde.
Definition at line 116 of file phys_param.C.
References Lorene::Metric::con(), Lorene::Time_slice_conf::ff, Lorene::Metric::get_mp(), Lorene::Metric_flat::get_triad(), met_gamt, Lorene::Vector::set(), Lorene::sqrt(), and Lorene::Vector::std_spectral_base().
Trace K of the extrinsic curvature at the current time step (jtime
)
Reimplemented from Lorene::Time_slice.
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 796 of file time_slice_conf.C.
References Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Time_slice_conf::ff, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice_conf::ln_psi(), Lorene::Time_slice_conf::nn(), Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi_evol, Lorene::Time_slice::scheme_order, Lorene::Time_slice::the_time, Lorene::Evolution< TyT >::time_derive(), Lorene::Time_slice::trk_evol, and Lorene::Evolution_std< TyT >::update().
void Lorene::Isol_hor::update_aa | ( | ) |
Conformal representation
Definition at line 842 of file isol_hor.C.
References aa_quad_evol, Lorene::Time_slice::beta(), Lorene::contract(), gamt_point, Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Time_slice::jtime, met_gamt, mp, Lorene::Time_slice_conf::nn(), Lorene::Vector::ope_killing_conf(), Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi4(), regul, regularise_one(), Lorene::Tensor::set(), Lorene::Time_slice_conf::set_hata(), Lorene::Time_slice::the_time, and Lorene::Tensor::up_down().
Vector
(see the documentation of hata_evol
)
Definition at line 825 of file time_slice_conf.C.
References Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::hata_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, and Lorene::Time_slice_conf::p_vec_X.
double Lorene::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
WARNING this should only be used for an isolated black hole.
Vector
Definition at line 1251 of file bound_hor.C.
References Lorene::Scalar::annule_hard(), b_tilde(), Lorene::Time_slice::beta(), Lorene::Vector::change_triad(), Lorene::Metric::con(), Lorene::contract(), Lorene::Map::cosp, Lorene::Map::cost, Lorene::Scalar::dec_dzpuis(), Lorene::Tensor::derive_cov(), Lorene::Tensor::down(), Lorene::Time_slice_conf::ff, Lorene::Time_slice::gam(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Map::get_rot_phi(), Lorene::Scalar::inc_dzpuis(), kappa_hor(), met_gamt, mp, Lorene::Time_slice_conf::nn(), Lorene::Metric::radial_vect(), Lorene::Vector::set(), Lorene::Valeur::set_base(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_spectral_va(), Lorene::Map::sinp, Lorene::Map::sint, Lorene::Mg3d::std_base_vect_cart(), Lorene::Time_slice_conf::trk(), Lorene::Map::xa, and Lorene::Map::ya.
Vector
Definition at line 1379 of file bound_hor.C.
References Lorene::Map::get_bvect_cart(), and mp.
Definition at line 881 of file isol_hor.h.
|
mutableprotectedinherited |
Potential A associated with the symmetric tensor
(see documentation of Sym_tensor::p_aaa
).
Definition at line 547 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the components
Definition at line 310 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the components
Definition at line 316 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the components
Definition at line 320 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the components
Definition at line 323 of file isol_hor.h.
|
mutableprotectedinherited |
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 233 of file time_slice.h.
|
mutableprotectedinherited |
Potential
(see documentation of Sym_tensor::p_tilde_b
).
Definition at line 552 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the shift function
Definition at line 301 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the shift function
Definition at line 304 of file isol_hor.h.
|
mutableprotectedinherited |
Values at successive time steps of the shift vector
Definition at line 219 of file time_slice.h.
|
protected |
Boost velocity in x-direction.
Definition at line 272 of file isol_hor.h.
|
protected |
Boost velocity in z-direction.
Definition at line 275 of file isol_hor.h.
|
protected |
Function used to construct
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 345 of file isol_hor.h.
|
protectedinherited |
Number of stored time slices.
Definition at line 179 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the covariant derivative of the lapse with respect to the flat metric
Definition at line 294 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the covariant derivative of the conformal factor
Definition at line 298 of file isol_hor.h.
|
protectedinherited |
Pointer on the flat metric
Definition at line 507 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the covariant components of the induced metric
Definition at line 198 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the contravariant components of the induced metric
Definition at line 203 of file time_slice.h.
|
protected |
Time derivative of the 3-metric tilde.
Definition at line 329 of file isol_hor.h.
|
mutableprotectedinherited |
Values at successive time steps of the components
It is the conformal representation of the traceless part of the extrinsic curvature:
Definition at line 542 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the components
It is the deviation of the conformal metric
Definition at line 530 of file time_slice.h.
|
protectedinherited |
Time step index of the latest slice.
Definition at line 190 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the covariant components of the extrinsic curvature tensor
Definition at line 208 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the contravariant components of the extrinsic curvature tensor
Definition at line 213 of file time_slice.h.
|
protected |
3 metric tilde
Definition at line 326 of file isol_hor.h.
|
protected |
Affine mapping.
Definition at line 260 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the lapse function
Definition at line 281 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the lapse function
Definition at line 284 of file isol_hor.h.
|
mutableprotectedinherited |
Values at successive time steps of the lapse function N.
Definition at line 216 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the factor
Definition at line 522 of file time_slice.h.
|
protected |
Number of zones.
Definition at line 263 of file isol_hor.h.
|
protected |
Angular velocity in LORENE's units.
Definition at line 269 of file isol_hor.h.
|
mutableprotectedinherited |
Pointer on the induced metric at the current time step (jtime
)
Definition at line 239 of file time_slice.h.
|
mutableprotectedinherited |
Pointer on the vector jtime
).
Definition at line 571 of file time_slice.h.
|
mutableprotectedinherited |
Pointer on the logarithm of jtime
)
Definition at line 566 of file time_slice.h.
|
mutableprotectedinherited |
Pointer on the factor jtime
)
Definition at line 563 of file time_slice.h.
|
mutableprotectedinherited |
Pointer on the conformal metric jtime
)
Definition at line 560 of file time_slice.h.
|
mutableprotectedinherited |
Pointer on the vector
(see the documentation of hata_evol
)
Definition at line 577 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the conformal factor
Definition at line 287 of file isol_hor.h.
|
mutableprotected |
Values at successive time steps of the lapse function
Definition at line 290 of file isol_hor.h.
|
mutableprotectedinherited |
Values at successive time steps of the conformal factor
Definition at line 517 of file time_slice.h.
|
protected |
Radius of the horizon in LORENE's units.
Definition at line 266 of file isol_hor.h.
|
protected |
Intensity of the correction on the shift vector.
Definition at line 278 of file isol_hor.h.
|
protectedinherited |
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 187 of file time_slice.h.
|
protectedinherited |
Time label of each slice.
Definition at line 193 of file time_slice.h.
|
protected |
Trace of the extrinsic curvature.
Definition at line 332 of file isol_hor.h.
|
mutableprotectedinherited |
Values at successive time steps of the trace K of the extrinsic curvature.
Definition at line 224 of file time_slice.h.
|
protected |
Time derivative of the trace of the extrinsic curvature.
Definition at line 335 of file isol_hor.h.