LORENE
Lorene::Tslice_dirac_max Class Reference

Spacelike time slice of a 3+1 spacetime with conformal decomposition in the maximal slicing and Dirac gauge. More...

#include <time_slice.h>

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

Public Member Functions

 Tslice_dirac_max (const Scalar &lapse_in, const Vector &shift_in, const Metric_flat &ff_in, const Scalar &psi_in, const Sym_tensor_trans &hh_in, const Sym_tensor &hata_in, int depth_in=3)
 Constructor from conformal decomposition. More...
 
 Tslice_dirac_max (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...
 
 Tslice_dirac_max (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...
 
 Tslice_dirac_max (const Star_rot_Dirac &star, double pdt, int depth_in=3)
 Construnction of a stationary slice from a rotating star. More...
 
 Tslice_dirac_max (const Tslice_dirac_max &)
 Copy constructor. More...
 
virtual ~Tslice_dirac_max ()
 Destructor. More...
 
void operator= (const Tslice_dirac_max &)
 Assignment to another Tslice_dirac_max. More...
 
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} $. 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 void set_khi_mu (const Scalar &khi_in, const Scalar &mu_in)
 Sets the potentials $\chi $ and $\mu$ of the TT part $ \bar{h}^{ij} $ of $ h^{ij} $ (see the documentation of Sym_tensor_tt for details). More...
 
virtual void set_AB_hh (const Scalar &A_in, const Scalar &B_in)
 Sets the potentials A and $\tilde{B}$ of the TT part $ \bar{h}^{ij} $ of $ h^{ij} $ (see the documentation of Sym_tensor for details). More...
 
virtual void set_trh (const Scalar &trh_in)
 Sets the trace, with respect to the flat metric ff , of $ h^{ij} $. More...
 
virtual Scalar solve_psi (const Scalar *ener_dens=0x0) const
 Solves the elliptic equation for the conformal factor $$ (Hamiltonian constraint). More...
 
virtual Scalar solve_npsi (const Scalar *ener_dens=0x0, const Scalar *trace_stress=0x0) const
 Solves the elliptic equation for $ N\Psi $ (maximal slicing condition + Hamiltonian constraint) More...
 
virtual Vector solve_beta (int method=6) const
 Solves the elliptic equation for the shift vector $\beta^i$ from $ A^{ij} $ (Eq. More...
 
void evolve (double pdt, int nb_time_steps, int niter_elliptic, double relax_elliptic, int check_mod, int save_mod, int method_poisson_vect=6, int nopause=1, const char *graph_device=0x0, bool verbose=true, const Scalar *ener_euler=0x0, const Vector *mom_euler=0x0, const Scalar *s_euler=0x0, const Sym_tensor *strain_euler=0x0)
 Time evolution by resolution of Einstein equations. More...
 
virtual double adm_mass () const
 Returns the ADM mass at (geometrical units) the current step. More...
 
virtual const Sym_tensorhh (Param *par_bc=0x0, Param *par_mat=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} $. More...
 
virtual const Scalartrk () const
 Trace K of the extrinsic curvature at the current time step (jtime ). More...
 
virtual const Vectorhdirac () const
 Vector $ H^i = {\cal D}_j \tilde\gamma^{ij} $ which vanishes in Dirac gauge. More...
 
virtual const ScalarA_hh () const
 Returns the potential A of $ \bar{h}^{ij} $. More...
 
virtual const ScalarB_hh () const
 Returns the potential $\tilde{B}$ of $ \bar{h}^{ij} $. More...
 
virtual const Scalartrh () const
 Computes the trace h, with respect to the flat metric ff , of $ h^{ij} $. More...
 
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} $. More...
 
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} $. More...
 
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$. More...
 
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. More...
 
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) $. More...
 
virtual void set_hata_TT (const Sym_tensor_tt &hata_tt)
 Sets the TT part of $ \hat{A}^{ij} $ (see member hata_evol ). More...
 
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 $. More...
 
virtual const Scalarnn () const
 Lapse function N 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 ScalarA_hata () const
 Returns the potential A of $ \hat{A}^{ij} $. More...
 
virtual const ScalarB_hata () const
 Returns the potential $\tilde{B}$ of $ \hat{A}^{ij} $. More...
 
virtual const Scalarpsi () const
 Conformal factor $ \Psi $ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $. More...
 
const Scalarpsi4 () const
 Factor $ \Psi^4 $ at the current time step (jtime ). More...
 
const Scalarln_psi () const
 Logarithm of $ \Psi $ at the current time step (jtime ). More...
 
virtual const Scalarnpsi () const
 Factor $ N\Psi $ at the current time step (jtime ). More...
 
virtual const Metrictgam () const
 Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $ Returns the value at the current time step (jtime ). More...
 
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) $. More...
 
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) $. More...
 
virtual const Vectorvec_X (int method_poisson=6) const
 Vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $. 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 $ X^i $ from the conformally-rescaled momentum $ \hat{S}^i = \Psi^6 S^i $, using the momentum constraint. More...
 
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). More...
 
void check_psi_dot (Tbl &tlnpsi_dot, Tbl &tdiff, Tbl &tdiff_rel) const
 Checks the $\frac{\partial}{\partial t} \ln\Psi $ 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 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...
 
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

void compute_sources (const Sym_tensor *strain_tensor=0x0) const
 Computes the sources source_A_XXX_evol and source_B_XXX_evol , for the solution of the evolution equation for $ h^{ij} $ and $ \hat{A}^{ij} $. More...
 
void initialize_sources_copy () const
 Copy the sources source_A_XXX_evol and source_B_XXX_evol to all time-steps. More...
 
void hh_det_one (int j, Param *par_bc=0x0, Param *par_mat=0x0) const
 Computes $ h^{ij} $ from the values of A and $\tilde{B}$ and using the condition $\det\tilde\gamma^{ij} = \det f^{ij} $, which fixes the trace of $ h^{ij} $. More...
 
void hh_det_one (const Sym_tensor_tt &hijtt, Param *par_mat=0x0) const
 Computes $ h^{ij} $ from the TT part using the condition $\det\tilde\gamma^{ij} = \det f^{ij} $, which fixes the trace of $ h^{ij} $. 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...
 
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...
 

Protected Attributes

Evolution_std< ScalarA_hh_evol
 The A potential of $ \bar{h}^{ij} $. More...
 
Evolution_std< ScalarB_hh_evol
 The $\tilde{B} $ potential of $ \bar{h}^{ij} $. More...
 
Evolution_std< Scalarsource_A_hh_evol
 The A potential of the source of equation for $ \bar{h}^{ij} $. More...
 
Evolution_std< Scalarsource_B_hh_evol
 The $\tilde{B} $ potential of the source of equation for $ \bar{h}^{ij} $. More...
 
Evolution_std< Scalarsource_A_hata_evol
 The potential A of the source of equation for $ \hat{A}^{ij} $. More...
 
Evolution_std< Scalarsource_B_hata_evol
 The potential $\tilde{B}$ of the source of equation for $ \hat{A}^{ij} $. More...
 
Evolution_std< Scalartrh_evol
 The trace, with respect to the flat metric ff , of $ h^{ij} $. More...
 
const Metric_flatff
 Pointer on the flat metric $ f_{ij} $ with respect to which the conformal decomposition is performed. More...
 
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} $. More...
 
Evolution_std< Scalarnpsi_evol
 Values at successive time steps of the factor $ N\Psi $. More...
 
Evolution_std< Sym_tensorhh_evol
 Values at successive time steps of the components $ h^{ij} $. More...
 
Evolution_std< Sym_tensorhata_evol
 Values at successive time steps of the components $ \hat{A}^{ij} $. More...
 
Evolution_std< ScalarA_hata_evol
 Potential A associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $. More...
 
Evolution_std< ScalarB_hata_evol
 Potential $ \tilde{B} $ associated with the symmetric tensor $ \hat{A}^{ij}_{TT} $. More...
 
Metricp_tgamma
 Pointer on the conformal metric $ \tilde\gamma_{ij} $ at the current time step (jtime) More...
 
Scalarp_psi4
 Pointer on the factor $ \Psi^4 $ at the current time step (jtime) More...
 
Scalarp_ln_psi
 Pointer on the logarithm of $ \Psi $ at the current time step (jtime) More...
 
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). More...
 
Vectorp_vec_X
 Pointer on the vector $ X^i $ representing the longitudinal part of $ \hat{A}^{ij} $. 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_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...
 

Detailed Description

Spacelike time slice of a 3+1 spacetime with conformal decomposition in the maximal slicing and Dirac gauge.

()

Definition at line 971 of file time_slice.h.

Constructor & Destructor Documentation

◆ Tslice_dirac_max() [1/5]

