LORENE
time_slice.h
1/*
2 * Definition of Lorene class Time_slice
3 *
4 */
5
6/*
7 * Copyright (c) 2004 Eric Gourgoulhon, Jose Luis Jaramillo & Jerome Novak
8 *
9 * This file is part of LORENE.
10 *
11 * LORENE is free software; you can redistribute it and/or modify
12 * it under the terms of the GNU General Public License version 2
13 * as published by the Free Software Foundation.
14 *
15 * LORENE is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with LORENE; if not, write to the Free Software
22 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23 *
24 */
25
26#ifndef __TIME_SLICE_H_
27#define __TIME_SLICE_H_
28
29/*
30 * $Id: time_slice.h,v 1.33 2019/12/06 13:56:54 j_novak Exp $
31 * $Log: time_slice.h,v $
32 * Revision 1.33 2019/12/06 13:56:54 j_novak
33 * Corrected mistake in a comment
34 *
35 * Revision 1.32 2014/10/13 08:52:37 j_novak
36 * Lorene classes and functions now belong to the namespace Lorene.
37 *
38 * Revision 1.31 2012/02/06 12:59:07 j_novak
39 * Correction of some errors.
40 *
41 * Revision 1.30 2010/10/20 07:58:09 j_novak
42 * Better implementation of the explicit time-integration. Not fully-tested yet.
43 *
44 * Revision 1.29 2008/12/04 18:22:49 j_novak
45 * Enhancement of the dzpuis treatment + various bug fixes.
46 *
47 * Revision 1.28 2008/12/02 15:02:21 j_novak
48 * Implementation of the new constrained formalism, following Cordero et al. 2009
49 * paper. The evolution eqs. are solved as a first-order system. Not tested yet!
50 *
51 * Revision 1.27 2008/08/19 06:41:59 j_novak
52 * Minor modifications to avoid warnings with gcc 4.3. Most of them concern
53 * cast-type operations, and constant strings that must be defined as const char*
54 *
55 * Revision 1.26 2007/11/06 14:47:06 j_novak
56 * New constructor from a rotating star in Dirac gauge (class Star_rot_Dirac).
57 * Evolution can take into account matter terms.
58 *
59 * Revision 1.25 2007/04/25 15:20:59 j_novak
60 * Corrected an error in the initialization of tildeB in
61 * Tslice_dirac_max::initial_dat_cts. + New method for solve_hij_AB.
62 *
63 * Revision 1.24 2007/03/21 14:51:48 j_novak
64 * Introduction of potentials A and tilde(B) of h^{ij} into Tslice_dirac_max.
65 *
66 * Revision 1.23 2005/03/28 19:44:00 f_limousin
67 * Function tgam() is now virtual.
68 *
69 * Revision 1.22 2004/06/24 20:36:07 e_gourgoulhon
70 * Class Time_slice_conf: added method check_psi_dot.
71 *
72 * Revision 1.21 2004/06/14 20:46:35 e_gourgoulhon
73 * Added argument method_poisson to Tslice_dirac_max::solve_hij.
74 *
75 * Revision 1.20 2004/05/31 20:28:20 e_gourgoulhon
76 * -- Class Time_slice : added inline functions get_latest_j() and
77 * get_time()
78 * -- Class Tslice_dirac_max: method hh_det_one takes now a time step
79 * as argument, to compute h^{ij} from khi and mu at an arbitrary
80 * time step and not only the latest one.
81 *
82 * Revision 1.19 2004/05/27 15:22:28 e_gourgoulhon
83 * Added functions save and sauve, as well as constructors from file.
84 *
85 * Revision 1.18 2004/05/20 20:30:37 e_gourgoulhon
86 * Added arguments check_mod and save_mod to method Tsclice_dirac_max::evolve.
87 *
88 * Revision 1.17 2004/05/17 19:52:16 e_gourgoulhon
89 * -- Method initial_data_cts: added arguments graph_device and
90 * method_poisson_vect.
91 * -- Method Tslice_dirac_max::solve_beta : added argument method
92 * -- Method Tslice_dirac_max::solve_hij : added argument graph_device
93 * -- Method Tslice_dirac_max::evolve : added arguments
94 * method_poisson_vect, nopause and graph_device.
95 *
96 * Revision 1.16 2004/05/12 15:16:25 e_gourgoulhon
97 * Added #include "metric.h" before #include "evolution.h".
98 *
99 * Revision 1.15 2004/05/09 20:56:09 e_gourgoulhon
100 * Added member adm_mass_evol and corresponding virtual method adm_mass().
101 *
102 * Revision 1.14 2004/05/06 15:23:10 e_gourgoulhon
103 * initial_data_cts is know a virtual function of class Time_slice_conf
104 * and is implemented also for class Tslice_dirac_max.
105 *
106 * Revision 1.13 2004/05/03 14:46:11 e_gourgoulhon
107 * Class Tslice_dirac_max: -- changed prototype of method solve_hij
108 * -- added new method evolve
109 *
110 * Revision 1.12 2004/04/30 14:36:15 j_novak
111 * Added the method Tslice_dirac_max::solve_hij(...)
112 * NOT READY YET!!!
113 *
114 * Revision 1.11 2004/04/30 10:51:38 e_gourgoulhon
115 * Class Tslice_dirac_max: added methods solve_n, solve_q and solve_beta
116 * for resolution of the elliptic part of Einstein equations.
117 *
118 * Revision 1.10 2004/04/29 17:07:27 e_gourgoulhon
119 * Added argument pdt to Time_slice_conf::initial_data_cts.
120 *
121 * Revision 1.9 2004/04/08 16:42:11 e_gourgoulhon
122 * Many changes:
123 * -- class Time_slice_conf: added methods set_*, changed argument list of
124 * method initial_data_cts.
125 * -- class Tslice_dirac_max: added methods set_* and hh_det_one().
126 *
127 * Revision 1.8 2004/04/05 21:21:51 e_gourgoulhon
128 * class Time_slice_conf: added method initial_data_cts (computation of
129 * initial data from conformally thin sandwich method).
130 * classes Time_slice_conf and Tslice_dirac_max: added constructor as
131 * standard time slice of Minkowski spacetime.
132 *
133 * Revision 1.7 2004/04/01 16:09:01 j_novak
134 * Trace of K_ij is now member of Time_slice (it was member of Time_slice_conf).
135 * Added new methods for checking 3+1 Einstein equations (preliminary).
136 *
137 * Revision 1.6 2004/03/30 14:00:30 j_novak
138 * New class Tslide_dirac_max (first version).
139 *
140 * Revision 1.5 2004/03/29 11:58:53 e_gourgoulhon
141 * Many modif. to class Time_slice_conf.
142 * Minor modif. to class Time_slice.
143 *
144 * Revision 1.4 2004/03/28 21:33:14 e_gourgoulhon
145 * Constructor Time_slice::Time_slice(int depth_in) declared "explicit".
146 *
147 * Revision 1.3 2004/03/28 21:27:57 e_gourgoulhon
148 * Class Time_slice: - renamed the Evolution_std with suffix "_evol".
149 * - added protected constructor for derived classes
150 * Added class Time_slice_conf.
151 *
152 * Revision 1.2 2004/03/26 13:33:02 j_novak
153 * New methods for accessing/updating members (nn(), beta(), gam_uu(), k_uu(), ...)
154 *
155 * Revision 1.1 2004/03/24 14:56:18 e_gourgoulhon
156 * First version
157 *
158 *
159 * $Header: /cvsroot/Lorene/C++/Include/time_slice.h,v 1.33 2019/12/06 13:56:54 j_novak Exp $
160 *
161 */
162
163#include "star_rot_dirac.h"
164#include "evolution.h"
165
166 //---------------------------//
167 // class Time_slice //
168 //---------------------------//
169
170namespace Lorene {
177
178 // Data :
179 // -----
180 protected:
182 int depth ;
183
191
193 int jtime ;
194
197
202
207
212
217
220
223
228
237
238 // Derived data :
239 // ------------
240 protected:
242 mutable Metric* p_gamma ;
243
244 // Constructors - Destructor
245 // -------------------------
246 public:
247
259 Time_slice(const Scalar& lapse_in, const Vector& shift_in,
260 const Sym_tensor& gamma_in, const Sym_tensor& kk_in,
261 int depth_in = 3) ;
262
272 Time_slice(const Scalar& lapse_in, const Vector& shift_in,
273 const Evolution_std<Sym_tensor>& gamma_in) ;
274
285 Time_slice(const Map& mp, const Base_vect& triad, int depth_in = 3) ;
286
301 Time_slice(const Map& mp, const Base_vect& triad, FILE* fich,
302 bool partial_read, int depth_in = 3) ;
303
304 Time_slice(const Time_slice& ) ;
305
306 protected:
310 explicit Time_slice(int depth_in) ;
311
312 public:
313 virtual ~Time_slice() ;
314
315
316 // Memory management
317 // -----------------
318 protected:
319
321 virtual void del_deriv() const ;
322
324 void set_der_0x0() const ;
325
326
327 // Mutators / assignment
328 // ---------------------
329 public:
331 void operator=(const Time_slice&) ;
332
334 void set_scheme_order(int ord) {
335 assert ((0<= ord)&&(ord < 4)) ;
336 scheme_order = ord ; } ;
337
338 // Accessors
339 // ---------
340 public:
341
343 int get_scheme_order() const { return scheme_order ; } ;
344
346 int get_latest_j() const {return jtime; } ;
347
349 const Evolution_std<double>& get_time() const {return the_time; } ;
350
352 virtual const Scalar& nn() const ;
353
355 virtual const Vector& beta() const ;
356
358 const Metric& gam() const ;
359
363 virtual const Sym_tensor& gam_dd() const ;
364
368 virtual const Sym_tensor& gam_uu() const ;
369
373 virtual const Sym_tensor& k_dd() const ;
374
378 virtual const Sym_tensor& k_uu() const ;
379
383 virtual const Scalar& trk() const ;
384
385
386 // Computational functions
387 // -----------------------
388 public:
404 Tbl check_hamiltonian_constraint(const Scalar* energy_density = 0x0,
405 ostream& ost = cout, bool verb=true) const ;
406
422 Tbl check_momentum_constraint(const Vector* momentum_density = 0x0,
423 ostream& ost = cout, bool verb=true) const ;
424
447 Tbl check_dynamical_equations(const Sym_tensor* strain_tensor = 0x0,
448 const Scalar* energy_density = 0x0,
449 ostream& ost = cout, bool verb=true) const ;
450
455 virtual double adm_mass() const ;
456
457 // Outputs
458 // -------
459 protected:
461 virtual ostream& operator>>(ostream& ) const ;
462
464 friend ostream& operator<<(ostream& , const Time_slice& ) ;
465
466 public:
473 void save(const char* rootname) const ;
474
475 protected:
484 virtual void sauve(FILE* fich, bool partial_save) const ;
485
486};
487
488ostream& operator<<(ostream& , const Time_slice& ) ;
489
490
491
492 //---------------------------//
493 // class Time_slice_conf //
494 //---------------------------//
495
502
503 // Data :
504 // -----
505 protected:
506
510 const Metric_flat& ff ;
511
521
526
527
534
546
551
556
557 // Derived data :
558 // ------------
559 protected:
563 mutable Metric* p_tgamma ;
564
566 mutable Scalar* p_psi4 ;
567
569 mutable Scalar* p_ln_psi ;
570
574 mutable Vector* p_hdirac ;
575
580 mutable Vector* p_vec_X ;
581
582 // Constructors - Destructor
583 // -------------------------
584 public:
585
610 Time_slice_conf(const Scalar& lapse_in, const Vector& shift_in,
611 const Metric_flat& ff_in, const Scalar& psi_in,
612 const Sym_tensor& hh_in, const Sym_tensor& hata_in,
613 const Scalar& trk_in, int depth_in = 3) ;
614
615
630 Time_slice_conf(const Scalar& lapse_in, const Vector& shift_in,
631 const Sym_tensor& gamma_in, const Sym_tensor& kk_in,
632 const Metric_flat& ff_in, int depth_in = 3) ;
633
646 Time_slice_conf(const Map& mp, const Base_vect& triad,
647 const Metric_flat& ff_in, int depth_in = 3) ;
648
665 Time_slice_conf(const Map& mp, const Base_vect& triad,
666 const Metric_flat& ff_in, FILE* fich,
667 bool partial_read, int depth_in = 3) ;
668
670
671 virtual ~Time_slice_conf() ;
672
673
674 // Memory management
675 // -----------------
676 protected:
677
679 virtual void del_deriv() const ;
680
682 void set_der_0x0() const ;
683
684
685 // Mutators / assignment
686 // ---------------------
687 public:
689 void operator=(const Time_slice_conf&) ;
690
692 void operator=(const Time_slice&) ;
693
704 virtual void set_psi_del_npsi(const Scalar& psi_in) ;
705
716 virtual void set_psi_del_n(const Scalar& psi_in) ;
717
722 virtual void set_npsi_del_psi(const Scalar& npsi_in) ;
723
728 virtual void set_npsi_del_n(const Scalar& npsi_in) ;
729
738 virtual void set_hh(const Sym_tensor& hh_in) ;
739
746 virtual void set_hata(const Sym_tensor& hata_in) ;
747
751 virtual void set_hata_TT(const Sym_tensor_tt& hata_tt) ;
752
759 virtual void set_hata_from_XAB(Param* par_bc=0x0, Param* par_mat=0x0) ;
760
761 // Accessors
762 // ---------
763 public:
764
765 // Virtual functions from base class Time_slice:
766 // ---------------------------------------------
767
769 virtual const Scalar& nn() const ;
770
774 virtual const Sym_tensor& gam_dd() const ;
775
779 virtual const Sym_tensor& gam_uu() const ;
780
784 virtual const Sym_tensor& k_dd() const ;
785
789 virtual const Sym_tensor& k_uu() const ;
790
791 // Virtual functions from this class:
792 // ----------------------------------
793
798 virtual const Scalar& A_hata() const ;
799
804 virtual const Scalar& B_hata() const ;
805
813 virtual const Scalar& psi() const ;
814
816 const Scalar& psi4() const ;
817
819 const Scalar& ln_psi() const ;
820
823 virtual const Scalar& npsi() const ;
824
829 virtual const Metric& tgam() const ;
830
837 virtual const Sym_tensor& hh(Param* = 0x0, Param* = 0x0) const ;
838
844 virtual const Sym_tensor& hata() const ;
845
851 virtual Sym_tensor aa() const ;
852
856 virtual const Scalar& trk() const ;
857
861 virtual const Vector& hdirac() const ;
862
866 virtual const Vector& vec_X(int method_poisson=6) const ;
867
868 // Computational methods
869 // ---------------------
870 public:
875 (const Vector& hat_S, const Sym_tensor_tt& hata_tt,
876 int iter_max = 200, double precis = 1.e-12,
877 double relax = 0.8, int methode_poisson = 6) ;
878
884 virtual void set_AB_hata(const Scalar& A_in, const Scalar& B_in) ;
885
919 virtual void initial_data_cts(const Sym_tensor& uu, const Scalar& trk_in,
920 const Scalar& trk_point, double pdt, double precis = 1.e-12,
921 int method_poisson_vect = 6, const char* graph_device = 0x0,
922 const Scalar* ener_dens = 0x0, const Vector* mom_dens = 0x0,
923 const Scalar* trace_stress = 0x0 ) ;
924
929 virtual double adm_mass() const ;
930
942 void check_psi_dot(Tbl& tlnpsi_dot, Tbl& tdiff, Tbl& tdiff_rel) const ;
943
944 // Outputs
945 // -------
946 protected:
948 virtual ostream& operator>>(ostream& ) const ;
949
958 virtual void sauve(FILE* fich, bool partial_save) const ;
959
960} ;
961 //----------------------------//
962 // class Tslice_dirac_max //
963 //----------------------------//
964
972
973 // Data :
974 // -----
975 protected:
981
987
993
999
1005
1011
1014
1015
1016 // Constructors - Destructor
1017 // -------------------------
1018 public:
1043 Tslice_dirac_max(const Scalar& lapse_in, const Vector& shift_in,
1044 const Metric_flat& ff_in, const Scalar& psi_in,
1045 const Sym_tensor_trans& hh_in, const Sym_tensor& hata_in,
1046 int depth_in = 3) ;
1047
1060 Tslice_dirac_max(const Map& mp, const Base_vect& triad,
1061 const Metric_flat& ff_in, int depth_in = 3) ;
1062
1079 Tslice_dirac_max(const Map& mp, const Base_vect& triad,
1080 const Metric_flat& ff_in, FILE* fich,
1081 bool partial_read, int depth_in = 3) ;
1082
1084 Tslice_dirac_max(const Star_rot_Dirac& star, double pdt, int depth_in = 3) ;
1085
1087
1088 virtual ~Tslice_dirac_max() ;
1089
1090
1091 // Mutators / assignment
1092 // ---------------------
1093 public:
1095 void operator=(const Tslice_dirac_max&) ;
1096
1097 // Virtual functions from base class Time_slice_conf:
1098 // -------------------------------------------------
1099
1108 virtual void set_hh(const Sym_tensor& hh_in) ;
1109
1143 virtual void initial_data_cts(const Sym_tensor& uu, const Scalar& trk_in,
1144 const Scalar& trk_point, double pdt, double precis = 1.e-12,
1145 int method_poisson_vect = 6, const char* graph_device = 0x0,
1146 const Scalar* ener_dens = 0x0, const Vector* mom_dens = 0x0,
1147 const Scalar* trace_stress = 0x0 ) ;
1148
1149 // Virtual functions from this class:
1150 // ----------------------------------
1151
1159 virtual void set_khi_mu(const Scalar& khi_in, const Scalar& mu_in) ;
1160
1167 virtual void set_AB_hh(const Scalar& A_in, const Scalar& B_in) ;
1168
1175 virtual void set_trh(const Scalar& trh_in) ;
1176
1187 virtual Scalar solve_psi(const Scalar* ener_dens=0x0) const ;
1188
1204 virtual Scalar solve_npsi(const Scalar* ener_dens=0x0,
1205 const Scalar* trace_stress=0x0) const ;
1206
1217 virtual Vector solve_beta(int method = 6) const ;
1218
1240 void evolve(double pdt, int nb_time_steps, int niter_elliptic,
1241 double relax_elliptic, int check_mod, int save_mod,
1242 int method_poisson_vect = 6, int nopause = 1,
1243 const char* graph_device = 0x0, bool verbose=true,
1244 const Scalar* ener_euler = 0x0,
1245 const Vector* mom_euler = 0x0, const Scalar* s_euler = 0x0,
1246 const Sym_tensor* strain_euler = 0x0) ;
1247
1252 virtual double adm_mass() const ;
1253
1254 protected:
1263 void compute_sources(const Sym_tensor* strain_tensor = 0x0) const ;
1264
1266 void initialize_sources_copy() const ;
1267
1275 void hh_det_one(int j, Param* par_bc = 0x0, Param* par_mat = 0x0) const ;
1276
1282 void hh_det_one(const Sym_tensor_tt& hijtt, Param* par_mat = 0x0) const ;
1283
1284 // Accessors
1285 // ---------
1286 public:
1287 // Virtual functions from base class Time_slice_conf:
1288 // -------------------------------------------------
1289
1296 virtual const Sym_tensor& hh(Param* par_bc = 0x0, Param* par_mat = 0x0) const ;
1297
1302 virtual const Scalar& trk() const ;
1303
1308 virtual const Vector& hdirac() const ;
1309
1310 // Virtual functions from this class:
1311 // ----------------------------------
1312
1317 virtual const Scalar& A_hh() const ;
1318
1323 virtual const Scalar& B_hh() const ;
1324
1329 virtual const Scalar& trh() const ;
1330
1331
1332 // Outputs
1333 // -------
1334 protected:
1336 virtual ostream& operator>>(ostream& ) const ;
1337
1346 virtual void sauve(FILE* fich, bool partial_save) const ;
1347
1348};
1349}
1350#endif
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition base_vect.h:105
Time evolution with full storage (*** under development ***).
Definition evolution.h:270
Time evolution with partial storage (*** under development ***).
Definition evolution.h:371
Flat metric for tensor calculation.
Definition metric.h:261
Metric for tensor calculation.
Definition metric.h:90
Parameter storage.
Definition param.h:125
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
Class for relativistic rotating stars in Dirac gauge and maximal slicing.
Transverse symmetric tensors of rank 2.
Definition sym_tensor.h:611
Transverse and traceless symmetric tensors of rank 2.
Definition sym_tensor.h:933
Class intended to describe valence-2 symmetric tensors.
Definition sym_tensor.h:226
Basic array class.
Definition tbl.h:161
virtual const Scalar & npsi() const
Factor at the current time step (jtime ).
Evolution_std< Sym_tensor > hh_evol
Values at successive time steps of the components .
Definition time_slice.h:533
virtual ~Time_slice_conf()
Destructor.
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
virtual const Sym_tensor & gam_dd() const
Induced metric (covariant components ) at the current time step (jtime ).
virtual const Sym_tensor & hh(Param *=0x0, Param *=0x0) const
Deviation of the conformal metric from the flat metric : .
virtual void set_psi_del_n(const Scalar &psi_in)
Sets the conformal factor relating the physical metric to the conformal one: .
virtual void set_hata_TT(const Sym_tensor_tt &hata_tt)
Sets the TT part of (see member hata_evol ).
virtual void sauve(FILE *fich, bool partial_save) const
Total or partial saves in a binary file.
virtual const Sym_tensor & k_dd() const
Extrinsic curvature tensor (covariant components ) at the current time step (jtime ).
virtual void set_npsi_del_psi(const Scalar &npsi_in)
Sets the factor at the current time step (jtime ) and deletes the value of .
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: .
virtual const Sym_tensor & hata() const
Conformal representation of the traceless part of the extrinsic curvature: .
virtual const Sym_tensor & k_uu() const
Extrinsic curvature tensor (contravariant components ) at the current time step (jtime ).
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 .
Definition time_slice.h:545
virtual const Vector & hdirac() const
Vector which vanishes in Dirac gauge.
virtual const Scalar & A_hata() const
Returns the potential A of .
virtual const Vector & vec_X(int method_poisson=6) const
Vector representing the longitudinal part of .
void operator=(const Time_slice_conf &)
Assignment to another Time_slice_conf.
virtual const Scalar & psi() const
Conformal factor relating the physical metric to the conformal one: .
Metric * p_tgamma
Pointer on the conformal metric at the current time step (jtime).
Definition time_slice.h:563
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 ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator<<).
virtual void del_deriv() const
Deletes all the derived quantities.
void check_psi_dot(Tbl &tlnpsi_dot, Tbl &tdiff, Tbl &tdiff_rel) const
Checks the relation.
Evolution_std< Scalar > npsi_evol
Values at successive time steps of the factor .
Definition time_slice.h:525
Vector * p_vec_X
Pointer on the vector representing the longitudinal part of .
Definition time_slice.h:580
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: .
Vector * p_hdirac
Pointer on the vector (which vanishes in Dirac gauge), at the current time step (jtime).
Definition time_slice.h:574
const Scalar & psi4() const
Factor at the current time step (jtime ).
virtual void set_npsi_del_n(const Scalar &npsi_in)
Sets the factor at the current time step (jtime ) and deletes the value of N.
Evolution_std< Scalar > psi_evol
Values at successive time steps of the conformal factor relating the physical metric to the conform...
Definition time_slice.h:520
virtual void set_AB_hata(const Scalar &A_in, const Scalar &B_in)
Sets the potentials A and of the TT part (see the documentation of Sym_tensor for details).
virtual void set_hata_from_XAB(Param *par_bc=0x0, Param *par_mat=0x0)
Sets the conformal representation of the traceless part of the extrinsic curvature from its potentia...
virtual const Metric & tgam() const
Conformal metric Returns the value at the current time step (jtime ).
Time_slice_conf(const Scalar &lapse_in, const Vector &shift_in, const Metric_flat &ff_in, const Scalar &psi_in, const Sym_tensor &hh_in, const Sym_tensor &hata_in, const Scalar &trk_in, int depth_in=3)
Constructor from conformal decomposition.
virtual const Sym_tensor & gam_uu() const
Induced metric (contravariant components ) at the current time step (jtime ).
Evolution_std< Scalar > A_hata_evol
Potential A associated with the symmetric tensor .
Definition time_slice.h:550
virtual double adm_mass() const
Returns the ADM mass (geometrical units) at the current step.
const Metric_flat & ff
Pointer on the flat metric with respect to which the conformal decomposition is performed.
Definition time_slice.h:510
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 from the conformally-rescaled momentum , using the momentum constraint.
Scalar * p_ln_psi
Pointer on the logarithm of at the current time step (jtime).
Definition time_slice.h:569
Evolution_std< Scalar > B_hata_evol
Potential associated with the symmetric tensor .
Definition time_slice.h:555
virtual void set_hh(const Sym_tensor &hh_in)
Sets the deviation of the conformal metric from the flat metric : .
virtual const Scalar & trk() const
Trace K of the extrinsic curvature at the current time step (jtime ).
Scalar * p_psi4
Pointer on the factor at the current time step (jtime).
Definition time_slice.h:566
Spacelike time slice of a 3+1 spacetime.
Definition time_slice.h:176
virtual const Sym_tensor & k_dd() const
Extrinsic curvature tensor (covariant components ) at the current time step (jtime ).
int jtime
Time step index of the latest slice.
Definition time_slice.h:193
Evolution_std< Scalar > trk_evol
Values at successive time steps of the trace K of the extrinsic curvature.
Definition time_slice.h:227
Metric * p_gamma
Pointer on the induced metric at the current time step (jtime).
Definition time_slice.h:242
Evolution_std< Sym_tensor > k_dd_evol
Values at successive time steps of the covariant components of the extrinsic curvature tensor .
Definition time_slice.h:211
virtual void sauve(FILE *fich, bool partial_save) const
Total or partial saves in a binary file.
Definition time_slice.C:510
void operator=(const Time_slice &)
Assignment to another Time_slice.
Definition time_slice.C:391
Evolution_full< Tbl > adm_mass_evol
ADM mass at each time step, since the creation of the slice.
Definition time_slice.h:236
virtual const Scalar & trk() const
Trace K of the extrinsic curvature at the current time step (jtime ).
virtual const Scalar & nn() const
Lapse function N at the current time step (jtime ).
void save(const char *rootname) const
Saves in a binary file.
Definition time_slice.C:464
virtual const Sym_tensor & gam_dd() const
Induced metric (covariant components ) at the current time step (jtime ).
void set_der_0x0() const
Sets to 0x0 all the pointers on derived quantities.
Definition time_slice.C:380
Evolution_std< Sym_tensor > gam_dd_evol
Values at successive time steps of the covariant components of the induced metric .
Definition time_slice.h:201
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.
Evolution_std< Sym_tensor > gam_uu_evol
Values at successive time steps of the contravariant components of the induced metric .
Definition time_slice.h:206
int scheme_order
Order of the finite-differences scheme for the computation of time derivatives.
Definition time_slice.h:190
virtual const Vector & beta() const
shift vector at the current time step (jtime )
virtual ~Time_slice()
Destructor.
Definition time_slice.C:360
virtual void del_deriv() const
Deletes all the derived quantities.
Definition time_slice.C:370
int get_latest_j() const
Gets the latest value of time step index.
Definition time_slice.h:346
const Metric & gam() const
Induced metric at the current time step (jtime ).
virtual const Sym_tensor & k_uu() const
Extrinsic curvature tensor (contravariant components ) at the current time step (jtime ).
int depth
Number of stored time slices.
Definition time_slice.h:182
const Evolution_std< double > & get_time() const
Gets the time coordinate t at successive time steps.
Definition time_slice.h:349
Evolution_std< Scalar > n_evol
Values at successive time steps of the lapse function N.
Definition time_slice.h:219
virtual const Sym_tensor & gam_uu() const
Induced metric (contravariant components ) at the current time step (jtime ).
virtual double adm_mass() const
Returns the ADM mass (geometrical units) at the current step.
Time_slice(const Scalar &lapse_in, const Vector &shift_in, const Sym_tensor &gamma_in, const Sym_tensor &kk_in, int depth_in=3)
General constructor (Hamiltonian-like).
Definition time_slice.C:115
int get_scheme_order() const
Gets the order of the finite-differences scheme.
Definition time_slice.h:343
Evolution_std< double > the_time
Time label of each slice.
Definition time_slice.h:196
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator<<).
Definition time_slice.C:414
Evolution_std< Sym_tensor > k_uu_evol
Values at successive time steps of the contravariant components of the extrinsic curvature tensor .
Definition time_slice.h:216
friend ostream & operator<<(ostream &, const Time_slice &)
Display.
Definition time_slice.C:456
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 set_scheme_order(int ord)
Sets the order of the finite-differences scheme.
Definition time_slice.h:334
Evolution_std< Vector > beta_evol
Values at successive time steps of the shift vector .
Definition time_slice.h:222
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.
virtual const Scalar & trh() const
Computes the trace h, with respect to the flat metric ff , of .
virtual double adm_mass() const
Returns the ADM mass at (geometrical units) the current step.
virtual const Sym_tensor & hh(Param *par_bc=0x0, Param *par_mat=0x0) const
Deviation of the conformal metric from the flat metric : .
virtual Scalar solve_psi(const Scalar *ener_dens=0x0) const
Solves the elliptic equation for the conformal factor $\Psi$ (Hamiltonian constraint).
Evolution_std< Scalar > source_B_hh_evol
The potential of the source of equation for .
Definition time_slice.h:998
Evolution_std< Scalar > source_A_hh_evol
The A potential of the source of equation for .
Definition time_slice.h:992
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 const Scalar & B_hh() 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...
Evolution_std< Scalar > A_hh_evol
The A potential of .
Definition time_slice.h:980
virtual ~Tslice_dirac_max()
Destructor.
void operator=(const Tslice_dirac_max &)
Assignment to another Tslice_dirac_max.
virtual const Scalar & trk() const
Trace K of the extrinsic curvature at the current time step (jtime ).
void hh_det_one(int j, Param *par_bc=0x0, Param *par_mat=0x0) const
Computes from the values of A and and using the condition , which fixes the trace of .
Evolution_std< Scalar > source_B_hata_evol
The potential of the source of equation for .
virtual void set_hh(const Sym_tensor &hh_in)
Sets the deviation of the conformal metric from the flat metric : .
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 equa...
virtual Scalar solve_npsi(const Scalar *ener_dens=0x0, const Scalar *trace_stress=0x0) const
Solves the elliptic equation for (maximal slicing condition + Hamiltonian constraint).
virtual void set_khi_mu(const Scalar &khi_in, const Scalar &mu_in)
Sets the potentials and of the TT part of (see the documentation of Sym_tensor_tt for details).
virtual void sauve(FILE *fich, bool partial_save) const
Total or partial saves in a binary file.
virtual void set_trh(const Scalar &trh_in)
Sets the trace, with respect to the flat metric ff , of .
virtual const Scalar & A_hh() const
Returns the potential A of .
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator<<).
virtual const Vector & hdirac() const
Vector which vanishes in Dirac gauge.
virtual Vector solve_beta(int method=6) const
Solves the elliptic equation for the shift vector from (Eq.
Evolution_std< Scalar > trh_evol
The trace, with respect to the flat metric ff , of .
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.
void initialize_sources_copy() const
Copy the sources source_A_XXX_evol and source_B_XXX_evol to all time-steps.
virtual void set_AB_hh(const Scalar &A_in, const Scalar &B_in)
Sets the potentials A and of the TT part of (see the documentation of Sym_tensor for details).
Evolution_std< Scalar > B_hh_evol
The potential of .
Definition time_slice.h:986
Evolution_std< Scalar > source_A_hata_evol
The potential A of the source of equation for .
Tensor field of valence 1.
Definition vector.h:188
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142