LORENE
Lorene::Time_slice Class Reference

Spacelike time slice of a 3+1 spacetime. More...

#include <time_slice.h>

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

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). More...
 
 Time_slice (const Scalar &lapse_in, const Vector &shift_in, const Evolution_std< Sym_tensor > &gamma_in)
 General constructor (Lagrangian-like). More...
 
 Time_slice (const Map &mp, const Base_vect &triad, int depth_in=3)
 Constructor as standard time slice of flat spacetime (Minkowski). More...
 
 Time_slice (const Map &mp, const Base_vect &triad, FILE *fich, bool partial_read, int depth_in=3)
 Constructor from binary file. More...
 
 Time_slice (const Time_slice &)
 Copy constructor. More...
 
virtual ~Time_slice ()
 Destructor. More...
 
void operator= (const Time_slice &)
 Assignment to another Time_slice. 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 Scalarnn () const
 Lapse function N at the current time step (jtime ) More...
 
virtual const Vectorbeta () const
 shift vector $ \beta^i $ at the current time step (jtime ) More...
 
const Metricgam () const
 Induced metric $ \mathbf{\gamma} $ at the current time step (jtime ) More...
 
virtual const Sym_tensorgam_dd () const
 Induced metric (covariant components $ \gamma_{ij} $) at the current time step (jtime ) More...
 
virtual const Sym_tensorgam_uu () const
 Induced metric (contravariant components $ \gamma^{ij} $) at the current time step (jtime ) More...
 
virtual const Sym_tensork_dd () const
 Extrinsic curvature tensor (covariant components $ K_{ij} $) at the current time step (jtime ) More...
 
virtual const Sym_tensork_uu () const
 Extrinsic curvature tensor (contravariant components $ K^{ij} $) at the current time step (jtime ) More...
 
virtual const Scalartrk () const
 Trace K of the extrinsic curvature 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...
 
virtual double adm_mass () const
 Returns the ADM mass (geometrical units) at the current step. More...
 
void save (const char *rootname) const
 Saves in a binary file. More...
 

Protected Member Functions

 Time_slice (int depth_in)
 Special constructor for derived classes. More...
 
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

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_tensorgam_dd_evol
 Values at successive time steps of the covariant components of the induced metric $ \gamma_{ij} $. More...
 
Evolution_std< Sym_tensorgam_uu_evol
 Values at successive time steps of the contravariant components of the induced metric $ \gamma^{ij} $. More...
 
Evolution_std< Sym_tensork_dd_evol
 Values at successive time steps of the covariant components of the extrinsic curvature tensor $ K_{ij} $. More...
 
Evolution_std< Sym_tensork_uu_evol
 Values at successive time steps of the contravariant components of the extrinsic curvature tensor $ K^{ij} $. More...
 
Evolution_std< Scalarn_evol
 Values at successive time steps of the lapse function N. More...
 
Evolution_std< Vectorbeta_evol
 Values at successive time steps of the shift vector $ \beta^i $. More...
 
Evolution_std< Scalartrk_evol
 Values at successive time steps of the trace K of the extrinsic curvature. More...
 
Evolution_full< Tbladm_mass_evol
 ADM mass at each time step, since the creation of the slice. More...
 
Metricp_gamma
 Pointer on the induced metric at the current time step (jtime) More...
 

Friends

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

Detailed Description

Spacelike time slice of a 3+1 spacetime.

()

Definition at line 176 of file time_slice.h.

Constructor & Destructor Documentation

◆ Time_slice() [1/6]

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

Parameters
lapse_inlapse function N
shift_inshift vector
gamma_ininduced metric (covariant or contravariant components)
kk_inextrinsic curvature (covariant or contravariant components)
depth_innumber 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 115 of file time_slice.C.

References gam_dd_evol, gam_uu_evol, Lorene::Tensor::get_index_type(), jtime, k_dd_evol, k_uu_evol, p_gamma, set_der_0x0(), the_time, Lorene::Tensor::trace(), and trk_evol.

◆ Time_slice() [2/6]

Lorene::Time_slice::Time_slice ( const Scalar lapse_in,
const Vector shift_in,
const Evolution_std< Sym_tensor > &  gamma_in 
)

General constructor (Lagrangian-like).

Parameters
lapse_inlapse function N
shift_inshift vector
gamma_ininduced 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 158 of file time_slice.C.

References set_der_0x0().

◆ Time_slice() [3/6]

Lorene::Time_slice::Time_slice ( const Map mp,
const Base_vect triad,
int  depth_in = 3 
)

Constructor as standard time slice of flat spacetime (Minkowski).

Parameters
mpMapping on which the various Lorene fields will be constructed
triadvector basis with respect to which the various tensor components will be defined
depth_innumber 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 183 of file time_slice.C.