Lorene::Tslice_dirac_max::Tslice_dirac_max ( const Scalar lapse_in,
const Vector shift_in,
const Metric_flat ff_in,
const Scalar psi_in,
const Sym_tensor_trans hh_in,
const Sym_tensor hata_in,
int  depth_in = 3 
)

Constructor from conformal decomposition.

Parameters
lapse_inlapse function N
shift_inshift vector
ff_inreference flat metric with respect to which the conformal decomposition is performed
psi_inconformal factor $\Psi$ relating the physical metric $ \gamma_{ij} $ to the conformal one: $ \gamma_{ij} = \Psi^4 \tilde\gamma_{ij} $
hh_indeviation $ 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_inconformal 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) $, with K = 0 in the present case
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 152 of file tslice_dirac_max.C.

◆ Tslice_dirac_max() [2/5]

Lorene::Tslice_dirac_max::Tslice_dirac_max ( 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
mpMapping on which the various Lorene fields will be constructed
triadvector basis with respect to which the various tensor components will be defined
ff_inreference flat metric with respect to which the conformal decomposition is performed
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 166 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Lorene::Time_slice::jtime, Lorene::Scalar::set_etat_zero(), source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, Lorene::Time_slice::the_time, and trh_evol.

◆ Tslice_dirac_max() [3/5]

Lorene::Tslice_dirac_max::Tslice_dirac_max ( 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
mpMapping on which the various Lorene fields will be constructed
triadvector basis with respect to which the various tensor components will be defined
ff_inreference flat metric with respect to which the conformal decomposition is performed
fichfile containing the saved Tslice_dirac_max
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 must coincide with that stored in the file.

Definition at line 196 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Lorene::Time_slice::depth, Lorene::fread_be(), Lorene::Map::get_mg(), Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice::the_time.

◆ Tslice_dirac_max() [4/5]

◆ Tslice_dirac_max() [5/5]

Lorene::Tslice_dirac_max::Tslice_dirac_max ( const Tslice_dirac_max tin)

Copy constructor.

Definition at line 302 of file tslice_dirac_max.C.

◆ ~Tslice_dirac_max()

Lorene::Tslice_dirac_max::~Tslice_dirac_max ( )
virtual

Destructor.

Definition at line 317 of file tslice_dirac_max.C.

Member Function Documentation

◆ A_hata()

const Scalar & Lorene::Time_slice_conf::A_hata ( ) const
virtualinherited

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 Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice::the_time.

◆ A_hh()

const Scalar & Lorene::Tslice_dirac_max::A_hh ( ) const
virtual

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

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

Definition at line 528 of file tslice_dirac_max.C.

References A_hh_evol, Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice::the_time.

◆ aa()

Sym_tensor Lorene::Time_slice_conf::aa ( ) const
virtualinherited

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 Lorene::Time_slice_conf::hata(), Lorene::Time_slice_conf::psi(), and Lorene::Time_slice_conf::psi4().

◆ adm_mass()

◆ B_hata()

const Scalar & Lorene::Time_slice_conf::B_hata ( ) const
virtualinherited

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 Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice::the_time.

◆ B_hh()

const Scalar & Lorene::Tslice_dirac_max::B_hh ( ) const
virtual

Returns the potential $\tilde{B}$ of $ \bar{h}^{ij} $.

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

Definition at line 539 of file tslice_dirac_max.C.

References B_hh_evol, Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice::the_time.

◆ beta()

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

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, and Lorene::Time_slice::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
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::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().

◆ check_hamiltonian_constraint()

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

◆ check_momentum_constraint()

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

◆ check_psi_dot()

void Lorene::Time_slice_conf::check_psi_dot ( Tbl tlnpsi_dot,
Tbl tdiff,
Tbl tdiff_rel 
) const
inherited

◆ compute_sources()

void Lorene::Tslice_dirac_max::compute_sources ( const Sym_tensor strain_tensor = 0x0) const
protected

Computes the sources source_A_XXX_evol and source_B_XXX_evol , for the solution of the evolution equation for $ h^{ij} $ and $ \hat{A}^{ij} $.

Parameters
strain_tensor[input] : 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 (default value), it is assumed that $ S_{ij} $ = 0 (vacuum).

Definition at line 222 of file tslice_dirac_max_setAB.C.

References Lorene::Time_slice_conf::A_hata_evol, A_hh_evol, Lorene::Time_slice_conf::aa(), Lorene::Tensor::annule_domain(), Lorene::Time_slice_conf::B_hata_evol, B_hh_evol, Lorene::Time_slice::beta(), Lorene::Sym_tensor::compute_A(), Lorene::Sym_tensor::compute_tilde_B_tt(), Lorene::Metric::con(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Scalar::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::Sym_tensor::derive_lie(), Lorene::Vector::divergence(), Lorene::Time_slice_conf::ff, Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::get_triad(), Lorene::Time_slice_conf::hata(), hh(), Lorene::Tensor::inc_dzpuis(), Lorene::Time_slice::jtime, Lorene::Time_slice_conf::ln_psi(), Lorene::log(), Lorene::maxabs(), Lorene::Time_slice_conf::nn(), Lorene::Vector::ope_killing_conf(), Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi4(), Lorene::Tensor::set(), Lorene::Sym_tensor_tt::set_A_tildeB(), Lorene::Scalar::set_dzpuis(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::set_spectral_va(), source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, Lorene::Scalar::std_spectral_base(), Lorene::Time_slice_conf::tgam(), Lorene::Time_slice::the_time, Lorene::Tensor::trace(), and Lorene::Valeur::ylm().

◆ compute_X_from_momentum_constraint()

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

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

◆ del_deriv()

void Lorene::Time_slice_conf::del_deriv ( ) const
protectedvirtualinherited

◆ evolve()

void Lorene::Tslice_dirac_max::evolve ( double  pdt,
int  nb_time_steps,
int  niter_elliptic,
double  relax_elliptic,
int  check_mod,
int  save_mod,
int  method_poisson_vect = 6,
int  nopause = 1,
const char *  graph_device = 0x0,
bool  verbose = true,
const Scalar ener_euler = 0x0,
const Vector mom_euler = 0x0,
const Scalar s_euler = 0x0,
const Sym_tensor strain_euler = 0x0 
)

Time evolution by resolution of Einstein equations.

Parameters
pdttime step dt.
nb_time_stepsnumber of time steps for the evolution
niter_ellipticnumber of iterations if the resolution of elliptic equations
relax_ellipticrelaxation factor for the elliptic equations
check_moddetermines the frequency of check of the constraint equations: they are checked every check_mod time step
save_moddetermines the frequency of writing to file the monotoring quantities: they are written to file every save_mod time step
methodmethod_poisson_vect to be used for solving vector Poisson equation (for the shift), see Vector::poisson(double, const Metric_flat&, int) const.
nopause= 1 if no pause between each time step, 0 otherwise
graph_devicename of type of graphical device: 0x0 (default value) will result in interactive choice; "/xwin" in X-Window display and "/n" in no output.

Definition at line 123 of file tslice_dirac_max_evolve.C.

References Lorene::Time_slice_conf::A_hata_evol, A_hh_evol, Lorene::Time_slice_conf::B_hata_evol, B_hh_evol, Lorene::Time_slice::beta(), Lorene::Time_slice::depth, Lorene::Time_slice_conf::ff, Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_base(), Lorene::Tensor::get_triad(), Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, Lorene::Base_val::mult_cost(), Lorene::Base_val::mult_x(), Lorene::Time_slice_conf::nn(), Lorene::Sym_tensor_tt::set_A_tildeB(), and Lorene::Scalar::std_spectral_base().

◆ gam()

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.

◆ gam_dd()

const Sym_tensor & Lorene::Time_slice_conf::gam_dd ( ) const
virtualinherited

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::Time_slice::jtime, Lorene::Time_slice_conf::psi4(), Lorene::Time_slice_conf::tgam(), and Lorene::Time_slice::the_time.

◆ gam_uu()

const Sym_tensor & Lorene::Time_slice_conf::gam_uu ( ) const
virtualinherited

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::Time_slice::jtime, Lorene::Time_slice_conf::psi4(), Lorene::Time_slice_conf::tgam(), and Lorene::Time_slice::the_time.

◆ get_latest_j()

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

Gets the latest value of time step index.

Definition at line 346 of file time_slice.h.

References Lorene::Time_slice::jtime.

◆ get_scheme_order()

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

Gets the order of the finite-differences scheme.

Definition at line 343 of file time_slice.h.

References Lorene::Time_slice::scheme_order.

◆ get_time()

const Evolution_std<double>& Lorene::Time_slice::get_time ( ) const
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.

◆ hata()

◆ hdirac()

const Vector & Lorene::Tslice_dirac_max::hdirac ( ) const
virtual

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

It is null in the present case...

Reimplemented from Lorene::Time_slice_conf.

Definition at line 510 of file tslice_dirac_max.C.

References Lorene::Time_slice_conf::ff, Lorene::Metric::get_mp(), Lorene::Metric_flat::get_triad(), Lorene::Time_slice_conf::p_hdirac, and Lorene::Tensor::set_etat_zero().

◆ hh()

const Sym_tensor & Lorene::Tslice_dirac_max::hh ( Param par_bc = 0x0,
Param par_mat = 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 from Lorene::Time_slice_conf.

Definition at line 477 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, hh_det_one(), Lorene::Time_slice_conf::hh_evol, and Lorene::Time_slice::jtime.

◆ hh_det_one() [1/2]

◆ hh_det_one() [2/2]

◆ initial_data_cts()

void Lorene::Tslice_dirac_max::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
uuvalue 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_invalue of $ K = K_i^{\ i} $ (freely specifiable data)
trk_pointvalue of $ \partial K / \partial t $ (freely specifiable data)
pdttime step, to be used in order to fill depth slices
precisconvergence threshold required to stop the iteration
method_poisson_vectmethod to be used for solving vector Poisson equation (for the shift), see Vector::poisson(double, const Metric_flat&, int) const.
graph_devicename of type of graphical device: 0x0 (default value) will result in interactive choice; "/xwin" in X-Window display and "/n" in no output.
ener_densmatter 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_densmatter 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_stresstrace 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 from Lorene::Time_slice_conf.

Definition at line 355 of file tslice_dirac_max.C.

References Lorene::Time_slice::depth, Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice_conf::initial_data_cts(), and Lorene::Time_slice::jtime.

◆ initialize_sources_copy()

void Lorene::Tslice_dirac_max::initialize_sources_copy ( ) const
protected

Copy the sources source_A_XXX_evol and source_B_XXX_evol to all time-steps.

Definition at line 414 of file tslice_dirac_max_setAB.C.

References Lorene::Time_slice::depth, Lorene::Time_slice::jtime, source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, and Lorene::Time_slice::the_time.

◆ k_dd()

const Sym_tensor & Lorene::Time_slice_conf::k_dd ( ) const
virtualinherited

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::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice_conf::k_uu(), and Lorene::Time_slice::the_time.

◆ k_uu()

const Sym_tensor & Lorene::Time_slice_conf::k_uu ( ) const
virtualinherited

◆ ln_psi()

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

Logarithm of $ \Psi $ at the current time step (jtime ).

Definition at line 722 of file time_slice_conf.C.

References Lorene::log(), Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::psi(), and Lorene::Scalar::std_spectral_base().

◆ nn()

const Scalar & Lorene::Time_slice_conf::nn ( ) const
virtualinherited

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, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice_conf::psi_evol, and Lorene::Time_slice::the_time.

◆ npsi()

const Scalar & Lorene::Time_slice_conf::npsi ( ) const
virtualinherited

◆ operator=()

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

◆ operator>>()

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

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

Reimplemented from Lorene::Time_slice_conf.

Definition at line 569 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Lorene::Time_slice::jtime, Lorene::maxabs(), Lorene::Time_slice_conf::operator>>(), and trh_evol.

◆ psi()

const Scalar & Lorene::Time_slice_conf::psi ( ) const
virtualinherited

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::Time_slice::jtime, Lorene::Time_slice::n_evol, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice_conf::psi_evol, and Lorene::Time_slice::the_time.

◆ psi4()

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

Factor $ \Psi^4 $ at the current time step (jtime ).

Definition at line 710 of file time_slice_conf.C.

References Lorene::Time_slice_conf::p_psi4, Lorene::pow(), Lorene::Time_slice_conf::psi(), and Lorene::Scalar::std_spectral_base().

◆ sauve()

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

Definition at line 590 of file tslice_dirac_max.C.

References A_hh(), A_hh_evol, B_hh(), B_hh_evol, Lorene::Time_slice::depth, Lorene::fwrite_be(), Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice_conf::sauve().

◆ save()

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

◆ set_AB_hata()

void Lorene::Time_slice_conf::set_AB_hata ( const Scalar A_in,
const Scalar B_in 
)
virtualinherited

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 Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, and Lorene::Time_slice::the_time.

◆ set_AB_hh()

void Lorene::Tslice_dirac_max::set_AB_hh ( const Scalar A_in,
const Scalar B_in 
)
virtual

◆ set_der_0x0()

void Lorene::Time_slice_conf::set_der_0x0 ( ) const
protectedinherited

◆ set_hata()

void Lorene::Time_slice_conf::set_hata ( const Sym_tensor hata_in)
virtualinherited

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 Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, and Lorene::Time_slice::the_time.

◆ set_hata_from_XAB()

void Lorene::Time_slice_conf::set_hata_from_XAB ( Param par_bc = 0x0,
Param par_mat = 0x0 
)
virtualinherited

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 Lorene::Time_slice_conf::A_hata_evol, Lorene::Time_slice_conf::B_hata_evol, Lorene::Time_slice_conf::ff, Lorene::Tensor::get_mp(), Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::k_dd_evol, Lorene::Time_slice::k_uu_evol, Lorene::Vector::ope_killing_conf(), Lorene::Time_slice_conf::p_vec_X, Lorene::Sym_tensor_tt::set_A_tildeB(), and Lorene::Time_slice::the_time.

◆ set_hata_TT()

void Lorene::Time_slice_conf::set_hata_TT ( const Sym_tensor_tt hata_tt)
virtualinherited

◆ set_hh()

void Lorene::Tslice_dirac_max::set_hh ( const Sym_tensor hh_in)
virtual

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} $.

$ h^{ij} $ must be such that $\det\tilde\gamma^{ij} = f^{-1} $. Sets the value at the current time step (jtime ).

Reimplemented from Lorene::Time_slice_conf.

Definition at line 339 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Lorene::Time_slice::jtime, Lorene::Time_slice_conf::set_hh(), source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, and trh_evol.

◆ set_khi_mu()

void Lorene::Tslice_dirac_max::set_khi_mu ( const Scalar khi_in,
const Scalar mu_in 
)
virtual

◆ set_npsi_del_n()

void Lorene::Time_slice_conf::set_npsi_del_n ( const Scalar npsi_in)
virtualinherited

Sets the factor $ N\Psi $ 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, Lorene::Time_slice_conf::npsi_evol, and Lorene::Time_slice::the_time.

◆ set_npsi_del_psi()

void Lorene::Time_slice_conf::set_npsi_del_psi ( const Scalar npsi_in)
virtualinherited

◆ set_psi_del_n()

void Lorene::Time_slice_conf::set_psi_del_n ( const Scalar psi_in)
virtualinherited

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::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::n_evol, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::psi_evol, and Lorene::Time_slice::the_time.

◆ set_psi_del_npsi()

void Lorene::Time_slice_conf::set_psi_del_npsi ( const Scalar psi_in)
virtualinherited

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::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice::jtime, Lorene::Time_slice_conf::npsi_evol, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_ln_psi, Lorene::Time_slice_conf::p_psi4, Lorene::Time_slice_conf::psi_evol, and Lorene::Time_slice::the_time.

◆ set_scheme_order()

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

Sets the order of the finite-differences scheme.

Definition at line 334 of file time_slice.h.

References Lorene::Time_slice::scheme_order.

◆ set_trh()

void Lorene::Tslice_dirac_max::set_trh ( const Scalar trh_in)
virtual

Sets the trace, with respect to the flat metric ff , of $ h^{ij} $.

Sets the value at the current time step (jtime ). Note that this method does not ensure that the conformal metric is unimodular.

Definition at line 440 of file tslice_dirac_max.C.

References Lorene::Time_slice::adm_mass_evol, Lorene::Time_slice::gam_dd_evol, Lorene::Time_slice::gam_uu_evol, Lorene::Time_slice_conf::hh_evol, Lorene::Time_slice::jtime, Lorene::Time_slice::p_gamma, Lorene::Time_slice_conf::p_hdirac, Lorene::Time_slice_conf::p_tgamma, source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, Lorene::Time_slice::the_time, and trh_evol.

◆ solve_beta()

Vector Lorene::Tslice_dirac_max::solve_beta ( int  method = 6) const
virtual

Solves the elliptic equation for the shift vector $\beta^i$ from $ A^{ij} $ (Eq.

(73) of Bonazzola et al. 2004).

Parameters
methodmethod to be used for solving vector Poisson equation (for the shift), see Vector::poisson(double, const Metric_flat&, int) const.
Returns
solution $\beta^i_{\rm new}$ of the elliptic equation (flat vector Laplacian) for the shift with the source computed from the quantities at the current time step.

Definition at line 227 of file tslice_dirac_max_solve.C.

References Lorene::Time_slice_conf::aa(), Lorene::Time_slice::beta(), Lorene::contract(), Lorene::Tensor::derive_con(), Lorene::Sym_tensor::divergence(), Lorene::Vector::divergence(), Lorene::Time_slice_conf::ff, hh(), Lorene::Tensor::inc_dzpuis(), Lorene::maxabs(), Lorene::Time_slice_conf::nn(), and Lorene::Vector::poisson().

◆ solve_npsi()

Scalar Lorene::Tslice_dirac_max::solve_npsi ( const Scalar ener_dens = 0x0,
const Scalar trace_stress = 0x0 
) const
virtual

Solves the elliptic equation for $ N\Psi $ (maximal slicing condition + Hamiltonian constraint)

Parameters
ener_densconformal matter energy density $ \hat{E} = \Psi^6 E $, where E is measured by the Eulerian observer; this quantity is passed as a pointer, the null value of which (default) meaning E=0.
trace_stresstrace of the conformal matter stress $ S^* = \Psi^6 S $, where S is measured by the Eulerian observer; this quantity is passed as a pointer, the null value of which (default) meaning S=0.
Returns
solution $(N\Psi)_{\rm new}$ of the elliptic equation (flat Laplacian) for $ N\Psi $ with the source computed from the quantities at the current time step.

Definition at line 121 of file tslice_dirac_max_solve.C.

References Lorene::Time_slice_conf::aa(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), hh(), Lorene::Scalar::inc_dzpuis(), Lorene::Time_slice_conf::npsi(), Lorene::Scalar::poisson(), Lorene::Time_slice_conf::psi(), Lorene::Time_slice_conf::psi4(), Lorene::Scalar::set_etat_zero(), Lorene::Time_slice_conf::tgam(), and Lorene::Tensor::up_down().

◆ solve_psi()

Scalar Lorene::Tslice_dirac_max::solve_psi ( const Scalar ener_dens = 0x0) const
virtual

Solves the elliptic equation for the conformal factor $$ (Hamiltonian constraint).

Parameters
ener_densconformal matter energy density $ \hat{E} = \Psi^6 E $, where E is measured by the Eulerian observer; this quantity is passed as a pointer, the null value of which (default) meaning E=0.
Returns
solution $\Psi_{\rm new}$ of the elliptic equation (flat Laplacian) for the lapse with the source computed from the quantities at the current time step.

Definition at line 175 of file tslice_dirac_max_solve.C.

References Lorene::contract(), Lorene::Metric::cov(), Lorene::Tensor::derive_cov(), Lorene::Scalar::derive_cov(), Lorene::Tensor_sym::derive_cov(), Lorene::Time_slice_conf::ff, Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Time_slice_conf::hata(), hh(), Lorene::Scalar::inc_dzpuis(), Lorene::Scalar::poisson(), Lorene::pow(), Lorene::Time_slice_conf::psi(), Lorene::Scalar::set_etat_zero(), Lorene::Time_slice_conf::tgam(), and Lorene::Tensor::up_down().

◆ tgam()

const Metric & Lorene::Time_slice_conf::tgam ( ) const
virtualinherited

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(), Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::hh(), and Lorene::Time_slice_conf::p_tgamma.

◆ trh()

const Scalar & Lorene::Tslice_dirac_max::trh ( ) const
virtual

Computes the trace h, with respect to the flat metric ff , of $ h^{ij} $.

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

Definition at line 550 of file tslice_dirac_max.C.

References hh_det_one(), Lorene::Time_slice::jtime, and trh_evol.

◆ trk()

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

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

It is null in the present case (maximal slicing)

Reimplemented from Lorene::Time_slice_conf.

Definition at line 494 of file tslice_dirac_max.C.

References Lorene::Time_slice_conf::ff, Lorene::Metric::get_mp(), Lorene::Time_slice::jtime, Lorene::Scalar::set_etat_zero(), Lorene::Time_slice::the_time, and Lorene::Time_slice::trk_evol.

◆ vec_X()

const Vector & Lorene::Time_slice_conf::vec_X ( int  method_poisson = 6) const
virtualinherited

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 Lorene::Time_slice_conf::ff, Lorene::Time_slice_conf::hata_evol, Lorene::Time_slice::jtime, and Lorene::Time_slice_conf::p_vec_X.

Member Data Documentation

◆ A_hata_evol

Evolution_std<Scalar> Lorene::Time_slice_conf::A_hata_evol
mutableprotectedinherited

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.

◆ A_hh_evol

Evolution_std<Scalar> Lorene::Tslice_dirac_max::A_hh_evol
mutableprotected

The A potential of $ \bar{h}^{ij} $.

(see the documentation of Sym_tensor::p_aaa for details).

Definition at line 980 of file time_slice.h.

◆ adm_mass_evol

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

◆ B_hata_evol

Evolution_std<Scalar> Lorene::Time_slice_conf::B_hata_evol
mutableprotectedinherited

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.

◆ B_hh_evol

Evolution_std<Scalar> Lorene::Tslice_dirac_max::B_hh_evol
mutableprotected

The $\tilde{B} $ potential of $ \bar{h}^{ij} $.

(see the documentation of Sym_tensor::p_tilde_b for details).

Definition at line 986 of file time_slice.h.

◆ beta_evol

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

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
protectedinherited

Number of stored time slices.

Definition at line 182 of file time_slice.h.

◆ ff

const Metric_flat& Lorene::Time_slice_conf::ff
protectedinherited

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.

◆ gam_dd_evol

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

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
mutableprotectedinherited

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

Definition at line 206 of file time_slice.h.

◆ hata_evol

Evolution_std<Sym_tensor> Lorene::Time_slice_conf::hata_evol
mutableprotectedinherited

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.

◆ hh_evol

Evolution_std<Sym_tensor> Lorene::Time_slice_conf::hh_evol
mutableprotectedinherited

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.

◆ jtime

int Lorene::Time_slice::jtime
protectedinherited

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
mutableprotectedinherited

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
mutableprotectedinherited

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
mutableprotectedinherited

Values at successive time steps of the lapse function N.

Definition at line 219 of file time_slice.h.

◆ npsi_evol

Evolution_std<Scalar> Lorene::Time_slice_conf::npsi_evol
mutableprotectedinherited

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

Definition at line 525 of file time_slice.h.

◆ p_gamma

Metric* Lorene::Time_slice::p_gamma
mutableprotectedinherited

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

Definition at line 242 of file time_slice.h.

◆ p_hdirac

Vector* Lorene::Time_slice_conf::p_hdirac
mutableprotectedinherited

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.

◆ p_ln_psi

Scalar* Lorene::Time_slice_conf::p_ln_psi
mutableprotectedinherited

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

Definition at line 569 of file time_slice.h.

◆ p_psi4

Scalar* Lorene::Time_slice_conf::p_psi4
mutableprotectedinherited

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

Definition at line 566 of file time_slice.h.

◆ p_tgamma

Metric* Lorene::Time_slice_conf::p_tgamma
mutableprotectedinherited

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

Definition at line 563 of file time_slice.h.

◆ p_vec_X

Vector* Lorene::Time_slice_conf::p_vec_X
mutableprotectedinherited

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.

◆ psi_evol

Evolution_std<Scalar> Lorene::Time_slice_conf::psi_evol
mutableprotectedinherited

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.

◆ scheme_order

int Lorene::Time_slice::scheme_order
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.

◆ source_A_hata_evol

Evolution_std<Scalar> Lorene::Tslice_dirac_max::source_A_hata_evol
mutableprotected

The potential A of the source of equation for $ \hat{A}^{ij} $.

(see the documentation of Sym_tensor::p_aaa for details).

Definition at line 1004 of file time_slice.h.

◆ source_A_hh_evol

Evolution_std<Scalar> Lorene::Tslice_dirac_max::source_A_hh_evol
mutableprotected

The A potential of the source of equation for $ \bar{h}^{ij} $.

(see the documentation of Sym_tensor::p_aaa for details).

Definition at line 992 of file time_slice.h.

◆ source_B_hata_evol

Evolution_std<Scalar> Lorene::Tslice_dirac_max::source_B_hata_evol
mutableprotected

The potential $\tilde{B}$ of the source of equation for $ \hat{A}^{ij} $.

(see the documentation of Sym_tensor::p_tilde_b for details).

Definition at line 1010 of file time_slice.h.

◆ source_B_hh_evol

Evolution_std<Scalar> Lorene::Tslice_dirac_max::source_B_hh_evol
mutableprotected

The $\tilde{B} $ potential of the source of equation for $ \bar{h}^{ij} $.

(see the documentation of Sym_tensor::p_tilde_b for details).

Definition at line 998 of file time_slice.h.

◆ the_time

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

Time label of each slice.

Definition at line 196 of file time_slice.h.

◆ trh_evol

Evolution_std<Scalar> Lorene::Tslice_dirac_max::trh_evol
mutableprotected

The trace, with respect to the flat metric ff , of $ h^{ij} $.

Definition at line 1013 of file time_slice.h.

◆ trk_evol

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

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: