Lorene::Time_slice_conf Class Reference
[Time evolution (under development)]

Spacelike time slice of a 3+1 spacetime with conformal decomposition. More...

#include <time_slice.h>

Inheritance diagram for Lorene::Time_slice_conf:
Lorene::Time_slice Lorene::Isol_hor Lorene::Tslice_dirac_max

List of all members.

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.
 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.
 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).
 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.
 Time_slice_conf (const Time_slice_conf &)
 Copy constructor.
virtual ~Time_slice_conf ()
 Destructor.
void operator= (const Time_slice_conf &)
 Assignment to another Time_slice_conf.
void operator= (const Time_slice &)
 Assignment to a Time_slice.
virtual void set_psi_del_npsi (const Scalar &psi_in)
 Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
virtual void set_psi_del_n (const Scalar &psi_in)
 Sets the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
virtual void set_npsi_del_psi (const Scalar &npsi_in)
 Sets the factor $ N\Psi $ at the current time step (jtime ) and deletes the value of $\Psi$.
virtual void set_npsi_del_n (const Scalar &npsi_in)
 Sets the factor $ N\Psi $ at the current time step (jtime ) and deletes the value of N.
virtual void set_hh (const Sym_tensor &hh_in)
 Sets the deviation $ h^{ij} $ of the conformal metric $ \tilde\gamma^{ij} $ from the flat metric $ f^{ij} $: $\tilde\gamma^{ij} = f^{ij} + h^{ij} $.
virtual void set_hata (const Sym_tensor &hata_in)
 Sets the conformal representation $ \hat{A}{ij} $ of the traceless part of the extrinsic curvature: $ \hat{A}^{ij} = \Psi^{10} \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.
virtual void set_hata_TT (const Sym_tensor_tt &hata_tt)
 Sets the TT part of $ \hat{A}^{ij} $ (see member hata_evol ).
virtual void set_hata_from_XAB (Param *par_bc=0x0, Param *par_mat=0x0)
 Sets the conformal representation $ \hat{A}{ij} $ of the traceless part of the extrinsic curvature from its potentials A, $ \tilde{B} $ and $ X^i $.
virtual const Scalarnn () const
 Lapse function N at the current time step (jtime ).
virtual const Sym_tensorgam_dd () const
 Induced metric (covariant components $ \gamma_{ij} $) at the current time step (jtime ).
virtual const Sym_tensorgam_uu () const
 Induced metric (contravariant components $ \gamma^{ij} $) at the current time step (jtime ).
virtual const Sym_tensork_dd () const
 Extrinsic curvature tensor (covariant components $ K_{ij} $) at the current time step (jtime ).
virtual const Sym_tensork_uu () const
 Extrinsic curvature tensor (contravariant components $ K^{ij} $) at the current time step (jtime ).
virtual const ScalarA_hata () const
 Returns the potential A of $ \hat{A}^{ij} $.
virtual const ScalarB_hata () const
 Returns the potential $\tilde{B}$ of $ \hat{A}^{ij} $.
virtual const Scalarpsi () const
 Conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
const Scalarpsi4 () const
 Factor $ \Psi^4 $ at the current time step (jtime ).
const Scalarln_psi () const
 Logarithm of $ \Psi $ at the current time step (jtime ).
virtual const Scalarnpsi () const
 Factor $ N\Psi $ at the current time step (jtime ).
virtual const Metrictgam () const
 Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $ Returns the value at the current time step (jtime ).
virtual const Sym_tensorhh (Param *=0x0, Param *=0x0) const
 Deviation $ h^{ij} $ of the conformal metric $ \tilde\gamma^{ij} $ from the flat metric $ f^{ij} $: $\tilde\gamma^{ij} = f^{ij} + h^{ij} $.
virtual const Sym_tensorhata () const
 Conformal representation $ \hat{A}^{ij} $ of the traceless part of the extrinsic curvature: $ \hat{A}^{ij} = \Psi^{10} \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.
virtual Sym_tensor aa () const
 Conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature: $ A^{ij} = \Psi^4 \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $.
virtual const Scalartrk () const
 Trace K of the extrinsic curvature at the current time step (jtime ).
virtual const Vectorhdirac () const
 Vector $ H^i = {\cal D}_j \tilde\gamma^{ij} $ which vanishes in Dirac gauge.
virtual const Vectorvec_X (int method_poisson=6) const
 Vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $.
void compute_X_from_momentum_constraint (const Vector &hat_S, const Sym_tensor_tt &hata_tt, int iter_max=200, double precis=1.e-12, double relax=0.8, int methode_poisson=6)
 Computes the vector $ X^i $ from the conformally-rescaled momentum $ \hat{S}^i = \Psi^6 S^i $, using the momentum constraint.
