Spacelike time slice of a 3+1 spacetime. More...
#include <time_slice.h>
Public Member Functions | |
| Time_slice (const Scalar &lapse_in, const Vector &shift_in, const Sym_tensor &gamma_in, const Sym_tensor &kk_in, int depth_in=3) | |
| General constructor (Hamiltonian-like). | |
| Time_slice (const Scalar &lapse_in, const Vector &shift_in, const Evolution_std< Sym_tensor > &gamma_in) | |
| General constructor (Lagrangian-like). | |
| Time_slice (const Map &mp, const Base_vect &triad, int depth_in=3) | |
| Constructor as standard time slice of flat spacetime (Minkowski). | |
| Time_slice (const Map &mp, const Base_vect &triad, FILE *fich, bool partial_read, int depth_in=3) | |
| Constructor from binary file. | |
| Time_slice (const Time_slice &) | |
| Copy constructor. | |
| virtual | ~Time_slice () |
| Destructor. | |
| void | operator= (const Time_slice &) |
Assignment to another Time_slice. | |
| 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 Scalar & | nn () const |
Lapse function N at the current time step (jtime ). | |
| virtual const Vector & | beta () const |
shift vector at the current time step (jtime ) | |
| const Metric & | gam () const |
Induced metric at the current time step (jtime ). | |
| virtual const Sym_tensor & | gam_dd () const |
Induced metric (covariant components ) at the current time step (jtime ). | |
| virtual const Sym_tensor & | gam_uu () const |
Induced metric (contravariant components ) at the current time step (jtime ). | |
| virtual const Sym_tensor & | k_dd () const |
Extrinsic curvature tensor (covariant components ) at the current time step (jtime ). | |
| virtual const Sym_tensor & | k_uu () const |
Extrinsic curvature tensor (contravariant components ) at the current time step (jtime ). | |
| virtual const Scalar & | trk () const |
Trace K of the extrinsic curvature 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 | |
| Time_slice (int depth_in) | |
| Special constructor for derived classes. | |
| virtual void | del_deriv () const |
| Deletes all the derived quantities. | |
| void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. | |
| virtual ostream & | operator>> (ostream &) const |
| Operator >> (virtual function called by the operator<<). | |
| virtual void | sauve (FILE *fich, bool partial_save) const |
| Total or partial saves in a binary file. | |
Protected Attributes | |
| 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 | |
| ostream & | operator<< (ostream &, const Time_slice &) |
| Display. | |
Spacelike time slice of a 3+1 spacetime.
()
Definition at line 169 of file time_slice.h.
| Time_slice::Time_slice | ( | const Scalar & | lapse_in, | |
| const Vector & | shift_in, | |||
| const Sym_tensor & | gamma_in, | |||
| const Sym_tensor & | kk_in, | |||
| int | depth_in = 3 | |||
| ) |
General constructor (Hamiltonian-like).
| lapse_in | lapse function N | |
| shift_in | shift vector | |
| gamma_in | induced metric (covariant or contravariant components) | |
| kk_in | extrinsic curvature (covariant or contravariant components) | |
| 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 108 of file time_slice.C.
References gam_dd_evol, gam_uu_evol, Tensor::get_index_type(), jtime, k_dd_evol, k_uu_evol, p_gamma, set_der_0x0(), the_time, Tensor::trace(), trk_evol, and Evolution_std< TyT >::update().
| Time_slice::Time_slice | ( | const Scalar & | lapse_in, | |
| const Vector & | shift_in, | |||
| const Evolution_std< Sym_tensor > & | gamma_in | |||
| ) |
General constructor (Lagrangian-like).
| lapse_in | lapse function N | |
| shift_in | shift vector | |
| gamma_in | induced metric (covariant or contravariant components) at various time steps; note that the scheme_order member is set to gamma_in.size - 1. scheme_order can be changed afterwards by the method set_scheme_order(int). |
Definition at line 151 of file time_slice.C.
References set_der_0x0().
Constructor as standard time slice of flat spacetime (Minkowski).
| mp | Mapping on which the various Lorene fields will be constructed | |
| triad | vector basis with respect to which the various tensor components will be defined | |
| 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 time_slice.C.
References beta_evol, Metric_flat::cov(), Map::flat_met_cart(), Map::flat_met_spher(), gam_dd_evol, jtime, k_dd_evol, n_evol, set_der_0x0(), Scalar::set_etat_one(), Scalar::set_etat_zero(), Tensor::set_etat_zero(), Scalar::std_spectral_base(), the_time, trk_evol, and Evolution_std< TyT >::update().
| Time_slice::Time_slice | ( | const Map & | mp, | |
| const Base_vect & | triad, | |||
| FILE * | fich, | |||
| bool | partial_read, | |||
| int | depth_in = 3 | |||
| ) |
Constructor from binary file.
The binary file must have been created by method save.
| mp | Mapping on which the various Lorene fields will be constructed | |
| triad | vector basis with respect to which the various tensor components will be defined | |
| fich | file containing the saved Time_slice | |
| 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; the given value must coincide with that stored in the file. |
Definition at line 231 of file time_slice.C.
References beta_evol, depth, fread_be(), Map::get_mg(), jtime, n_evol, scheme_order, set_der_0x0(), the_time, and Evolution_std< TyT >::update().
| Time_slice::Time_slice | ( | const Time_slice & | tin | ) |
| Time_slice::Time_slice | ( | int | depth_in | ) | [explicit, protected] |
Special constructor for derived classes.
Definition at line 329 of file time_slice.C.
References set_der_0x0().
| Time_slice::~Time_slice | ( | ) | [virtual] |
| double Time_slice::adm_mass | ( | ) | const [virtual] |
Returns the ADM mass (geometrical units) at the current step.
Moreover this method updates adm_mass_evol if necessary.
Reimplemented in Tslice_dirac_max.
Definition at line 70 of file tslice_adm_mass.C.
References adm_mass_evol, Tensor::derive_con(), Tensor_sym::derive_con(), Vector::flux(), gam_dd(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Tbl::get_taille(), jtime, Tbl::set(), Tbl::set_etat_qcq(), the_time, Tensor::trace(), Tensor::up(), Evolution_full< TyT >::update(), and Map::val_r().
| const Vector & Time_slice::beta | ( | ) | const [virtual] |
shift vector
at the current time step (jtime )
Definition at line 83 of file time_slice_access.C.
References beta_evol, Evolution< TyT >::is_known(), and jtime.
| Tbl Time_slice::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.
| strain_tensor | : a pointer on the strain_tensor measured by the Eulerian observer of 4-velocity ; if this is the null pointer, it is assumed that = 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 |
= 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(), beta(), Metric::con(), contract(), Tensor::derive_con(), Scalar::derive_con(), Sym_tensor::derive_lie(), gam(), Tbl::get_dim(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), jtime, k_dd(), k_dd_evol, k_uu(), maxabs(), nn(), Tensor_sym::position(), Metric::ricci(), scheme_order, Itbl::set(), Evolution< TyT >::time_derive(), Tensor::trace(), trk(), and Tensor::up_down().
| Tbl Time_slice::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.
| energy_density | : a pointer on the energy density E measured by the Eulerian observer of 4-velocity ; if this is the null pointer, it is assumed that E = 0 (vacuum). | |
| ost | : output stream for a formatted output of the result |
Definition at line 75 of file tslice_check_einstein.C.
References contract(), Scalar::dec_dzpuis(), gam(), k_dd(), k_uu(), maxabs(), Metric::ricci_scal(), and trk().
| Tbl Time_slice::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.
| momentum_density | : a pointer on the momentum density measured by the Eulerian observer of 4-velocity ; if this is the null pointer, it is assumed that = 0 (vacuum). | |
| ost | : output stream for a formatted output of the result |
= 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(), gam(), Tensor::get_index_type(), k_uu(), maxabs(), and trk().
| void Time_slice::del_deriv | ( | ) | const [protected, virtual] |
Deletes all the derived quantities.
Reimplemented in Time_slice_conf.
Definition at line 363 of file time_slice.C.
References adm_mass_evol, Evolution< TyT >::downdate(), jtime, p_gamma, and set_der_0x0().
| const Metric & Time_slice::gam | ( | ) | const |
Induced metric
at the current time step (jtime ).
Definition at line 91 of file time_slice_access.C.
| const Sym_tensor & Time_slice::gam_dd | ( | ) | const [virtual] |
Induced metric (covariant components
) at the current time step (jtime ).
Reimplemented in Time_slice_conf.
Definition at line 103 of file time_slice_access.C.
References Metric::cov(), gam_dd_evol, gam_uu_evol, Evolution< TyT >::is_known(), jtime, p_gamma, the_time, and Evolution_std< TyT >::update().
| const Sym_tensor & Time_slice::gam_uu | ( | ) | const [virtual] |
Induced metric (contravariant components
) at the current time step (jtime ).
Reimplemented in Time_slice_conf.
Definition at line 118 of file time_slice_access.C.
References gam(), gam_dd_evol, gam_uu_evol, Evolution< TyT >::is_known(), jtime, the_time, and Evolution_std< TyT >::update().
| int Time_slice::get_latest_j | ( | ) | const [inline] |
Gets the latest value of time step index.
Definition at line 339 of file time_slice.h.
References jtime.
| int Time_slice::get_scheme_order | ( | ) | const [inline] |
Gets the order of the finite-differences scheme.
Definition at line 336 of file time_slice.h.
References scheme_order.
| const Evolution_std<double>& Time_slice::get_time | ( | ) | const [inline] |
Gets the time coordinate t at successive time steps.
Definition at line 342 of file time_slice.h.
References the_time.
| const Sym_tensor & Time_slice::k_dd | ( | ) | const [virtual] |
Extrinsic curvature tensor (covariant components
) at the current time step (jtime ).
Reimplemented in Time_slice_conf.
Definition at line 131 of file time_slice_access.C.
References beta(), Tensor::down(), gam(), gam_dd(), gam_dd_evol, Evolution< TyT >::is_known(), jtime, k_dd_evol, nn(), Vector::ope_killing(), scheme_order, the_time, Evolution< TyT >::time_derive(), and Evolution_std< TyT >::update().
| const Sym_tensor & Time_slice::k_uu | ( | ) | const [virtual] |
Extrinsic curvature tensor (contravariant components
) at the current time step (jtime ).
Reimplemented in Time_slice_conf.
Definition at line 153 of file time_slice_access.C.
References beta(), gam(), gam_uu(), gam_uu_evol, Evolution< TyT >::is_known(), jtime, k_uu_evol, nn(), Vector::ope_killing(), scheme_order, the_time, Evolution< TyT >::time_derive(), and Evolution_std< TyT >::update().
| const Scalar & Time_slice::nn | ( | ) | const [virtual] |
Lapse function N at the current time step (jtime ).
Reimplemented in Time_slice_conf.
Definition at line 75 of file time_slice_access.C.
References Evolution< TyT >::is_known(), jtime, and n_evol.
| void Time_slice::operator= | ( | const Time_slice & | tin | ) |
Assignment to another Time_slice.
Reimplemented in Isol_hor, Time_slice_conf, and Tslice_dirac_max.
Definition at line 384 of file time_slice.C.
References beta_evol, del_deriv(), depth, gam_dd_evol, gam_uu_evol, jtime, k_dd_evol, k_uu_evol, n_evol, scheme_order, the_time, and trk_evol.
| ostream & Time_slice::operator>> | ( | ostream & | flux | ) | const [protected, virtual] |
Operator >> (virtual function called by the operator<<).
Reimplemented in Isol_hor, and Tslice_dirac_max.
Definition at line 407 of file time_slice.C.
References adm_mass(), adm_mass_evol, beta_evol, depth, gam_dd_evol, gam_uu_evol, Evolution< TyT >::is_known(), jtime, k_dd_evol, k_uu_evol, maxabs(), n_evol, p_gamma, scheme_order, the_time, and trk_evol.
| void Time_slice::sauve | ( | FILE * | fich, | |
| bool | partial_save | |||
| ) | const [protected, virtual] |
Total or partial saves in a binary file.
This protected method is to be called either from public method save or from method sauve of a derived class.
| fich | binary file | |
| partial_save | indicates whether the whole object must be saved. |
Reimplemented in Isol_hor, and Tslice_dirac_max.
Definition at line 503 of file time_slice.C.
References beta(), beta_evol, depth, fwrite_be(), Evolution< TyT >::is_known(), jtime, n_evol, nn(), scheme_order, and the_time.
| void Time_slice::save | ( | const char * | rootname | ) | const |
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 457 of file time_slice.C.
References beta(), depth, fwrite_be(), Map::get_mg(), Tensor::get_mp(), Tensor::get_triad(), jtime, nn(), sauve(), Base_vect::sauve(), Map::sauve(), and Mg3d::sauve().
| void Time_slice::set_der_0x0 | ( | ) | const [protected] |
Sets to 0x0 all the pointers on derived quantities.
Reimplemented in Time_slice_conf.
Definition at line 373 of file time_slice.C.
References p_gamma.
| void Time_slice::set_scheme_order | ( | int | ord | ) | [inline] |
Sets the order of the finite-differences scheme.
Definition at line 327 of file time_slice.h.
References scheme_order.
| const Scalar & Time_slice::trk | ( | ) | const [virtual] |
Trace K of the extrinsic curvature at the current time step (jtime ).
Reimplemented in Time_slice_conf, and Tslice_dirac_max.
Definition at line 173 of file time_slice_access.C.
References gam(), Evolution< TyT >::is_known(), jtime, k_dd(), k_uu(), k_uu_evol, the_time, trk_evol, and Evolution_std< TyT >::update().
| ostream& operator<< | ( | ostream & | , | |
| const Time_slice & | ||||
| ) | [friend] |
Display.
Evolution_full<Tbl> Time_slice::adm_mass_evol [mutable, protected] |
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<Vector> Time_slice::beta_evol [mutable, protected] |
Values at successive time steps of the shift vector
.
Definition at line 215 of file time_slice.h.
int Time_slice::depth [protected] |
Number of stored time slices.
Definition at line 175 of file time_slice.h.
Evolution_std<Sym_tensor> Time_slice::gam_dd_evol [mutable, protected] |
Values at successive time steps of the covariant components of the induced metric
.
Definition at line 194 of file time_slice.h.
Evolution_std<Sym_tensor> Time_slice::gam_uu_evol [mutable, protected] |
Values at successive time steps of the contravariant components of the induced metric
.
Definition at line 199 of file time_slice.h.
int Time_slice::jtime [protected] |
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] |
Values at successive time steps of the covariant components of the extrinsic curvature tensor
.
Definition at line 204 of file time_slice.h.
Evolution_std<Sym_tensor> Time_slice::k_uu_evol [mutable, protected] |
Values at successive time steps of the contravariant components of the extrinsic curvature tensor
.
Definition at line 209 of file time_slice.h.
Evolution_std<Scalar> Time_slice::n_evol [mutable, protected] |
Values at successive time steps of the lapse function N.
Definition at line 212 of file time_slice.h.
Metric* Time_slice::p_gamma [mutable, protected] |
Pointer on the induced metric at the current time step (jtime).
Definition at line 235 of file time_slice.h.
int Time_slice::scheme_order [protected] |
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] |
Time label of each slice.
Definition at line 189 of file time_slice.h.
Evolution_std<Scalar> Time_slice::trk_evol [mutable, protected] |
Values at successive time steps of the trace K of the extrinsic curvature.
Definition at line 220 of file time_slice.h.
1.6.1