References beta_evol, Lorene::Metric_flat::cov(), Lorene::Map::flat_met_cart(), Lorene::Map::flat_met_spher(), gam_dd_evol, jtime, k_dd_evol, n_evol, set_der_0x0(), Lorene::Scalar::set_etat_one(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::std_spectral_base(), the_time, and trk_evol.

◆ Time_slice() [4/6]

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

Parameters
mpMapping on which the various Lorene fields will be constructed
triadvector basis with respect to which the various tensor components will be defined
fichfile containing the saved Time_slice
partial_readindicates whether the full object must be read in file or whether the final construction is devoted to a constructor of a derived class
depth_innumber of stored time slices; the given value must coincide with that stored in the file.

Definition at line 238 of file time_slice.C.

References beta_evol, depth, Lorene::fread_be(), Lorene::Map::get_mg(), jtime, n_evol, scheme_order, set_der_0x0(), the_time, and Lorene::Evolution_std< TyT >::update().

◆ Time_slice() [5/6]

Lorene::Time_slice::Time_slice ( const Time_slice tin)

Copy constructor.

Definition at line 317 of file time_slice.C.

References set_der_0x0().

◆ Time_slice() [6/6]

Lorene::Time_slice::Time_slice ( int  depth_in)
explicitprotected

Special constructor for derived classes.

Definition at line 336 of file time_slice.C.

References set_der_0x0().

◆ ~Time_slice()

Lorene::Time_slice::~Time_slice ( )
virtual

Destructor.

Definition at line 360 of file time_slice.C.

References del_deriv().

Member Function Documentation

◆ adm_mass()

double Lorene::Time_slice::adm_mass ( ) const
virtual

◆ beta()

const Vector & Lorene::Time_slice::beta ( ) const
virtual

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

Definition at line 90 of file time_slice_access.C.

References beta_evol, and jtime.

◆ check_dynamical_equations()

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

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::Tbl::annule_hard(), beta(), Lorene::Metric::con(), Lorene::contract(), Lorene::Tensor::derive_con(), Lorene::Scalar::derive_con(), Lorene::Sym_tensor::derive_lie(), gam(), Lorene::Tbl::get_dim(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), jtime, k_dd(), k_dd_evol, k_uu(), Lorene::maxabs(), nn(), Lorene::Tensor_sym::position(), Lorene::Metric::ricci(), scheme_order, Lorene::Itbl::set(), Lorene::Tensor::trace(), trk(), and Lorene::Tensor::up_down().

◆ check_hamiltonian_constraint()

Tbl Lorene::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.

\[ 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(), gam(), k_dd(), k_uu(), Lorene::maxabs(), Lorene::Metric::ricci_scal(), and trk().

◆ check_momentum_constraint()

Tbl Lorene::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.

\[ 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(), gam(), Lorene::Tensor::get_index_type(), k_uu(), Lorene::maxabs(), and trk().

◆ del_deriv()

void Lorene::Time_slice::del_deriv ( ) const
protectedvirtual

Deletes all the derived quantities.

Reimplemented in Lorene::Time_slice_conf.

Definition at line 370 of file time_slice.C.

References adm_mass_evol, jtime, p_gamma, and set_der_0x0().

◆ gam()

const Metric & Lorene::Time_slice::gam ( ) const

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

Definition at line 98 of file time_slice_access.C.

References gam_dd(), and p_gamma.

◆ gam_dd()

const Sym_tensor & Lorene::Time_slice::gam_dd ( ) const
virtual

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

Reimplemented in Lorene::Time_slice_conf.

Definition at line 110 of file time_slice_access.C.

References Lorene::Metric::cov(), gam_dd_evol, gam_uu_evol, jtime, p_gamma, and the_time.

◆ gam_uu()

const Sym_tensor & Lorene::Time_slice::gam_uu ( ) const
virtual

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

Reimplemented in Lorene::Time_slice_conf.

Definition at line 125 of file time_slice_access.C.

References gam(), gam_dd_evol, gam_uu_evol, jtime, and the_time.

◆ get_latest_j()

int Lorene::Time_slice::get_latest_j ( ) const
inline

Gets the latest value of time step index.

Definition at line 346 of file time_slice.h.

References jtime.

◆ get_scheme_order()

int Lorene::Time_slice::get_scheme_order ( ) const
inline

Gets the order of the finite-differences scheme.

Definition at line 343 of file time_slice.h.

References scheme_order.

◆ get_time()

const Evolution_std<double>& Lorene::Time_slice::get_time ( ) const
inline

Gets the time coordinate t at successive time steps.

Definition at line 349 of file time_slice.h.

References the_time.

◆ k_dd()

const Sym_tensor & Lorene::Time_slice::k_dd ( ) const
virtual

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

Reimplemented in Lorene::Time_slice_conf.

Definition at line 138 of file time_slice_access.C.

References beta(), Lorene::Tensor::down(), gam(), gam_dd(), gam_dd_evol, jtime, k_dd_evol, nn(), Lorene::Vector::ope_killing(), scheme_order, and the_time.

◆ k_uu()

const Sym_tensor & Lorene::Time_slice::k_uu ( ) const
virtual

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

Reimplemented in Lorene::Time_slice_conf.

Definition at line 160 of file time_slice_access.C.

References beta(), gam(), gam_uu(), gam_uu_evol, jtime, k_uu_evol, nn(), Lorene::Vector::ope_killing(), scheme_order, and the_time.

◆ nn()

const Scalar & Lorene::Time_slice::nn ( ) const
virtual

Lapse function N at the current time step (jtime )

Reimplemented in Lorene::Time_slice_conf.

Definition at line 82 of file time_slice_access.C.

References jtime, and n_evol.

◆ operator=()

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

Assignment to another Time_slice.

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

◆ operator>>()

ostream & Lorene::Time_slice::operator>> ( ostream &  flux) const
protectedvirtual

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

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

Definition at line 414 of file time_slice.C.

References adm_mass(), adm_mass_evol, beta_evol, depth, gam_dd_evol, gam_uu_evol, jtime, k_dd_evol, k_uu_evol, Lorene::maxabs(), n_evol, p_gamma, scheme_order, the_time, and trk_evol.

◆ sauve()

void Lorene::Time_slice::sauve ( FILE *  fich,
bool  partial_save 
) const
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.

Parameters
fichbinary file
partial_saveindicates whether the whole object must be saved.

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

Definition at line 510 of file time_slice.C.

References beta(), beta_evol, depth, Lorene::fwrite_be(), Lorene::Evolution< TyT >::is_known(), jtime, n_evol, nn(), scheme_order, and the_time.

◆ save()

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

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

Definition at line 464 of file time_slice.C.

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

◆ set_der_0x0()

void Lorene::Time_slice::set_der_0x0 ( ) const
protected

Sets to 0x0 all the pointers on derived quantities.

Definition at line 380 of file time_slice.C.

References p_gamma.

◆ set_scheme_order()

void Lorene::Time_slice::set_scheme_order ( int  ord)
inline

Sets the order of the finite-differences scheme.

Definition at line 334 of file time_slice.h.

References scheme_order.

◆ trk()

const Scalar & Lorene::Time_slice::trk ( ) const
virtual

Trace K of the extrinsic curvature at the current time step (jtime )

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

Definition at line 180 of file time_slice_access.C.

References gam(), jtime, k_dd(), k_uu(), k_uu_evol, the_time, and trk_evol.

Friends And Related Function Documentation

◆ operator<<

ostream& operator<< ( ostream &  flux,
const Time_slice sigma 
)
friend

Display.

Definition at line 456 of file time_slice.C.

Member Data Documentation

◆ adm_mass_evol

Evolution_full<Tbl> Lorene::Time_slice::adm_mass_evol
mutableprotected

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.

◆ beta_evol

Evolution_std<Vector> Lorene::Time_slice::beta_evol
mutableprotected

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

Definition at line 222 of file time_slice.h.

◆ depth

int Lorene::Time_slice::depth
protected

Number of stored time slices.

Definition at line 182 of file time_slice.h.

◆ gam_dd_evol

Evolution_std<Sym_tensor> Lorene::Time_slice::gam_dd_evol
mutableprotected

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

Definition at line 201 of file time_slice.h.

◆ gam_uu_evol

Evolution_std<Sym_tensor> Lorene::Time_slice::gam_uu_evol
mutableprotected

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

Definition at line 206 of file time_slice.h.

◆ jtime

int Lorene::Time_slice::jtime
protected

Time step index of the latest slice.

Definition at line 193 of file time_slice.h.

◆ k_dd_evol

Evolution_std<Sym_tensor> Lorene::Time_slice::k_dd_evol
mutableprotected

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.

◆ k_uu_evol

Evolution_std<Sym_tensor> Lorene::Time_slice::k_uu_evol
mutableprotected

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.

◆ n_evol

Evolution_std<Scalar> Lorene::Time_slice::n_evol
mutableprotected

Values at successive time steps of the lapse function N.

Definition at line 219 of file time_slice.h.

◆ p_gamma

Metric* Lorene::Time_slice::p_gamma
mutableprotected

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

Definition at line 242 of file time_slice.h.

◆ scheme_order

int Lorene::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 190 of file time_slice.h.

◆ the_time

Evolution_std<double> Lorene::Time_slice::the_time
protected

Time label of each slice.

Definition at line 196 of file time_slice.h.

◆ trk_evol

Evolution_std<Scalar> Lorene::Time_slice::trk_evol
mutableprotected

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: