Tslice_dirac_max Class Reference
[Time evolution (under development)]

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 Tslice_dirac_max:
Time_slice_conf Time_slice

List of all members.

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.
 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).
 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.
 Tslice_dirac_max (const Star_rot_Dirac &star, double pdt, int depth_in=3)
 Construnction of a stationary slice from a rotating star.
 Tslice_dirac_max (const Tslice_dirac_max &)
 Copy constructor.
virtual ~Tslice_dirac_max ()
 Destructor.
void operator= (const Tslice_dirac_max &)
 Assignment to another Tslice_dirac_max.
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 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 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).
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).
virtual void set_trh (const Scalar &trh_in)
 Sets the trace, with respect to the flat metric ff , of $ h^{ij} $.
virtual Scalar solve_psi (const Scalar *ener_dens=0x0) const
 Solves the elliptic equation for the conformal factor $$ (Hamiltonian constraint).
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).
virtual Vector solve_beta (int method=6) const
 Solves the elliptic equation for the shift vector $\beta^i$ from $ A^{ij} $ (Eq.
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.
virtual double adm_mass () const
 Returns the ADM mass at (geometrical units) the current step.
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} $.
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 ScalarA_hh () const
 Returns the potential A of $ \bar{h}^{ij} $.
virtual const ScalarB_hh () const
 Returns the potential $\tilde{B}$ of $ \bar{h}^{ij} $.
virtual const Scalartrh () const
 Computes the trace h, with respect to the flat metric ff , of $ h^{ij} $.
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_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_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 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.
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

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} $.
void initialize_sources_copy () const
 Copy the sources source_A_XXX_evol and source_B_XXX_evol to all time-steps.
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} $.
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} $.
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.
virtual void del_deriv () const
 Deletes all the derived quantities.
void set_der_0x0 () const
 Sets to 0x0 all the pointers on derived quantities.

Protected Attributes

Evolution_std< ScalarA_hh_evol
 The A potential of $ \bar{h}^{ij} $.
Evolution_std< ScalarB_hh_evol
 The $\tilde{B} $ potential of $ \bar{h}^{ij} $.
Evolution_std< Scalarsource_A_hh_evol
 The A potential of the source of equation for $ \bar{h}^{ij} $.
Evolution_std< Scalarsource_B_hh_evol
 The $\tilde{B} $ potential of the source of equation for $ \bar{h}^{ij} $.
Evolution_std< Scalarsource_A_hata_evol
 The potential A of the source of equation for $ \hat{A}^{ij} $.
Evolution_std< Scalarsource_B_hata_evol
 The potential $\tilde{B}$ of the source of equation for $ \hat{A}^{ij} $.
Evolution_std< Scalartrh_evol
 The trace, with respect to the flat metric ff , of $ h^{ij} $.
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 in the maximal slicing and Dirac gauge.

()

Definition at line 964 of file time_slice.h.


Constructor & Destructor Documentation

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_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 $ \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_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 145 of file tslice_dirac_max.C.

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:
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 159 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Time_slice::jtime, Scalar::set_etat_zero(), source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, Time_slice::the_time, trh_evol, and Evolution_std< TyT >::update().

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:
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 Tslice_dirac_max
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 189 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Time_slice::depth, fread_be(), Map::get_mg(), Time_slice_conf::hh_evol, Time_slice::jtime, Time_slice::the_time, and Evolution_std< TyT >::update().

Tslice_dirac_max::Tslice_dirac_max ( const Star_rot_Dirac star,
double  pdt,
int  depth_in = 3 
)
Tslice_dirac_max::Tslice_dirac_max ( const Tslice_dirac_max tin  ) 

Copy constructor.

Definition at line 295 of file tslice_dirac_max.C.

Tslice_dirac_max::~Tslice_dirac_max (  )  [virtual]

Destructor.

Definition at line 310 of file tslice_dirac_max.C.


Member Function Documentation