virtual void set_AB_hata (const Scalar &A_in, const Scalar &B_in)
 Sets the potentials A and $\tilde{B}$ of the TT part $ \hat{A}^{ij} $ (see the documentation of Sym_tensor for details).
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.
virtual double adm_mass () const
 Returns the ADM mass (geometrical units) at the current step.
void check_psi_dot (Tbl &tlnpsi_dot, Tbl &tdiff, Tbl &tdiff_rel) const
 Checks the $\frac{\partial}{\partial t} \ln\Psi $ relation.
void set_scheme_order (int ord)
 Sets the order of the finite-differences scheme.
int get_scheme_order () const
 Gets the order of the finite-differences scheme.
int get_latest_j () const
 Gets the latest value of time step index.
const Evolution_std< double > & get_time () const
 Gets the time coordinate t at successive time steps.
virtual const Vectorbeta () const
 shift vector $ \beta^i $ at the current time step (jtime )
const Metricgam () const
 Induced metric $ \mathbf{\gamma} $ at the current time step (jtime ).
Tbl check_hamiltonian_constraint (const Scalar *energy_density=0x0, ostream &ost=cout, bool verb=true) const
 Checks the level at which the hamiltonian constraint is verified.
Tbl check_momentum_constraint (const Vector *momentum_density=0x0, ostream &ost=cout, bool verb=true) const
 Checks the level at which the momentum constraints are verified.
Tbl check_dynamical_equations (const Sym_tensor *strain_tensor=0x0, const Scalar *energy_density=0x0, ostream &ost=cout, bool verb=true) const
 Checks the level at which the dynamical equations are verified.
void save (const char *rootname) const
 Saves in a binary file.

Protected Member Functions

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

const Metric_flatff
 Pointer on the flat metric $ f_{ij} $ with respect to which the conformal decomposition is performed.
Evolution_std< Scalarpsi_evol
 Values at successive time steps of the conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $.
Evolution_std< Scalarnpsi_evol
 Values at successive time steps of the factor $ N\Psi $.
Evolution_std< Sym_tensorhh_evol
 Values at successive time steps of the components $ h^{ij} $.
Evolution_std< Sym_tensorhata_evol
 Values at successive time steps of the components $ \hat{A}^{ij} $.
Evolution_std< ScalarA_hata_evol
 Potential A associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $.
Evolution_std< ScalarB_hata_evol
 Potential $ \tilde{B} $ associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $.
Metricp_tgamma
 Pointer on the conformal metric $ \tilde\gamma_{ij} $ at the current time step (jtime).
Scalarp_psi4
 Pointer on the factor $ \Psi^4 $ at the current time step (jtime).
Scalarp_ln_psi
 Pointer on the logarithm of $ \Psi $ at the current time step (jtime).
Vectorp_hdirac
 Pointer on the vector $ H^i = {\cal D}_j \tilde\gamma^{ij} $ (which vanishes in Dirac gauge), at the current time step (jtime).
Vectorp_vec_X
 Pointer on the vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $.
int depth
 Number of stored time slices.
int scheme_order
 Order of the finite-differences scheme for the computation of time derivatives.
int jtime
 Time step index of the latest slice.
Evolution_std< double > the_time
 Time label of each slice.
Evolution_std< Sym_tensorgam_dd_evol
 Values at successive time steps of the covariant components of the induced metric $ \gamma_{ij} $.
Evolution_std< Sym_tensorgam_uu_evol
 Values at successive time steps of the contravariant components of the induced metric $ \gamma^{ij} $.
Evolution_std< Sym_tensork_dd_evol
 Values at successive time steps of the covariant components of the extrinsic curvature tensor $ K_{ij} $.
Evolution_std< Sym_tensork_uu_evol
 Values at successive time steps of the contravariant components of the extrinsic curvature tensor $ K^{ij} $.
Evolution_std< Scalarn_evol
 Values at successive time steps of the lapse function N.
Evolution_std< Vectorbeta_evol
 Values at successive time steps of the shift vector $ \beta^i $.
Evolution_std< Scalartrk_evol
 Values at successive time steps of the trace K of the extrinsic curvature.
Evolution_full< Tbladm_mass_evol
 ADM mass at each time step, since the creation of the slice.
Metricp_gamma
 Pointer on the induced metric at the current time step (jtime).

Friends

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

Detailed Description

Spacelike time slice of a 3+1 spacetime with conformal decomposition.

()

Definition at line 501 of file time_slice.h.


Constructor & Destructor Documentation

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.

Parameters:
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 $\Psi$ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $
hh_in deviation $ h^{ij} $ of the conformal metric $ \tilde\gamma^{ij} $ from the flat metric $ f^{ij} $: $\tilde\gamma^{ij} = f^{ij} + h^{ij} $. $ h^{ij} $ must be such that $\det\tilde\gamma^{ij} = f^{-1} $.
hata_in conformal representation $ A^{ij} $ of the traceless part of the extrinsic curvature: $ \hat{A}^{ij} = \Psi^{10} \left( K^{ij} - \frac{1}{3} K \gamma^{ij} \right) $
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, Lorene::Time_slice::trk_evol, and Lorene::Evolution_std< TyT >::update().

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.

Parameters:
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_flat::con(), Lorene::Metric::con(), Lorene::Metric_flat::determinant(), Lorene::Metric::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, Lorene::Time_slice::trk_evol, and Lorene::Evolution_std< TyT >::update().

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

Parameters:
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(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

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.

Parameters:
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(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

Lorene::Time_slice_conf::Time_slice_conf ( const Time_slice_conf tin  ) 

Copy constructor.

Definition at line 323 of file time_slice_conf.C.

References set_der_0x0().

Lorene::Time_slice_conf::~Time_slice_conf (  )  [virtual]

Destructor.

Definition at line 340 of file time_slice_conf.C.

References del_deriv().


Member Function Documentation

const Scalar & Lorene::Time_slice_conf::A_hata (  )  const [virtual]

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

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

Definition at line 667 of file time_slice_conf.C.

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

Sym_tensor Lorene::Time_slice_conf::aa (  )  const [virtual]

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

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

Definition at line 768 of file time_slice_conf.C.

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

double Lorene::Time_slice_conf::adm_mass (  )  const [virtual]
const Scalar & Lorene::Time_slice_conf::B_hata (  )  const [virtual]

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

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

Definition at line 681 of file time_slice_conf.C.

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

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

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

Definition at line 90 of file time_slice_access.C.

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

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

Checks the level at which the dynamical equations are verified.

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

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

Definition at line 142 of file tslice_check_einstein.C.

References Lorene::Tensor::annule_domain(), 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::Evolution< TyT >::time_derive(), Lorene::Tensor::trace(), Lorene::Time_slice::trk(), and Lorene::Tensor::up_down().

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

Checks the level at which the hamiltonian constraint is verified.

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

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

Definition at line 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().

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

Checks the level at which the momentum constraints are verified.

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

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

Definition at line 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 $\frac{\partial}{\partial t} \ln\Psi $ relation.

Parameters:
tlnpsi_dot [output] maximun value in each domain of $ \left| \frac{\partial}{\partial t} \ln\Psi \right| $
tdiff [output] maximum value in each domain of $ \left| \frac{\partial}{\partial t} \ln\Psi - \beta^i {\cal D}_i \ln \Psi - \frac{1}{6} ( {\cal D}_i \beta^i - N K) \right| $
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::Evolution< TyT >::is_known(), 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, Lorene::Evolution< TyT >::time_derive(), 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 $ X^i $ from the conformally-rescaled momentum $ \hat{S}^i = \Psi^6 S^i $, 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().

void Lorene::Time_slice_conf::del_deriv (  )  const [protected, virtual]

Deletes all the derived quantities.

Reimplemented from Lorene::Time_slice.

Definition at line 350 of file time_slice_conf.C.

References p_hdirac, p_ln_psi, p_psi4, p_tgamma, p_vec_X, and set_der_0x0().

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

Induced metric $ \mathbf{\gamma} $ 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.

const Sym_tensor & Lorene::Time_slice_conf::gam_dd (  )  const [virtual]

Induced metric (covariant components $ \gamma_{ij} $) 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::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, psi4(), tgam(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

const Sym_tensor & Lorene::Time_slice_conf::gam_uu (  )  const [virtual]

Induced metric (contravariant components $ \gamma^{ij} $) 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::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, psi4(), tgam(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

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

Gets the latest value of time step index.

Definition at line 346 of file time_slice.h.

References Lorene::Time_slice::jtime.

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

Gets the order of the finite-differences scheme.

Definition at line 343 of file time_slice.h.

References Lorene::Time_slice::scheme_order.

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

Gets the time coordinate t at successive time steps.

Definition at line 349 of file time_slice.h.

References Lorene::Time_slice::the_time.

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

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

Reimplemented in Lorene::Tslice_dirac_max.

Definition at line 818 of file time_slice_conf.C.

References ff, hh(), and p_hdirac.

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

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

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

Reimplemented in Lorene::Tslice_dirac_max.

Definition at line 761 of file time_slice_conf.C.

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

void Lorene::Time_slice_conf::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 
) [virtual]

Computes valid initial data by solving the constraint equations in the conformal thin-sandwich approach.

Parameters:
uu value of $ {\tilde u}^{ij} = \partial h^{ij} /\partial t $ (freely specifiable data). This quantity must be trace-free with respect to the conformal metric $\tilde\gamma_{ij}$, reflecting the unimodular character of $\tilde\gamma_{ij}$.
trk_in value of $ K = K_i^{\ i} $ (freely specifiable data)
trk_point value of $ \partial K / \partial t $ (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 A_hata(), A_hata_evol, aa(), B_hata(), B_hata_evol, Lorene::Time_slice::beta(), Lorene::Time_slice::beta_evol, Lorene::Scalar::check_dzpuis(), Lorene::contract(), Lorene::Tensor::dec_dzpuis(), del_deriv(), Lorene::Time_slice::depth, Lorene::Tensor::derive_con(), Lorene::Scalar::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Vector::derive_lie(), des_meridian(), Lorene::diffrel(), Lorene::Vector::divergence(), Lorene::Sym_tensor::divergence(), Lorene::Evolution< TyT >::downdate(), ff, Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::get_triad(), hata_evol, hdirac(), hh(), hh_evol, Lorene::Tensor::inc_dzpuis(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Scalar::laplacian(), ln_psi(), Lorene::max(), Lorene::maxabs(), Lorene::Time_slice::n_evol, nn(), npsi_evol, Lorene::Vector::ope_killing_conf(), Lorene::Vector::poisson(), Lorene::Scalar::poisson(), psi(), psi4(), psi_evol, Lorene::Metric::ricci_scal(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::set_etat_zero(), set_hata(), set_psi_del_npsi(), Lorene::Scalar::std_spectral_base(), tgam(), Lorene::Time_slice::the_time, Lorene::Tensor::trace(), trk(), Lorene::Time_slice::trk_evol, Lorene::Tensor::up_down(), Lorene::Evolution_std< TyT >::update(), and Lorene::Map::val_r().

const Sym_tensor & Lorene::Time_slice_conf::k_dd (  )  const [virtual]

Extrinsic curvature tensor (covariant components $ K_{ij} $) 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::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, k_uu(), Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

const Sym_tensor & Lorene::Time_slice_conf::k_uu (  )  const [virtual]

Extrinsic curvature tensor (contravariant components $ K^{ij} $) 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::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::Time_slice::k_uu_evol, psi(), psi4(), Lorene::Time_slice::the_time, trk(), and Lorene::Evolution_std< TyT >::update().

const Scalar & Lorene::Time_slice_conf::ln_psi (  )  const

Logarithm of $ \Psi $ 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().

const Scalar & Lorene::Time_slice_conf::nn (  )  const [virtual]
const Scalar & Lorene::Time_slice_conf::npsi (  )  const [virtual]
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 operator=().

void Lorene::Time_slice_conf::operator= ( const Time_slice_conf tin  ) 

Assignment to another Time_slice_conf.

Reimplemented from Lorene::Time_slice.

Reimplemented in Lorene::Isol_hor, and Lorene::Tslice_dirac_max.

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, and psi_evol.

ostream & Lorene::Time_slice_conf::operator>> ( ostream &  flux  )  const [protected, virtual]

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

Reimplemented from Lorene::Time_slice.

Reimplemented in Lorene::Isol_hor, and Lorene::Tslice_dirac_max.

Definition at line 944 of file time_slice_conf.C.

References ff, Lorene::Metric_flat::get_triad(), hata_evol, hh_evol, Lorene::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, Lorene::maxabs(), npsi_evol, p_hdirac, p_ln_psi, p_psi4, p_tgamma, p_vec_X, and psi_evol.

const Scalar & Lorene::Time_slice_conf::psi (  )  const [virtual]

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

$ \Psi $ is defined by

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

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

Definition at line 696 of file time_slice_conf.C.

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

const Scalar & Lorene::Time_slice_conf::psi4 (  )  const

Factor $ \Psi^4 $ 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().

void Lorene::Time_slice_conf::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.

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

Reimplemented from Lorene::Time_slice.

Reimplemented in Lorene::Isol_hor, and Lorene::Tslice_dirac_max.

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::Evolution< TyT >::is_known(), Lorene::Time_slice::jtime, psi(), and psi_evol.

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

Saves in a binary file.

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

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

Definition at line 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::Time_slice::sauve(), Lorene::Base_vect::sauve(), Lorene::Map::sauve(), and Lorene::Mg3d::sauve().

void Lorene::Time_slice_conf::set_AB_hata ( const Scalar A_in,
const Scalar B_in 
) [virtual]

Sets the potentials A and $\tilde{B}$ of the TT part $ \hat{A}^{ij} $ (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, Lorene::Evolution< TyT >::downdate(), 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().

void Lorene::Time_slice_conf::set_der_0x0 (  )  const [protected]

Sets to 0x0 all the pointers on derived quantities.

Reimplemented from Lorene::Time_slice.

Definition at line 364 of file time_slice_conf.C.

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

void Lorene::Time_slice_conf::set_hata ( const Sym_tensor hata_in  )  [virtual]

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

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

Definition at line 537 of file time_slice_conf.C.

References A_hata_evol, B_hata_evol, Lorene::Evolution< TyT >::downdate(), 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().

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

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

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

Definition at line 569 of file time_slice_conf.C.

References A_hata_evol, B_hata_evol, Lorene::Evolution< TyT >::downdate(), ff, Lorene::Tensor::get_mp(), hata_evol, Lorene::Map::inc_dzpuis(), 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(), p_vec_X, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

void Lorene::Time_slice_conf::set_hata_TT ( const Sym_tensor_tt hata_tt  )  [virtual]
void Lorene::Time_slice_conf::set_hh ( const Sym_tensor hh_in  )  [virtual]
void Lorene::Time_slice_conf::set_npsi_del_n ( const Scalar npsi_in  )  [virtual]
void Lorene::Time_slice_conf::set_npsi_del_psi ( const Scalar npsi_in  )  [virtual]
void Lorene::Time_slice_conf::set_psi_del_n ( const Scalar psi_in  )  [virtual]

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

$ \Psi $ is defined by

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

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

Definition at line 431 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, p_ln_psi, p_psi4, psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

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

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

$ \Psi $ is defined by

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

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

Definition at line 407 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, npsi_evol, Lorene::Time_slice::p_gamma, p_ln_psi, p_psi4, psi_evol, Lorene::Time_slice::the_time, and Lorene::Evolution_std< TyT >::update().

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

Sets the order of the finite-differences scheme.

Definition at line 334 of file time_slice.h.

References Lorene::Time_slice::scheme_order.

const Metric & Lorene::Time_slice_conf::tgam (  )  const [virtual]

Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $ 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.

const Scalar & Lorene::Time_slice_conf::trk (  )  const [virtual]
const Vector & Lorene::Time_slice_conf::vec_X ( int  method_poisson = 6  )  const [virtual]

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

(see the documentation of hata_evol)

Definition at line 828 of file time_slice_conf.C.

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


Friends And Related Function Documentation

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

Display.


Member Data Documentation

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

(see documentation of Sym_tensor::p_aaa).

Definition at line 550 of file time_slice.h.

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

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

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

Definition at line 236 of file time_slice.h.

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

(see documentation of Sym_tensor::p_tilde_b).

Definition at line 555 of file time_slice.h.

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

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

Definition at line 222 of file time_slice.h.

int Lorene::Time_slice::depth [protected, inherited]

Number of stored time slices.

Definition at line 182 of file time_slice.h.

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

Definition at line 510 of file time_slice.h.

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

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

Definition at line 201 of file time_slice.h.

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

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

Definition at line 206 of file time_slice.h.

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

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

Definition at line 545 of file time_slice.h.

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

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

Definition at line 533 of file time_slice.h.

int Lorene::Time_slice::jtime [protected, inherited]

Time step index of the latest slice.

Definition at line 193 of file time_slice.h.

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

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

Definition at line 211 of file time_slice.h.

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

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

Definition at line 216 of file time_slice.h.

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

Values at successive time steps of the lapse function N.

Definition at line 219 of file time_slice.h.

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

Definition at line 525 of file time_slice.h.

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

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

Definition at line 242 of file time_slice.h.

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

Definition at line 574 of file time_slice.h.

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

Definition at line 569 of file time_slice.h.

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

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

Definition at line 566 of file time_slice.h.

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

Definition at line 563 of file time_slice.h.

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

(see the documentation of hata_evol)

Definition at line 580 of file time_slice.h.

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

$ \Psi $ is defined by

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

Definition at line 520 of file time_slice.h.

int Lorene::Time_slice::scheme_order [protected, inherited]

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

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

Definition at line 190 of file time_slice.h.

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

Time label of each slice.

Definition at line 196 of file time_slice.h.

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

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

Definition at line 227 of file time_slice.h.


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

Generated on 7 Dec 2019 for LORENE by  doxygen 1.6.1