88#include "time_slice.h"
91#include "utilitaires.h"
96 double pdt,
double precis,
int method_poisson_vect,
97 const char* graph_device,
const Scalar* p_ener_dens,
98 const Vector* p_mom_dens,
const Scalar* p_trace_stress) {
107 "Time_slice_conf::initial_data_cts : the trace of u^{ij} with respect\n"
108 <<
" to the conformal metric tgam_{ij} is not zero !\n"
109 <<
" error = " << tr_uu << endl ;
133 int nz = map.get_mg()->get_nzone() ;
134 double ray_des = 1.25 * map.val_r(nz-2, 1., 0., 0.) ;
138 if (p_ener_dens != 0x0) ener_dens = *(p_ener_dens) ;
141 Vector mom_dens(map, CON, triad) ;
142 if (p_mom_dens != 0x0) mom_dens = *(p_mom_dens) ;
145 Scalar trace_stress(map) ;
146 if (p_trace_stress != 0x0) trace_stress = *(p_trace_stress) ;
152 Vector source_beta(map, CON, triad) ;
157 for (
int i=0; i<imax; i++) {
178 source_psi = tmp -
psi()*
psi4()* ( 0.5*qpig* ener_dens
180 - 8.33333333333333e-2*
trk()*
trk() ) ;
185 source_nn =
psi4()*(
nn()*( qpig* (ener_dens + trace_stress) + aa_quad
186 - 0.3333333333333333*
trk()*
trk() )
203 dnn - 6.*
nn() * dln_psi, 0) ;
205 source_beta += 2.*
nn() * ( 2.*qpig*
psi4() * mom_dens
211 beta().derive_cov(
ff).derive_cov(
ff), 1, 2)
212 + 0.3333333333333333*
218 source_beta -= vtmp ;
236 "Absolute error in the resolution of the equation for Psi") ;
238 des_meridian(psi_jp1, 0., ray_des,
"Psi", ngraph0, graph_device) ;
249 "Absolute error in the resolution of the equation for N") ;
251 des_meridian(nn_jp1, 0., ray_des,
"N", ngraph0+1, graph_device) ;
257 method_poisson_vect) ;
259 des_meridian(beta_jp1(1), 0., ray_des,
"\\gb\\ur\\d", ngraph0+2,
261 des_meridian(beta_jp1(2), 0., ray_des,
"\\gb\\u\\gh\\d", ngraph0+3,
263 des_meridian(beta_jp1(3), 0., ray_des,
"\\gb\\u\\gf\\d", ngraph0+4,
270 maxabs(test_beta - source_beta,
271 "Absolute error in the resolution for beta") ;
281 cout <<
"step = " << i <<
" : diff_psi = " << diff_psi
282 <<
" diff_nn = " << diff_nn
283 <<
" diff_beta = " << diff_beta << endl ;
284 if ( (diff_psi < precis) && (diff_nn < precis) && (diff_beta < precis) )
316 double ttime1 = ttime ;
318 for (
int j=1; j <
depth; j++) {
330 the_time.update(ttime1, jtime1, ttime1) ;
339 for (
int j=1; j <
depth; j++) {
Vectorial bases (triads) with respect to which the tensorial components are defined.
virtual const Scalar & ricci_scal() const
Returns the Ricci scalar.
Tensor field of valence 0 (or component of a tensorial field).
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this.
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
const Scalar & laplacian(int ced_mult_r=4) const
Returns the Laplacian of *this.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values of the Scalar in the co...
bool check_dzpuis(int dzi) const
Returns false if the last domain is compactified and *this is not zero in this domain and dzpuis is n...
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Class intended to describe valence-2 symmetric tensors.
const Vector & divergence(const Metric &) const
Returns the divergence of this with respect to a Metric .
Evolution_std< Sym_tensor > hh_evol
Values at successive time steps of the components .
virtual const Sym_tensor & hh(Param *=0x0, Param *=0x0) const
Deviation of the conformal metric from the flat metric : .
virtual const Scalar & nn() const
Lapse function N at the current time step (jtime ).
virtual void set_psi_del_npsi(const Scalar &psi_in)
Sets the conformal factor relating the physical metric to the conformal one: .
const Scalar & ln_psi() const
Logarithm of at the current time step (jtime ).
Evolution_std< Sym_tensor > hata_evol
Values at successive time steps of the components .
virtual const Vector & hdirac() const
Vector which vanishes in Dirac gauge.
virtual const Scalar & A_hata() const
Returns the potential A of .
virtual const Scalar & psi() const
Conformal factor relating the physical metric to the conformal one: .
virtual const Scalar & B_hata() const
Returns the potential of .
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 approa...
virtual void del_deriv() const
Deletes all the derived quantities.
Evolution_std< Scalar > npsi_evol
Values at successive time steps of the factor .
virtual Sym_tensor aa() const
Conformal representation of the traceless part of the extrinsic curvature: .
virtual void set_hata(const Sym_tensor &hata_in)
Sets the conformal representation of the traceless part of the extrinsic curvature: .
const Scalar & psi4() const
Factor at the current time step (jtime ).
Evolution_std< Scalar > psi_evol
Values at successive time steps of the conformal factor relating the physical metric to the conform...
virtual const Metric & tgam() const
Conformal metric Returns the value at the current time step (jtime ).
Evolution_std< Scalar > A_hata_evol
Potential A associated with the symmetric tensor .
const Metric_flat & ff
Pointer on the flat metric with respect to which the conformal decomposition is performed.
Evolution_std< Scalar > B_hata_evol
Potential associated with the symmetric tensor .
virtual const Scalar & trk() const
Trace K of the extrinsic curvature at the current time step (jtime ).
int jtime
Time step index of the latest slice.
Evolution_std< Scalar > trk_evol
Values at successive time steps of the trace K of the extrinsic curvature.
Evolution_std< Sym_tensor > k_dd_evol
Values at successive time steps of the covariant components of the extrinsic curvature tensor .
virtual const Vector & beta() const
shift vector at the current time step (jtime )
int depth
Number of stored time slices.
Evolution_std< Scalar > n_evol
Values at successive time steps of the lapse function N.
Evolution_std< double > the_time
Time label of each slice.
Evolution_std< Sym_tensor > k_uu_evol
Values at successive time steps of the contravariant components of the extrinsic curvature tensor .
Evolution_std< Vector > beta_evol
Values at successive time steps of the shift vector .
Tensor field of valence 1.
Vector poisson(double lambda, int method=6) const
Solves the vector Poisson equation with *this as a source.
const Scalar & divergence(const Metric &) const
The divergence of this with respect to a Metric .
Sym_tensor ope_killing_conf(const Metric &gam) const
Computes the conformal Killing operator associated with a given metric.
Vector derive_lie(const Vector &v) const
Computes the Lie derivative of this with respect to some vector field v.
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
const Map & get_mp() const
Returns the mapping.
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
Tensor up_down(const Metric &gam) const
Computes a new tensor by raising or lowering all the indices of *this .
const Tensor & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified ext...
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
const Tensor & derive_cov(const Metric &gam) const
Returns the covariant derivative of this with respect to some metric .
Tensor trace(int ind1, int ind2) const
Trace on two different type indices.
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Tbl maxabs(const Tensor &aa, const char *comment=0x0, ostream &ost=cout, bool verb=true)
Maxima in each domain of the absolute values of the tensor components.
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Standard units of space, time and mass.