LORENE
|
Spacelike time slice of a 3+1 spacetime with conformal decomposition. More...
#include <time_slice.h>
Public Member Functions | |
Time_slice_conf (const Scalar &lapse_in, const Vector &shift_in, const Metric_flat &ff_in, const Scalar &psi_in, const Sym_tensor &hh_in, const Sym_tensor &hata_in, const Scalar &trk_in, int depth_in=3) | |
Constructor from conformal decomposition. More... | |
Time_slice_conf (const Scalar &lapse_in, const Vector &shift_in, const Sym_tensor &gamma_in, const Sym_tensor &kk_in, const Metric_flat &ff_in, int depth_in=3) | |
Constructor from physical metric. More... | |
Time_slice_conf (const Map &mp, const Base_vect &triad, const Metric_flat &ff_in, int depth_in=3) | |
Constructor as standard time slice of flat spacetime (Minkowski). More... | |
Time_slice_conf (const Map &mp, const Base_vect &triad, const Metric_flat &ff_in, FILE *fich, bool partial_read, int depth_in=3) | |
Constructor from binary file. More... | |
Time_slice_conf (const Time_slice_conf &) | |
Copy constructor. More... | |
virtual | ~Time_slice_conf () |
Destructor. More... | |
void | operator= (const Time_slice_conf &) |
Assignment to another Time_slice_conf . More... | |
void | operator= (const Time_slice &) |
Assignment to a Time_slice . More... | |
virtual void | set_psi_del_npsi (const Scalar &psi_in) |
Sets the conformal factor relating the physical metric to the conformal one: . More... | |
virtual void | set_psi_del_n (const Scalar &psi_in) |
Sets the conformal factor relating the physical metric to the conformal one: . More... | |
virtual void | set_npsi_del_psi (const Scalar &npsi_in) |
Sets the factor at the current time step (jtime ) and deletes the value of . More... | |
virtual void | set_npsi_del_n (const Scalar &npsi_in) |
Sets the factor at the current time step (jtime ) and deletes the value of N. More... | |
virtual void | set_hh (const Sym_tensor &hh_in) |
Sets the deviation of the conformal metric from the flat metric : . More... | |
virtual void | set_hata (const Sym_tensor &hata_in) |
Sets the conformal representation of the traceless part of the extrinsic curvature: . More... | |
virtual void | set_hata_TT (const Sym_tensor_tt &hata_tt) |
Sets the TT part of (see member hata_evol ). More... | |
virtual void | set_hata_from_XAB (Param *par_bc=0x0, Param *par_mat=0x0) |
Sets the conformal representation of the traceless part of the extrinsic curvature from its potentials A, and . More... | |
virtual const Scalar & | nn () const |
Lapse function N at the current time step (jtime ) More... | |
virtual const Sym_tensor & | gam_dd () const |
Induced metric (covariant components ) at the current time step (jtime ) More... | |
virtual const Sym_tensor & | gam_uu () const |
Induced metric (contravariant components ) at the current time step (jtime ) More... | |
virtual const Sym_tensor & | k_dd () const |
Extrinsic curvature tensor (covariant components ) at the current time step (jtime ) More... | |
virtual const Sym_tensor & | k_uu () const |
Extrinsic curvature tensor (contravariant components ) at the current time step (jtime ) More... | |
virtual const Scalar & | A_hata () const |
Returns the potential A of . More... | |
virtual const Scalar & | B_hata () const |
Returns the potential of . More... | |
virtual const Scalar & | psi () const |
Conformal factor relating the physical metric to the conformal one: . More... | |
const Scalar & | psi4 () const |
Factor at the current time step (jtime ). More... | |
const Scalar & | ln_psi () const |
Logarithm of at the current time step (jtime ). More... | |
virtual const Scalar & | npsi () const |
Factor at the current time step (jtime ). More... | |
virtual const Metric & | tgam () const |
Conformal metric Returns the value at the current time step (jtime ). More... | |
virtual const Sym_tensor & | hh (Param *=0x0, Param *=0x0) const |
Deviation of the conformal metric from the flat metric : . More... | |
virtual const Sym_tensor & | hata () const |
Conformal representation of the traceless part of the extrinsic curvature: . More... | |
virtual Sym_tensor | aa () const |
Conformal representation of the traceless part of the extrinsic curvature: . More... | |
virtual const Scalar & | trk () const |
Trace K of the extrinsic curvature at the current time step (jtime ) More... | |
virtual const Vector & | hdirac () const |
Vector which vanishes in Dirac gauge. More... | |
virtual const Vector & | vec_X (int method_poisson=6) const |
Vector representing the longitudinal part of . More... | |
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 from the conformally-rescaled momentum , using the momentum constraint. More... | |
virtual void | set_AB_hata (const Scalar &A_in, const Scalar &B_in) |
Sets the potentials A and of the TT part (see the documentation of Sym_tensor for details). More... | |
virtual void | initial_data_cts (const Sym_tensor &uu, const Scalar &trk_in, const Scalar &trk_point, double pdt, double precis=1.e-12, int method_poisson_vect=6, const char *graph_device=0x0, const Scalar *ener_dens=0x0, const Vector *mom_dens=0x0, const Scalar *trace_stress=0x0) |
Computes valid initial data by solving the constraint equations in the conformal thin-sandwich approach. More... | |
virtual double | adm_mass () const |
Returns the ADM mass (geometrical units) at the current step. More... | |
void | check_psi_dot (Tbl &tlnpsi_dot, Tbl &tdiff, Tbl &tdiff_rel) const |
Checks the relation. More... | |
void | set_scheme_order (int ord) |
Sets the order of the finite-differences scheme. More... | |
int | get_scheme_order () const |
Gets the order of the finite-differences scheme. More... | |
int | get_latest_j () const |
Gets the latest value of time step index. More... | |
const Evolution_std< double > & | get_time () const |
Gets the time coordinate t at successive time steps. More... | |
virtual const Vector & | beta () const |
shift vector at the current time step (jtime ) More... | |
const Metric & | gam () const |
Induced metric at the current time step (jtime ) More... | |
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. More... | |
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. More... | |
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. More... | |
void | save (const char *rootname) const |
Saves in a binary file. More... | |
Protected Member Functions | |
virtual void | del_deriv () const |
Deletes all the derived quantities. More... | |
void | set_der_0x0 () const |
Sets to 0x0 all the pointers on derived quantities. More... | |
virtual ostream & | operator>> (ostream &) const |
Operator >> (virtual function called by the operator<<). More... | |
virtual void | sauve (FILE *fich, bool partial_save) const |
Total or partial saves in a binary file. More... | |
Protected Attributes | |
const Metric_flat & | ff |
Pointer on the flat metric with respect to which the conformal decomposition is performed. More... | |
Evolution_std< Scalar > | psi_evol |
Values at successive time steps of the conformal factor relating the physical metric to the conformal one: . More... | |
Evolution_std< Scalar > | npsi_evol |
Values at successive time steps of the factor . More... | |
Evolution_std< Sym_tensor > | hh_evol |
Values at successive time steps of the components . More... | |
Evolution_std< Sym_tensor > | hata_evol |
Values at successive time steps of the components . More... | |
Evolution_std< Scalar > | A_hata_evol |
Potential A associated with the symmetric tensor . More... | |
Evolution_std< Scalar > | B_hata_evol |
Potential associated with the symmetric tensor . More... | |
Metric * | p_tgamma |
Pointer on the conformal metric at the current time step (jtime ) More... | |
Scalar * | p_psi4 |
Pointer on the factor at the current time step (jtime ) More... | |
Scalar * | p_ln_psi |
Pointer on the logarithm of at the current time step (jtime ) More... | |
Vector * | p_hdirac |
Pointer on the vector (which vanishes in Dirac gauge), at the current time step (jtime ). More... | |
Vector * | p_vec_X |
Pointer on the vector representing the longitudinal part of . More... | |
int | depth |
Number of stored time slices. More... | |
int | scheme_order |
Order of the finite-differences scheme for the computation of time derivatives. More... | |
int | jtime |
Time step index of the latest slice. More... | |
Evolution_std< double > | the_time |
Time label of each slice. More... | |
Evolution_std< Sym_tensor > | gam_dd_evol |
Values at successive time steps of the covariant components of the induced metric . More... | |
Evolution_std< Sym_tensor > | gam_uu_evol |
Values at successive time steps of the contravariant components of the induced metric . More... | |
Evolution_std< Sym_tensor > | k_dd_evol |
Values at successive time steps of the covariant components of the extrinsic curvature tensor . More... | |
Evolution_std< Sym_tensor > | k_uu_evol |
Values at successive time steps of the contravariant components of the extrinsic curvature tensor . More... | |
Evolution_std< Scalar > | n_evol |
Values at successive time steps of the lapse function N. More... | |
Evolution_std< Vector > | beta_evol |
Values at successive time steps of the shift vector . More... | |
Evolution_std< Scalar > | trk_evol |
Values at successive time steps of the trace K of the extrinsic curvature. More... | |
Evolution_full< Tbl > | adm_mass_evol |
ADM mass at each time step, since the creation of the slice. More... | |
Metric * | p_gamma |
Pointer on the induced metric at the current time step (jtime ) More... | |
Spacelike time slice of a 3+1 spacetime with conformal decomposition.
()
Definition at line 501 of file time_slice.h.
Lorene::Time_slice_conf::Time_slice_conf | ( | const Scalar & | lapse_in, |
const Vector & | shift_in, | ||
const Metric_flat & | ff_in, | ||
const Scalar & | psi_in, | ||
const Sym_tensor & | hh_in, | ||
const Sym_tensor & | hata_in, | ||
const Scalar & | trk_in, | ||
int | depth_in = 3 |
||
) |
Constructor from conformal decomposition.
lapse_in | lapse function N |
shift_in | shift vector |
ff_in | reference flat metric with respect to which the conformal decomposition is performed |
psi_in | conformal factor relating the physical metric to the conformal one: |
hh_in | deviation of the conformal metric from the flat metric : . must be such that . |
hata_in | conformal representation of the traceless part of the extrinsic curvature: |
trk_in | trace K of the extrinsic curvature |
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 121 of file time_slice_conf.C.
References A_hata(), B_hata(), Lorene::Time_slice::beta_evol, Lorene::Metric_flat::con(), Lorene::Metric_flat::determinant(), ff, Lorene::Tensor::get_index_type(), Lorene::Time_slice::jtime, Lorene::max(), Lorene::maxabs(), Lorene::Time_slice::n_evol, set_der_0x0(), Lorene::Time_slice::the_time, and Lorene::Time_slice::trk_evol.
Lorene::Time_slice_conf::Time_slice_conf | ( | const Scalar & | lapse_in, |
const Vector & | shift_in, | ||
const Sym_tensor & | gamma_in, | ||
const Sym_tensor & | kk_in, | ||
const Metric_flat & | ff_in, | ||
int | depth_in = 3 |
||
) |
Constructor from physical metric.
The conformal decomposition is performed by the constructor.
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) |
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 175 of file time_slice_conf.C.
References A_hata(), B_hata(), Lorene::Metric::con(), Lorene::Metric_flat::con(), Lorene::Metric::determinant(), Lorene::Metric_flat::determinant(), ff, hata_evol, hh_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_uu(), Lorene::Time_slice::p_gamma, p_psi4, Lorene::pow(), psi_evol, set_der_0x0(), Lorene::Scalar::std_spectral_base(), Lorene::Time_slice::the_time, and Lorene::Time_slice::trk_evol.
Lorene::Time_slice_conf::Time_slice_conf | ( | const Map & | mp, |
const Base_vect & | triad, | ||
const Metric_flat & | ff_in, | ||
int | depth_in = 3 |
||
) |
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 |
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 212 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, hata_evol, hh_evol, Lorene::Time_slice::jtime, npsi_evol, psi_evol, set_der_0x0(), Lorene::Scalar::set_etat_one(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::std_spectral_base(), and Lorene::Time_slice::the_time.
Lorene::Time_slice_conf::Time_slice_conf | ( | const Map & | mp, |
const Base_vect & | triad, | ||
const Metric_flat & | ff_in, | ||
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 |
ff_in | reference flat metric with respect to which the conformal decomposition is performed |
fich | file containing the saved Time_slice_conf |
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 must coincide with that stored in the file. |
Definition at line 253 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, Lorene::Time_slice::depth, Lorene::fread_be(), Lorene::Map::get_mg(), hata_evol, Lorene::Time_slice::jtime, psi_evol, set_der_0x0(), and Lorene::Time_slice::the_time.
Lorene::Time_slice_conf::Time_slice_conf | ( | const Time_slice_conf & | tin | ) |
|
virtual |
|
virtual |
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 667 of file time_slice_conf.C.
References A_hata_evol, hata_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice::the_time.
|
virtual |
Conformal representation of the traceless part of the extrinsic curvature: .
Returns the value at the current time step (jtime
).
Definition at line 768 of file time_slice_conf.C.
|
virtual |
Returns the ADM mass (geometrical units) at the current step.
Moreover this method updates adm_mass_evol
if necessary.
Reimplemented from Lorene::Time_slice.
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 111 of file tslice_adm_mass.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Scalar::derive_con(), ff, Lorene::Vector::flux(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Tbl::get_taille(), hdirac(), hh(), Lorene::Time_slice::jtime, psi(), Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), Lorene::Time_slice::the_time, Lorene::Tensor::trace(), and Lorene::Map::val_r().
|
virtual |
Returns the potential of .
See the documentation of Sym_tensor_tt
for details. Returns the value at the current time step (jtime
).
Definition at line 681 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, hata_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice::the_time.
|
virtualinherited |
shift vector at the current time step (jtime
)
Definition at line 90 of file time_slice_access.C.
References Lorene::Time_slice::beta_evol, and Lorene::Time_slice::jtime.
|
inherited |
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 |
Definition at line 142 of file tslice_check_einstein.C.
References Lorene::Tensor::annule_domain(), Lorene::Tbl::annule_hard(), Lorene::Time_slice::beta(), Lorene::Metric::con(), Lorene::contract(), Lorene::Tensor::derive_con(), Lorene::Scalar::derive_con(), Lorene::Sym_tensor::derive_lie(), Lorene::Time_slice::gam(), Lorene::Tbl::get_dim(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), 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::Tensor_sym::position(), Lorene::Metric::ricci(), Lorene::Time_slice::scheme_order, Lorene::Itbl::set(), Lorene::Tensor::trace(), Lorene::Time_slice::trk(), and Lorene::Tensor::up_down().
|
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 ; 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 82 of file tslice_check_einstein.C.
References Lorene::contract(), Lorene::Scalar::dec_dzpuis(), Lorene::Time_slice::gam(), Lorene::Time_slice::k_dd(), Lorene::Time_slice::k_uu(), Lorene::maxabs(), Lorene::Metric::ricci_scal(), and Lorene::Time_slice::trk().
|
inherited |
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 |
Definition at line 112 of file tslice_check_einstein.C.
References Lorene::Scalar::derive_con(), Lorene::Sym_tensor::divergence(), Lorene::Time_slice::gam(), Lorene::Tensor::get_index_type(), Lorene::Time_slice::k_uu(), Lorene::maxabs(), and Lorene::Time_slice::trk().
void Lorene::Time_slice_conf::check_psi_dot | ( | Tbl & | tlnpsi_dot, |
Tbl & | tdiff, | ||
Tbl & | tdiff_rel | ||
) | const |
Checks the relation.
tlnpsi_dot | [output] maximun value in each domain of |
tdiff | [output] maximum value in each domain of |
tdiff_rel | [output] relative error on the above relation in each domain. |
Definition at line 907 of file time_slice_conf.C.
References Lorene::abs(), Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Scalar::dec_dzpuis(), Lorene::Vector::divergence(), ff, Lorene::Time_slice::jtime, ln_psi(), Lorene::max(), Lorene::Time_slice::n_evol, nn(), npsi(), npsi_evol, psi(), psi_evol, Lorene::Time_slice::scheme_order, and trk().
void Lorene::Time_slice_conf::compute_X_from_momentum_constraint | ( | const Vector & | hat_S, |
const Sym_tensor_tt & | hata_tt, | ||
int | iter_max = 200 , |
||
double | precis = 1.e-12 , |
||
double | relax = 0.8 , |
||
int | methode_poisson = 6 |
||
) |
Computes the vector from the conformally-rescaled momentum , using the momentum constraint.
Definition at line 840 of file time_slice_conf.C.
References Lorene::abs(), Lorene::contract(), Lorene::Tensor::get_index_type(), Lorene::max(), and Lorene::Vector::poisson().
|
protectedvirtual |
Deletes all the derived quantities.
Reimplemented from Lorene::Time_slice.
Definition at line 350 of file time_slice_conf.C.
References Lorene::Time_slice::del_deriv(), p_hdirac, p_ln_psi, p_psi4, p_tgamma, p_vec_X, and set_der_0x0().
|
inherited |
Induced metric at the current time step (jtime
)
Definition at line 98 of file time_slice_access.C.
References Lorene::Time_slice::gam_dd(), and Lorene::Time_slice::p_gamma.
|
virtual |
Induced metric (covariant components ) at the current time step (jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 611 of file time_slice_conf.C.
References Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::jtime, psi4(), tgam(), and Lorene::Time_slice::the_time.
|
virtual |
Induced metric (contravariant components ) at the current time step (jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 622 of file time_slice_conf.C.
References Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, psi4(), tgam(), and Lorene::Time_slice::the_time.
|
inlineinherited |
Gets the latest value of time step index.
Definition at line 346 of file time_slice.h.
References Lorene::Time_slice::jtime.
|
inlineinherited |
Gets the order of the finite-differences scheme.
Definition at line 343 of file time_slice.h.
References Lorene::Time_slice::scheme_order.
|
inlineinherited |
Gets the time coordinate t at successive time steps.
Definition at line 349 of file time_slice.h.
References Lorene::Time_slice::the_time.
|
virtual |
Conformal representation of the traceless part of the extrinsic curvature: .
Returns the value at the current time step (jtime
).
Definition at line 775 of file time_slice_conf.C.
References Lorene::Time_slice::beta(), Lorene::Sym_tensor::derive_lie(), Lorene::Vector::divergence(), ff, hata_evol, hh(), hh_evol, Lorene::Time_slice::jtime, nn(), Lorene::Vector::ope_killing_conf(), psi(), psi4(), Lorene::Time_slice::scheme_order, and Lorene::Time_slice::the_time.
|
virtual |
Vector which vanishes in Dirac gauge.
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 818 of file time_slice_conf.C.
|
virtual |
Deviation of the conformal metric from the flat metric : .
Returns the value at the current time step (jtime
).
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 761 of file time_slice_conf.C.
References hh_evol, and Lorene::Time_slice::jtime.
|
virtual |
Computes valid initial data by solving the constraint equations in the conformal thin-sandwich approach.
uu | value of (freely specifiable data). This quantity must be trace-free with respect to the conformal metric , reflecting the unimodular character of . |
trk_in | value of (freely specifiable data) |
trk_point | value of (freely specifiable data) |
pdt | time step, to be used in order to fill depth slices |
precis | convergence threshold required to stop the iteration |
method_poisson_vect | method to be used for solving vector Poisson equation (for the shift), see Vector::poisson(double, const Metric_flat&, int) const . |
graph_device | name of type of graphical device: 0x0 (default value) will result in interactive choice; "/xwin" in X-Window display and "/n" in no output. |
ener_dens | matter energy density E as measured by the Eulerian observer; this quantity is passed as a pointer, the null value of which (default) meaning E=0. |
mom_dens | matter momentum density J as measured by the Eulerian observer; this quantity is passed as a pointer, the null value of which (default) meaning J=0. |
trace_stress | trace of the matter stress S as measured by the Eulerian observer; this quantity is passed as a pointer, the null value of which (default) meaning S=0. |
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 94 of file tslice_conf_init.C.
References aa(), Lorene::Time_slice::beta(), Lorene::Scalar::check_dzpuis(), Lorene::contract(), Lorene::Scalar::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Vector::derive_lie(), Lorene::Sym_tensor::divergence(), Lorene::Vector::divergence(), ff, Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::get_triad(), hdirac(), hh(), Lorene::Tensor::inc_dzpuis(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, ln_psi(), Lorene::max(), Lorene::maxabs(), nn(), Lorene::Scalar::poisson(), psi(), psi4(), Lorene::Metric::ricci_scal(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), set_hata(), tgam(), Lorene::Time_slice::the_time, Lorene::Tensor::trace(), trk(), Lorene::Time_slice::trk_evol, Lorene::Tensor::up_down(), and Lorene::Map::val_r().
|
virtual |
Extrinsic curvature tensor (covariant components ) at the current time step (jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 633 of file time_slice_conf.C.
References Lorene::Time_slice::gam(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, k_uu(), and Lorene::Time_slice::the_time.
|
virtual |
Extrinsic curvature tensor (contravariant components ) at the current time step (jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 646 of file time_slice_conf.C.
References Lorene::Time_slice::gam(), hata(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_uu_evol, psi(), psi4(), Lorene::Time_slice::the_time, and trk().
const Scalar & Lorene::Time_slice_conf::ln_psi | ( | ) | const |
Logarithm of at the current time step (jtime
).
Definition at line 722 of file time_slice_conf.C.
References Lorene::log(), p_ln_psi, psi(), and Lorene::Scalar::std_spectral_base().
|
virtual |
Lapse function N at the current time step (jtime
)
Reimplemented from Lorene::Time_slice.
Definition at line 594 of file time_slice_conf.C.
References Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, npsi_evol, psi_evol, and Lorene::Time_slice::the_time.
|
virtual |
Factor at the current time step (jtime
).
Definition at line 735 of file time_slice_conf.C.
References Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, npsi_evol, psi_evol, and Lorene::Time_slice::the_time.
void Lorene::Time_slice_conf::operator= | ( | const Time_slice_conf & | tin | ) |
Assignment to another Time_slice_conf
.
Definition at line 379 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, del_deriv(), hata_evol, hh_evol, npsi_evol, Lorene::Time_slice::operator=(), and psi_evol.
void Lorene::Time_slice_conf::operator= | ( | const Time_slice & | tin | ) |
Assignment to a Time_slice
.
Definition at line 394 of file time_slice_conf.C.
References del_deriv(), and Lorene::Time_slice::operator=().
|
protectedvirtual |
Operator >> (virtual function called by the operator<<).
Reimplemented from Lorene::Time_slice.
Reimplemented in Lorene::Tslice_dirac_max, and Lorene::Isol_hor.
Definition at line 944 of file time_slice_conf.C.
References ff, Lorene::Metric_flat::get_triad(), hata_evol, hh_evol, Lorene::Time_slice::jtime, Lorene::maxabs(), npsi_evol, Lorene::Time_slice::operator>>(), p_hdirac, p_ln_psi, p_psi4, p_tgamma, p_vec_X, and psi_evol.
|
virtual |
Conformal factor relating the physical metric to the conformal one: .
is defined by
Returns the value at the current time step (jtime
).
Definition at line 696 of file time_slice_conf.C.
References Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, npsi_evol, psi_evol, and Lorene::Time_slice::the_time.
const Scalar & Lorene::Time_slice_conf::psi4 | ( | ) | const |
Factor at the current time step (jtime
).
Definition at line 710 of file time_slice_conf.C.
References p_psi4, Lorene::pow(), psi(), and Lorene::Scalar::std_spectral_base().
|
protectedvirtual |
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 from Lorene::Time_slice.
Reimplemented in Lorene::Tslice_dirac_max, and Lorene::Isol_hor.
Definition at line 977 of file time_slice_conf.C.
References A_hata(), A_hata_evol, B_hata(), B_hata_evol, Lorene::Time_slice::depth, Lorene::fwrite_be(), hata(), hata_evol, Lorene::Time_slice::jtime, psi(), psi_evol, and Lorene::Time_slice::sauve().
|
inherited |
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 464 of file time_slice.C.
References Lorene::Time_slice::beta(), Lorene::Time_slice::depth, Lorene::fwrite_be(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Tensor::get_triad(), Lorene::Time_slice::jtime, Lorene::Time_slice::nn(), Lorene::Base_vect::sauve(), Lorene::Time_slice::sauve(), Lorene::Mg3d::sauve(), and Lorene::Map::sauve().
Sets the potentials A and of the TT part (see the documentation of Sym_tensor
for details).
Sets the value at the current time step (jtime
).
Definition at line 894 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, and Lorene::Time_slice::the_time.
|
protected |
|
virtual |
Sets the conformal representation of the traceless part of the extrinsic curvature: .
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 537 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, and Lorene::Time_slice::the_time.
|
virtual |
Sets the conformal representation of the traceless part of the extrinsic curvature from its potentials A, and .
These potentials must be up-to-date. It sets the value at the current time step (jtime
).
Definition at line 569 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, ff, Lorene::Tensor::get_mp(), hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Vector::ope_killing_conf(), p_vec_X, Lorene::Sym_tensor_tt::set_A_tildeB(), and Lorene::Time_slice::the_time.
|
virtual |
Sets the TT part of (see member hata_evol
).
Sets the value at current time-step (jtime
) and updates the potentials A and .
Definition at line 549 of file time_slice_conf.C.
References A_hata_evol, B_hata_evol, Lorene::Sym_tensor::compute_A(), Lorene::Sym_tensor::compute_tilde_B_tt(), Lorene::Scalar::dec_dzpuis(), Lorene::Scalar::get_dzpuis(), hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, and Lorene::Time_slice::the_time.
|
virtual |
Sets the deviation of the conformal metric from the flat metric : .
must be such that . Sets the value at the current time step (jtime
).
Reimplemented in Lorene::Tslice_dirac_max.
Definition at line 492 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Metric_flat::con(), Lorene::Metric_flat::determinant(), ff, Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, hh_evol, Lorene::Time_slice::jtime, Lorene::max(), Lorene::maxabs(), Lorene::Time_slice::p_gamma, p_hdirac, p_tgamma, and Lorene::Time_slice::the_time.
|
virtual |
Sets the factor at the current time step (jtime
) and deletes the value of N.
Definition at line 481 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, npsi_evol, and Lorene::Time_slice::the_time.
|
virtual |
Sets the factor at the current time step (jtime
) and deletes the value of .
Definition at line 456 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, npsi_evol, Lorene::Time_slice::p_gamma, p_ln_psi, p_psi4, psi_evol, and Lorene::Time_slice::the_time.
|
virtual |
Sets the conformal factor relating the physical metric to the conformal one: .
is defined by
Sets the value at the current time step (jtime
) and deletes the value of N.
Definition at line 431 of file time_slice_conf.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::n_evol, Lorene::Time_slice::p_gamma, p_ln_psi, p_psi4, psi_evol, and Lorene::Time_slice::the_time.
|
virtual |
Sets the conformal factor relating the physical metric to the conformal one: .
is defined by
Sets the value at the current time step (jtime
) and deletes the value of .
Definition at line 407 of file time_slice_conf.C.
References Lorene::Time_slice::adm_mass_evol, Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, npsi_evol, Lorene::Time_slice::p_gamma, p_ln_psi, p_psi4, psi_evol, and Lorene::Time_slice::the_time.
|
inlineinherited |
Sets the order of the finite-differences scheme.
Definition at line 334 of file time_slice.h.
References Lorene::Time_slice::scheme_order.
|
virtual |
Conformal metric Returns the value at the current time step (jtime
).
Reimplemented in Lorene::Isol_hor.
Definition at line 750 of file time_slice_conf.C.
References Lorene::Metric_flat::con(), ff, hh(), and p_tgamma.
|
virtual |
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 799 of file time_slice_conf.C.
References Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Vector::divergence(), ff, Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice::jtime, ln_psi(), nn(), psi(), psi_evol, Lorene::Time_slice::scheme_order, Lorene::Time_slice::the_time, and Lorene::Time_slice::trk_evol.
|
virtual |
Vector representing the longitudinal part of .
(see the documentation of hata_evol
)
Definition at line 828 of file time_slice_conf.C.
References ff, hata_evol, Lorene::Time_slice::jtime, and p_vec_X.
|
mutableprotected |
Potential A associated with the symmetric tensor .
(see documentation of Sym_tensor::p_aaa
).
Definition at line 550 of file time_slice.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 236 of file time_slice.h.
|
mutableprotected |
Potential associated with the symmetric tensor .
(see documentation of Sym_tensor::p_tilde_b
).
Definition at line 555 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the shift vector .
Definition at line 222 of file time_slice.h.
|
protectedinherited |
Number of stored time slices.
Definition at line 182 of file time_slice.h.
|
protected |
Pointer on the flat metric with respect to which the conformal decomposition is performed.
Definition at line 510 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the covariant components of the induced metric .
Definition at line 201 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the contravariant components of the induced metric .
Definition at line 206 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the components .
It is the conformal representation of the traceless part of the extrinsic curvature: . One can uniquely (up to the boundary conditions) define the decomposition: , where represents the longitudinal part and is the transverse-traceless part.
Definition at line 545 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the components .
It is the deviation of the conformal metric from the flat metric : .
Definition at line 533 of file time_slice.h.
|
protectedinherited |
Time step index of the latest slice.
Definition at line 193 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the covariant components of the extrinsic curvature tensor .
Definition at line 211 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the contravariant components of the extrinsic curvature tensor .
Definition at line 216 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the lapse function N.
Definition at line 219 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the factor .
Definition at line 525 of file time_slice.h.
|
mutableprotectedinherited |
Pointer on the induced metric at the current time step (jtime
)
Definition at line 242 of file time_slice.h.
|
mutableprotected |
Pointer on the vector (which vanishes in Dirac gauge), at the current time step (jtime
).
Definition at line 574 of file time_slice.h.
|
mutableprotected |
Pointer on the logarithm of at the current time step (jtime
)
Definition at line 569 of file time_slice.h.
|
mutableprotected |
Pointer on the factor at the current time step (jtime
)
Definition at line 566 of file time_slice.h.
|
mutableprotected |
Pointer on the conformal metric at the current time step (jtime
)
Definition at line 563 of file time_slice.h.
|
mutableprotected |
Pointer on the vector representing the longitudinal part of .
(see the documentation of hata_evol
)
Definition at line 580 of file time_slice.h.
|
mutableprotected |
Values at successive time steps of the conformal factor relating the physical metric to the conformal one: .
is defined by
Definition at line 520 of file time_slice.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 190 of file time_slice.h.
|
protectedinherited |
Time label of each slice.
Definition at line 196 of file time_slice.h.
|
mutableprotectedinherited |
Values at successive time steps of the trace K of the extrinsic curvature.
Definition at line 227 of file time_slice.h.