const Scalar & Time_slice_conf::A_hata (  )  const [virtual, inherited]

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 660 of file time_slice_conf.C.

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

const Scalar & 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 521 of file tslice_dirac_max.C.

References A_hh_evol, Time_slice_conf::hh_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::the_time, and Evolution_std< TyT >::update().

Sym_tensor Time_slice_conf::aa (  )  const [virtual, inherited]

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 761 of file time_slice_conf.C.

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

double Tslice_dirac_max::adm_mass (  )  const [virtual]

Returns the ADM mass at (geometrical units) the current step.

Moreover this method updates adm_mass_evol if necessary.

Reimplemented from Time_slice.

Definition at line 139 of file tslice_adm_mass.C.

References Time_slice::adm_mass_evol, Scalar::derive_con(), Time_slice_conf::ff, Vector::flux(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Tbl::get_taille(), hh(), Time_slice::jtime, Time_slice_conf::psi(), Tbl::set(), Tbl::set_etat_qcq(), Time_slice::the_time, Tensor::trace(), Evolution_full< TyT >::update(), and Map::val_r().

const Scalar & Time_slice_conf::B_hata (  )  const [virtual, inherited]

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 674 of file time_slice_conf.C.

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

const Scalar & 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 532 of file tslice_dirac_max.C.

References B_hh_evol, Time_slice_conf::hh_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::the_time, and Evolution_std< TyT >::update().

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

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

Definition at line 83 of file time_slice_access.C.

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

Tbl Time_slice::check_dynamical_equations ( const Sym_tensor strain_tensor = 0x0,
const Scalar energy_density = 0x0,
ostream &  ost = cout,
bool  verb = true 
) const [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 135 of file tslice_check_einstein.C.

References Tensor::annule_domain(), Time_slice::beta(), Metric::con(), contract(), Tensor::derive_con(), Scalar::derive_con(), Sym_tensor::derive_lie(), Time_slice::gam(), Tbl::get_dim(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Time_slice::jtime, Time_slice::k_dd(), Time_slice::k_dd_evol, Time_slice::k_uu(), maxabs(), Time_slice::nn(), Tensor_sym::position(), Metric::ricci(), Time_slice::scheme_order, Itbl::set(), Evolution< TyT >::time_derive(), Tensor::trace(), Time_slice::trk(), and Tensor::up_down().

Tbl 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 75 of file tslice_check_einstein.C.

References contract(), Scalar::dec_dzpuis(), Time_slice::gam(), Time_slice::k_dd(), Time_slice::k_uu(), maxabs(), Metric::ricci_scal(), and Time_slice::trk().

Tbl 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 105 of file tslice_check_einstein.C.

References Scalar::derive_con(), Sym_tensor::divergence(), Time_slice::gam(), Tensor::get_index_type(), Time_slice::k_uu(), maxabs(), and Time_slice::trk().

void Tslice_dirac_max::compute_sources ( const Sym_tensor strain_tensor = 0x0  )  const [protected]
void 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 833 of file time_slice_conf.C.

References abs(), contract(), Tensor::get_index_type(), max(), and Vector::poisson().

void Time_slice_conf::del_deriv (  )  const [protected, virtual, inherited]
void 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:
pdt time step dt.
nb_time_steps number of time steps for the evolution
niter_elliptic number of iterations if the resolution of elliptic equations
relax_elliptic relaxation factor for the elliptic equations
check_mod determines the frequency of check of the constraint equations: they are checked every check_mod time step
save_mod determines the frequency of writing to file the monotoring quantities: they are written to file every save_mod time step
method method_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_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.

Definition at line 113 of file tslice_dirac_max_evolve.C.

References Time_slice_conf::A_hata_evol, A_hh(), A_hh_evol, Param::add_double_mod(), Param::add_int(), Param::add_tbl_mod(), adm_mass(), Tensor::annule(), Tbl::annule_hard(), Scalar::annule_hard(), arrete(), Time_slice_conf::B_hata_evol, B_hh(), B_hh_evol, Time_slice::beta(), Time_slice::beta_evol, central_value(), Time_slice::check_dynamical_equations(), Time_slice::check_hamiltonian_constraint(), Time_slice::check_momentum_constraint(), Param::clean_all(), compute_sources(), Time_slice_conf::compute_X_from_momentum_constraint(), Time_slice_conf::del_deriv(), Time_slice::depth, des_evol(), des_meridian(), Time_slice_conf::ff, Map::get_bvect_spher(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Tensor::get_triad(), Time_slice_conf::hata(), Time_slice_conf::hata_evol, hh(), hh_det_one(), Time_slice_conf::hh_evol, Tensor::inc_dzpuis(), Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice_conf::k_dd(), log(), maxabs_all_domains(), Sym_tensor::mu(), Base_val::mult_cost(), Scalar::mult_r(), Base_val::mult_x(), Time_slice_conf::nn(), Time_slice_conf::npsi(), Time_slice_conf::npsi_evol, Time_slice_conf::psi(), Time_slice_conf::psi_evol, Time_slice::save(), Evolution< TyT >::save(), Tbl::set(), Tensor::set(), Sym_tensor_tt::set_A_tildeB(), set_AB_hh(), Scalar::set_domain(), Scalar::set_dzpuis(), Tbl::set_etat_qcq(), Time_slice_conf::set_npsi_del_n(), Time_slice_conf::set_psi_del_npsi(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), solve_beta(), solve_npsi(), solve_psi(), source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, Scalar::std_spectral_base(), Time_slice::the_time, Evolution< TyT >::time_derive(), trh(), Evolution_full< TyT >::update(), Evolution_std< TyT >::update(), Map::val_r(), Time_slice_conf::vec_X(), Valeur::ylm(), and Valeur::ylm_i().

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

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

Definition at line 91 of file time_slice_access.C.

References Time_slice::gam_dd(), and Time_slice::p_gamma.

const Sym_tensor & Time_slice_conf::gam_dd (  )  const [virtual, inherited]

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

Reimplemented from Time_slice.

Definition at line 604 of file time_slice_conf.C.

References Time_slice::gam_dd_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice_conf::psi4(), Time_slice_conf::tgam(), Time_slice::the_time, and Evolution_std< TyT >::update().

const Sym_tensor & Time_slice_conf::gam_uu (  )  const [virtual, inherited]

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

Reimplemented from Time_slice.

Definition at line 615 of file time_slice_conf.C.

References Time_slice::gam_uu_evol, Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice_conf::psi4(), Time_slice_conf::tgam(), Time_slice::the_time, and Evolution_std< TyT >::update().

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

Gets the latest value of time step index.

Definition at line 339 of file time_slice.h.

References Time_slice::jtime.

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

Gets the order of the finite-differences scheme.

Definition at line 336 of file time_slice.h.

References Time_slice::scheme_order.

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

Gets the time coordinate t at successive time steps.

Definition at line 342 of file time_slice.h.

References Time_slice::the_time.

const Sym_tensor & Time_slice_conf::hata (  )  const [virtual, inherited]
const Vector & 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 Time_slice_conf.

Definition at line 503 of file tslice_dirac_max.C.

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

const Sym_tensor & 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 Time_slice_conf.

Definition at line 470 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, hh_det_one(), Time_slice_conf::hh_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

void Tslice_dirac_max::hh_det_one ( const Sym_tensor_tt hijtt,
Param par_mat = 0x0 
) const [protected]
void Tslice_dirac_max::hh_det_one ( int  j,
Param par_bc = 0x0,
Param par_mat = 0x0 
) const [protected]
void 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:
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.

Definition at line 348 of file tslice_dirac_max.C.

References Time_slice_conf::A_hata_evol, A_hh_evol, Tensor::annule_domain(), Time_slice_conf::B_hata_evol, B_hh_evol, compute_sources(), Time_slice_conf::del_deriv(), Time_slice::depth, Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Mg3d::get_type_r(), Time_slice_conf::hh_evol, initialize_sources_copy(), Time_slice::jtime, maxabs(), Scalar::set_dzpuis(), Time_slice::the_time, and Evolution_std< TyT >::update().

void Tslice_dirac_max::initialize_sources_copy (  )  const [protected]
const Sym_tensor & Time_slice_conf::k_dd (  )  const [virtual, inherited]

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

Reimplemented from Time_slice.

Definition at line 626 of file time_slice_conf.C.

References Time_slice::gam(), Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::k_dd_evol, Time_slice_conf::k_uu(), Time_slice::the_time, and Evolution_std< TyT >::update().

const Sym_tensor & Time_slice_conf::k_uu (  )  const [virtual, inherited]

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

Reimplemented from Time_slice.

Definition at line 639 of file time_slice_conf.C.

References Time_slice::gam(), Time_slice_conf::hata(), Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::k_uu_evol, Time_slice_conf::psi(), Time_slice_conf::psi4(), Time_slice::the_time, Time_slice_conf::trk(), and Evolution_std< TyT >::update().

const Scalar & Time_slice_conf::ln_psi (  )  const [inherited]

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

Definition at line 715 of file time_slice_conf.C.

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

const Scalar & Time_slice_conf::nn (  )  const [virtual, inherited]

Lapse function N at the current time step (jtime ).

Reimplemented from Time_slice.

Definition at line 587 of file time_slice_conf.C.

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

const Scalar & Time_slice_conf::npsi (  )  const [virtual, inherited]
void Tslice_dirac_max::operator= ( const Tslice_dirac_max tin  ) 

Assignment to another Tslice_dirac_max.

Reimplemented from Time_slice_conf.

Definition at line 317 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, and trh_evol.

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

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

Reimplemented from Time_slice.

Definition at line 562 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Evolution< TyT >::is_known(), Time_slice::jtime, maxabs(), and trh_evol.

const Scalar & Time_slice_conf::psi (  )  const [virtual, inherited]

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 689 of file time_slice_conf.C.

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

const Scalar & Time_slice_conf::psi4 (  )  const [inherited]

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

Definition at line 703 of file time_slice_conf.C.

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

void Tslice_dirac_max::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 Time_slice.

Definition at line 583 of file tslice_dirac_max.C.

References A_hh(), A_hh_evol, B_hh(), B_hh_evol, Time_slice::depth, fwrite_be(), Time_slice_conf::hh_evol, Evolution< TyT >::is_known(), and Time_slice::jtime.

void 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 457 of file time_slice.C.

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

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

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

$ h^{ij} $ is not modified. Sets the value at the current time step (jtime ).

Definition at line 80 of file tslice_dirac_max_setAB.C.

References A_hh_evol, Time_slice::adm_mass_evol, B_hh_evol, Evolution< TyT >::downdate(), Time_slice::gam_dd_evol, Time_slice::gam_uu_evol, Time_slice_conf::hh_evol, Time_slice::jtime, Time_slice::p_gamma, Time_slice_conf::p_hdirac, Time_slice_conf::p_tgamma, Time_slice::the_time, trh_evol, and Evolution_std< TyT >::update().

void Time_slice_conf::set_der_0x0 (  )  const [protected, inherited]

Sets to 0x0 all the pointers on derived quantities.

Reimplemented from Time_slice.

Definition at line 357 of file time_slice_conf.C.

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

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

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 530 of file time_slice_conf.C.

References Time_slice_conf::A_hata_evol, Time_slice_conf::B_hata_evol, Evolution< TyT >::downdate(), Time_slice_conf::hata_evol, Time_slice::jtime, Time_slice::k_dd_evol, Time_slice::k_uu_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

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

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 562 of file time_slice_conf.C.

References Time_slice_conf::A_hata_evol, Time_slice_conf::B_hata_evol, Evolution< TyT >::downdate(), Time_slice_conf::ff, Tensor::get_mp(), Time_slice_conf::hata_evol, Map::inc_dzpuis(), Evolution< TyT >::is_known(), Time_slice::jtime, Time_slice::k_dd_evol, Time_slice::k_uu_evol, Vector::ope_killing_conf(), Time_slice_conf::p_vec_X, Time_slice::the_time, and Evolution_std< TyT >::update().

void Time_slice_conf::set_hata_TT ( const Sym_tensor_tt hata_tt  )  [virtual, inherited]
void 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 Time_slice_conf.

Definition at line 332 of file tslice_dirac_max.C.

References A_hh_evol, B_hh_evol, Evolution< TyT >::downdate(), Time_slice::jtime, source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, and trh_evol.

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

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

The value of $ h^{ij} $ is then deduced from the unimodularity condition on the conformal metric. Sets the value at the current time step (jtime ).

Definition at line 404 of file tslice_dirac_max.C.

References Scalar::dec_dzpuis(), Time_slice_conf::ff, Map::flat_met_spher(), Map::get_bvect_spher(), Tensor::get_mp(), Metric_flat::get_triad(), Time_slice_conf::hh_evol, Time_slice::jtime, Tensor::set_etat_zero(), Sym_tensor::set_longit_trans(), Time_slice::the_time, trh_evol, and Evolution_std< TyT >::update().

void Time_slice_conf::set_npsi_del_n ( const Scalar npsi_in  )  [virtual, inherited]

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

Definition at line 474 of file time_slice_conf.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::jtime, Time_slice::n_evol, Time_slice_conf::npsi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

void Time_slice_conf::set_npsi_del_psi ( const Scalar npsi_in  )  [virtual, inherited]
void Time_slice_conf::set_psi_del_n ( const Scalar psi_in  )  [virtual, inherited]

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 424 of file time_slice_conf.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::gam_dd_evol, Time_slice::gam_uu_evol, Time_slice::jtime, Time_slice::n_evol, Time_slice::p_gamma, Time_slice_conf::p_ln_psi, Time_slice_conf::p_psi4, Time_slice_conf::psi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

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

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 400 of file time_slice_conf.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::gam_dd_evol, Time_slice::gam_uu_evol, Time_slice::jtime, Time_slice_conf::npsi_evol, Time_slice::p_gamma, Time_slice_conf::p_ln_psi, Time_slice_conf::p_psi4, Time_slice_conf::psi_evol, Time_slice::the_time, and Evolution_std< TyT >::update().

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

Sets the order of the finite-differences scheme.

Definition at line 327 of file time_slice.h.

References Time_slice::scheme_order.

void 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 433 of file tslice_dirac_max.C.

References Time_slice::adm_mass_evol, Evolution< TyT >::downdate(), Time_slice::gam_dd_evol, Time_slice::gam_uu_evol, Time_slice_conf::hh_evol, Time_slice::jtime, Time_slice::p_gamma, Time_slice_conf::p_hdirac, Time_slice_conf::p_tgamma, source_A_hata_evol, source_A_hh_evol, source_B_hata_evol, source_B_hh_evol, Time_slice::the_time, trh_evol, and Evolution_std< TyT >::update().

Vector 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:
method method 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 220 of file tslice_dirac_max_solve.C.

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

Scalar 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_dens conformal 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_stress trace 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 114 of file tslice_dirac_max_solve.C.

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

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

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

Parameters:
ener_dens conformal 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 168 of file tslice_dirac_max_solve.C.

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

const Metric & Time_slice_conf::tgam (  )  const [virtual, inherited]

Conformal metric $ \tilde\gamma_{ij} = \Psi^{-4} \gamma_{ij} $ Returns the value at the current time step (jtime ).

Reimplemented in Isol_hor.

Definition at line 743 of file time_slice_conf.C.

References Metric_flat::con(), Time_slice_conf::ff, Time_slice_conf::hh(), and Time_slice_conf::p_tgamma.

const Scalar & 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 543 of file tslice_dirac_max.C.

References hh_det_one(), Evolution< TyT >::is_known(), Time_slice::jtime, and trh_evol.

const Scalar & 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 Time_slice_conf.

Definition at line 487 of file tslice_dirac_max.C.

References Time_slice_conf::ff, Metric::get_mp(), Evolution< TyT >::is_known(), Time_slice::jtime, Scalar::set_etat_zero(), Time_slice::the_time, Time_slice::trk_evol, and Evolution_std< TyT >::update().

const Vector & Time_slice_conf::vec_X ( int  method_poisson = 6  )  const [virtual, inherited]

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

(see the documentation of hata_evol)

Definition at line 821 of file time_slice_conf.C.

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


Friends And Related Function Documentation

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

Display.


Member Data Documentation

Evolution_std<Scalar> Time_slice_conf::A_hata_evol [mutable, protected, inherited]

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

(see documentation of Sym_tensor::p_aaa).

Definition at line 543 of file time_slice.h.

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

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

Definition at line 973 of file time_slice.h.

Evolution_full<Tbl> 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 229 of file time_slice.h.

Evolution_std<Scalar> Time_slice_conf::B_hata_evol [mutable, protected, inherited]

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

(see documentation of Sym_tensor::p_tilde_b).

Definition at line 548 of file time_slice.h.

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

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

Definition at line 979 of file time_slice.h.

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

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

Definition at line 215 of file time_slice.h.

int Time_slice::depth [protected, inherited]

Number of stored time slices.

Definition at line 175 of file time_slice.h.

const Metric_flat& Time_slice_conf::ff [protected, inherited]

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

Definition at line 503 of file time_slice.h.

Evolution_std<Sym_tensor> 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 194 of file time_slice.h.

Evolution_std<Sym_tensor> 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 199 of file time_slice.h.

Evolution_std<Sym_tensor> Time_slice_conf::hata_evol [mutable, protected, inherited]

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 538 of file time_slice.h.

Evolution_std<Sym_tensor> Time_slice_conf::hh_evol [mutable, protected, inherited]

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 526 of file time_slice.h.

int Time_slice::jtime [protected, inherited]

Time step index of the latest slice.

Definition at line 186 of file time_slice.h.

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

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

Definition at line 204 of file time_slice.h.

Evolution_std<Sym_tensor> 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 209 of file time_slice.h.

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

Values at successive time steps of the lapse function N.

Definition at line 212 of file time_slice.h.

Evolution_std<Scalar> Time_slice_conf::npsi_evol [mutable, protected, inherited]

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

Definition at line 518 of file time_slice.h.

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

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

Definition at line 235 of file time_slice.h.

Vector* Time_slice_conf::p_hdirac [mutable, protected, inherited]

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 567 of file time_slice.h.

Scalar* Time_slice_conf::p_ln_psi [mutable, protected, inherited]

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

Definition at line 562 of file time_slice.h.

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

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

Definition at line 559 of file time_slice.h.

Metric* Time_slice_conf::p_tgamma [mutable, protected, inherited]

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

Definition at line 556 of file time_slice.h.

Vector* Time_slice_conf::p_vec_X [mutable, protected, inherited]

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

(see the documentation of hata_evol)

Definition at line 573 of file time_slice.h.

Evolution_std<Scalar> Time_slice_conf::psi_evol [mutable, protected, inherited]

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 513 of file time_slice.h.

int 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 183 of file time_slice.h.

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 997 of file time_slice.h.

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 985 of file time_slice.h.

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 1003 of file time_slice.h.

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 991 of file time_slice.h.

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

Time label of each slice.

Definition at line 189 of file time_slice.h.

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

Definition at line 1006 of file time_slice.h.

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

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

Definition at line 220 of file time_slice.h.


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

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1