164#include "utilitaires.h"
165#include "time_slice.h"
169#include "evolution.h"
203 ff_in.con(), psi_in*psi_in*psi_in*psi_in*psi_in*psi_in*aa_in,
257 bool partial_read,
int depth_in)
259 fich, partial_read, depth_in),
280 for (
int j=jmin; j<=
jtime; j++) {
281 fread_be(&indicator,
sizeof(
int), 1, fich) ;
282 if (indicator == 1) {
289 for (
int j=jmin; j<=
jtime; j++) {
290 fread_be(&indicator,
sizeof(
int), 1, fich) ;
291 if (indicator == 1) {
292 Scalar nn_auto_file(
mp, *(
mp.get_mg()), fich) ;
298 for (
int j=jmin; j<=
jtime; j++) {
299 fread_be(&indicator,
sizeof(
int), 1, fich) ;
300 if (indicator == 1) {
301 Scalar psi_auto_file(
mp, *(
mp.get_mg()), fich) ;
307 for (
int j=jmin; j<=
jtime; j++) {
308 fread_be(&indicator,
sizeof(
int), 1, fich) ;
309 if (indicator == 1) {
310 Vector beta_auto_file(
mp, mpi.get_bvect_spher(), fich) ;
326 Scalar trK_point_file (
mp, *(
mp.get_mg()), fich) ;
385 flux <<
'\n' <<
"radius of the horizon : " <<
radius <<
'\n' ;
386 flux <<
"boost in x-direction : " <<
boost_x <<
'\n' ;
387 flux <<
"boost in z-direction : " <<
boost_z <<
'\n' ;
388 flux <<
"angular velocity omega : " <<
omega_hor() <<
'\n' ;
389 flux <<
"area of the horizon : " <<
area_hor() <<
'\n' ;
390 flux <<
"ang. mom. of horizon : " <<
ang_mom_hor() <<
'\n' ;
391 flux <<
"ADM ang. mom. : " <<
ang_mom_adm() <<
'\n' ;
392 flux <<
"Mass of the horizon : " <<
mass_hor() <<
'\n' ;
393 flux <<
"ADM Mass : " <<
adm_mass() <<
'\n' ;
422 for (
int j=jmin; j<=
jtime; j++) {
423 int indicator = (
psi_evol.is_known(j)) ? 1 : 0 ;
424 fwrite_be(&indicator,
sizeof(
int), 1, fich) ;
425 if (indicator == 1)
psi_evol[j].sauve(fich) ;
429 for (
int j=jmin; j<=
jtime; j++) {
430 int indicator = (
n_auto_evol.is_known(j)) ? 1 : 0 ;
431 fwrite_be(&indicator,
sizeof(
int), 1, fich) ;
436 for (
int j=jmin; j<=
jtime; j++) {
438 fwrite_be(&indicator,
sizeof(
int), 1, fich) ;
443 for (
int j=jmin; j<=
jtime; j++) {
445 fwrite_be(&indicator,
sizeof(
int), 1, fich) ;
590 Vector dn_comp (
mp, COV,
mp.get_bvect_cart()) ;
595 for (
int i=1 ; i<=3 ; i++){
596 if (auxi(i).get_etat() != ETATZERO)
603 for (
int i=1 ; i<=3 ; i++){
682 Vector dpsi_comp (
mp, COV,
mp.get_bvect_cart()) ;
687 for (
int i=1 ; i<=3 ; i++){
688 if (auxi(i).get_etat() != ETATZERO)
695 for (
int i=1 ; i<=3 ; i++){
719 Vector tmp_vect (
mp, CON,
mp.get_bvect_cart()) ;
743 auxi = 0.5 - 0.5/
mp.r ;
774 Vector temp_vect1(
mp, CON,
mp.get_bvect_spher()) ;
775 temp_vect1.
set(1) = 0.0/
mp.r/
mp.r ;
776 temp_vect1.
set(2) = 0. ;
777 temp_vect1.
set(3) = 0. ;
780 Vector temp_vect2(
mp, CON,
mp.get_bvect_spher()) ;
794 trK.set_etat_zero() ;
837 Vector temp_vect(
mp, CON,
mp.get_bvect_spher()) ;
848 int nnt =
mp.get_mg()->get_nt(1) ;
849 int nnp =
mp.get_mg()->get_np(1) ;
853 for (
int k=0; k<nnp; k++)
854 for (
int j=0; j<nnt; j++){
855 if (
nn().val_grid_point(1, k, j , 0) < 1e-12){
866 cout <<
"regul = " <<
regul << endl ;
871 for (
int i=1 ; i<=3 ; i++)
872 for (
int j=i ; j<=3 ; j++) {
873 auxi = aa_new(i, j) ;
874 auxi = division_xpun (auxi, 0) ;
875 aa_new.
set(i,j) = auxi / nn_sxpun / 2. ;
906 Scalar psi_perturb (
pow(gamm(3,3), 0.25)) ;
910 cout <<
"met_gamt" << endl <<
norme(
met_gamt.cov()(1,1)) << endl
915 cout <<
"determinant" <<
norme(
met_gamt.determinant()) << endl ;
921void Isol_hor::aa_kerr_ww(
double mm,
double aaa) {
931 Scalar rbl = rr + mm + (mm*mm - aaa*aaa) / (4*rr) ;
939 ww_pert = - 1*(mm*aaa*aaa*aaa*
pow(
sint, 4.)*
cost) * sigma_inv ;
940 ww_pert.set_spectral_va().set_base_r(0,
R_CHEBPIM_P) ;
941 for (
int l=1; l<
nz-1; l++)
942 ww_pert.set_spectral_va().set_base_r(l,
R_CHEB) ;
943 ww_pert.set_spectral_va().set_base_r(
nz-1,
R_CHEBU) ;
944 ww_pert.set_spectral_va().set_base_t(
T_COSSIN_CI) ;
945 ww_pert.set_spectral_va().set_base_p(
P_COSSIN) ;
952 Vector dw_by (
mp, COV,
mp.get_bvect_spher()) ;
956 dw_by.set(2).set_spectral_va().set_base_r(0,
R_CHEBPIM_P) ;
957 for (
int l=1; l<
nz-1; l++)
958 dw_by.set(2).set_spectral_va().set_base_r(l,
R_CHEB) ;
959 dw_by.set(2).set_spectral_va().set_base_r(
nz-1,
R_CHEBU) ;
960 dw_by.set(2).set_spectral_va().set_base_t(
T_COSSIN_SI) ;
961 dw_by.set(2).set_spectral_va().set_base_p(
P_COSSIN) ;
971 Vector dw_pert(ww_pert.derive_con(ff)) ;
990 Scalar aquad = aquad_1 + aquad_2 + aquad_3 ;
1015 Scalar s_r (ww_pert.derive_cov(ff)(2)) ;
1026 Scalar s_t (ww_pert.derive_cov(ff)(1)) ;
1032 temp.div_rsint_dzpuis(2) ;
1048 aij.set(3,1) = s_r ;
1049 aij.set(3,1).div_rsint() ;
1050 aij.set(3,2) = s_t ;
1051 aij.set(3,2).div_rsint() ;
1053 Base_val base_31 (aij(3,1).get_spectral_va().get_base()) ;
1054 Base_val base_32 (aij(3,2).get_spectral_va().get_base()) ;
1056 aij.set(3,1) = aij(3,1) *
pow(
gam_dd()(1,1), 5./3.)
1058 aij.set(3,1) = aij(3,1) *
pow(
psi(), -6.) ;
1059 aij.set(3,1).set_spectral_va().set_base(base_31) ;
1060 aij.set(3,2) = aij(3,2) *
pow(
gam_dd()(1,1), 5./3.)
1062 aij.set(3,2) = aij(3,2) *
pow(
psi(), -6.) ;
1063 aij.set(3,2).set_spectral_va().set_base(base_32) ;
1083 kij = kij *
pow(
gam().determinant(), -1./3.) ;
1084 kij.std_spectral_base() ;
1112 double axibreak =
mp.integrale_surface(integrand, 1.0000000001)/
1119double fonc_expansion(
double rr,
const Param& par_expansion) {
1128void Isol_hor::adapt_hor(
double c_min,
double c_max) {
1131 Scalar app_hor(
mp) ;
1132 app_hor.annule_hard() ;
1136 double precis = 1.e-13 ;
1141 for (
int nt=0; nt<
mp.get_mg()->get_nt(1); nt++)
1142 for (
int np=0; np<
mp.get_mg()->get_np(1); np++) {
1144 double theta =
mp.get_mg()->get_grille3d(1)->tet[nt] ;
1145 double phi =
mp.get_mg()->get_grille3d(1)->phi[np] ;
1147 Param par_expansion ;
1149 par_expansion.add_double(theta, 0) ;
1150 par_expansion.add_double(
phi, 1) ;
1151 double r_app_hor =
zerosec_b(fonc_expansion, par_expansion, c_min, c_max,
1152 precis, nitmax, nit) ;
1154 app_hor.set_grid_point(1, np, nt, 0) = r_app_hor ;
1163 Scalar trans_11 (
mp) ;
1165 r_new.annule_hard() ;
1167 for (
int l=1; l<
nz; l++)
1168 for (
int nr=0; nr<
mp.get_mg()->get_nr(1); nr++)
1169 for (
int nt=0; nt<
mp.get_mg()->get_nt(1); nt++)
1170 for (
int np=0; np<
mp.get_mg()->get_np(1); np++) {
1171 r_new.set_grid_point(l, np, nt, nr) = rr.val_grid_point(l, np, nt, nr) -
1172 app_hor.val_grid_point(1, np, nt, 0) + 1 ;
1176 r_new.std_spectral_base() ;
1182 Scalar trans_12 (r_new.dsdt()) ;
1184 Scalar trans_13 (r_new.stdsdp()) ;
1186 for (
int nr=0; nr<
mp.get_mg()->get_nr(1); nr++)
1187 for (
int nt=0; nt<
mp.get_mg()->get_nt(1); nt++)
1188 for (
int np=0; np<
mp.get_mg()->get_np(1); np++) {
1189 trans_12.set_grid_point(
nz-1, np, nt, nr) = trans_12.val_grid_point(1, np, nt, nr) ;
1190 trans_13.set_grid_point(
nz-1, np, nt, nr) = trans_13.val_grid_point(1, np, nt, nr) ;
1194 Tensor trans (
mp, 2, comp,
mp.get_bvect_spher()) ;
1195 trans.set(1,1) = 1 ;
1198 trans.set(2,2) = 1. ;
1199 trans.set(3,3) = 1. ;
1200 trans.set(2,1) = 0. ;
1201 trans.set(3,1) = 0. ;
1202 trans.set(3,2) = 0. ;
1203 trans.set(2,3) = 0. ;
1204 trans.std_spectral_base() ;
1206 cout <<
"trans(1,3)" << endl <<
norme(trans(1,3)) << endl ;
1208 Sym_tensor gamma_uu (
gam().con()) ;
1209 Sym_tensor kk_uu (
k_uu()) ;
1211 gamma_uu =
contract(gamma_uu, 0, 1, trans * trans, 1, 3) ;
1212 kk_uu =
contract(kk_uu, 0, 1, trans * trans, 1, 3) ;
1214 Sym_tensor copie_gamma (gamma_uu) ;
1215 Sym_tensor copie_kk (kk_uu) ;
1218 kk_uu.dec_dzpuis(2) ;
1219 for(
int i=1; i<=3; i++)
1220 for(
int j=i; j<=3; j++){
1221 kk_uu.set(i,j).annule_hard() ;
1222 gamma_uu.set(i,j).annule_hard() ;
1225 copie_kk.dec_dzpuis(2) ;
1227 Scalar expa_trans(
mp) ;
1228 expa_trans.annule_hard() ;
1243 for(
int i=1; i<=3; i++)
1244 for(
int j=i; j<=3; j++)
1245 for (
int l=1; l<
nz; l++)
1246 for (
int nr=0; nr<
mp.get_mg()->get_nr(1); nr++)
1247 for (
int nt=0; nt<
mp.get_mg()->get_nt(1); nt++)
1248 for (
int np=0; np<
mp.get_mg()->get_np(1); np++) {
1250 double theta =
mp.get_mg()->get_grille3d(1)->tet[nt] ;
1251 double phi =
mp.get_mg()->get_grille3d(1)->phi[np] ;
1252 double r_inv = rr.val_grid_point(l, np, nt, nr) +
1253 app_hor.val_grid_point(1, np, nt, 0) - 1. ;
1255 gamma_uu.set(i,j).set_grid_point(l, np, nt, nr) =
1256 copie_gamma(i,j).val_point(r_inv, theta,
phi) ;
1257 kk_uu.set(i,j).set_grid_point(l, np, nt, nr) =
1258 copie_kk(i,j).val_point(r_inv, theta,
phi) ;
1260 expa_trans.set_grid_point(l, np, nt, nr) = expa.
val_point(r_inv, theta,
phi) ;
1262 kk_uu.std_spectral_base() ;
1263 gamma_uu.std_spectral_base() ;
1264 expa_trans.std_spectral_base() ;
1266 for (
int l=1; l<
nz; l++)
1267 for (
int nr=0; nr<
mp.get_mg()->get_nr(1); nr++)
1268 for (
int nt=0; nt<
mp.get_mg()->get_nt(1); nt++)
1269 for (
int np=0; np<
mp.get_mg()->get_np(1); np++) {
1270 gamma_uu.set(1,2).set_grid_point(l,np,nt,nr) = gamma_uu.set(1,2).val_grid_point(l,np,nt,nr)
1271 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1272 - 1 / rr.val_grid_point(l, np, nt, nr) ) ;
1273 gamma_uu.set(1,3).set_grid_point(l,np,nt,nr) = gamma_uu.set(1,3).val_grid_point(l,np,nt,nr)
1274 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1275 - 1 / rr.val_grid_point(l, np, nt, nr) ) ;
1276 gamma_uu.set(2,2).set_grid_point(l,np,nt,nr) = gamma_uu.set(2,2).val_grid_point(l,np,nt,nr)
1277 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1278 - 1 / rr.val_grid_point(l, np, nt, nr) )
1279 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1280 - 1 / rr.val_grid_point(l, np, nt, nr) ) ;
1281 gamma_uu.set(2,3).set_grid_point(l,np,nt,nr) = gamma_uu.set(2,3).val_grid_point(l,np,nt,nr)
1282 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1283 - 1 / rr.val_grid_point(l, np, nt, nr) )
1284 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1285 - 1 / rr.val_grid_point(l, np, nt, nr) ) ;
1286 gamma_uu.set(3,3).set_grid_point(l,np,nt,nr) = gamma_uu.set(3,3).val_grid_point(l,np,nt,nr)
1287 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1288 - 1 / rr.val_grid_point(l, np, nt, nr) )
1289 / (1 +app_hor.val_grid_point(1, np, nt, 0)/ rr.val_grid_point(l, np, nt, nr)
1290 - 1 / rr.val_grid_point(l, np, nt, nr) ) ;
1306 cout <<
"gamma_uu(1,1)" << endl <<
norme(gamma_uu(1,1)) << endl ;
1307 cout <<
"gamma_uu(2,1)" << endl <<
norme(gamma_uu(2,1)) << endl ;
1308 cout <<
"gamma_uu(3,1)" << endl <<
norme(gamma_uu(3,1)) << endl ;
1309 cout <<
"gamma_uu(2,2)" << endl <<
norme(gamma_uu(2,2)) << endl ;
1310 cout <<
"gamma_uu(3,2)" << endl <<
norme(gamma_uu(3,2)) << endl ;
1311 cout <<
"gamma_uu(3,3)" << endl <<
norme(gamma_uu(3,3)) << endl ;
1313 kk_uu.inc_dzpuis(2) ;
1314 expa_trans.inc_dzpuis(2) ;
1316 Metric gamm (gamma_uu) ;
Bases of the spectral expansions.
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
virtual const Vector & dnn() const
Covariant derivative of the lapse function at the current time step jtime.
Scalar expansion() const
Expansion of the outgoing null normal ( ).
void met_kerr_perturb()
Initialisation of the metric tilde from equation (15) of Dain (2002).
const Map_af & get_mp() const
Returns the mapping (readonly).
virtual const Scalar & aa_quad() const
Conformal representation .
Evolution_std< Sym_tensor > aa_nn
Values at successive time steps of the components .
void init_bhole()
Sets the values of the fields to :
Sym_tensor gamt_point
Time derivative of the 3-metric tilde.
double regularise_one()
Corrects the shift in the innermost shell, so that it remains and that equals zero on the horizon.
Evolution_std< Scalar > aa_quad_evol
Values at successive time steps of the components .
double omega
Angular velocity in LORENE's units.
double regul
Intensity of the correction on the shift vector.
virtual const Scalar & n_comp() const
Lapse function at the current time step jtime.
Scalar trK
Trace of the extrinsic curvature.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator<<).
const Scalar darea_hor() const
Element of area of the horizon.
Metric met_gamt
3 metric tilde
Evolution_std< Vector > dpsi_evol
Values at successive time steps of the covariant derivative of the conformal factor .
Evolution_std< Scalar > psi_auto_evol
Values at successive time steps of the conformal factor .
double ang_mom_hor() const
Angular momentum (modulo).
Evolution_std< Sym_tensor > aa_comp_evol
Values at successive time steps of the components of the conformal representation of the traceless p...
void set_psi(const Scalar &psi_in)
Sets the conformal factor relating the physical metric to the conformal one: .
virtual const Vector & beta_comp() const
Shift function at the current time step jtime.
void init_met_trK()
Sets the 3-metric tilde to the flat metric and gamt_point, trK and trK_point to zero.
Scalar trK_point
Time derivative of the trace of the extrinsic curvature.
Evolution_std< Vector > beta_auto_evol
Values at successive time steps of the shift function .
virtual const Vector & dpsi() const
Covariant derivative with respect to the flat metric of the conformal factor at the current time ste...
void set_nn(const Scalar &nn_in)
Sets the lapse.
double boost_z
Boost velocity in z-direction.
virtual const Scalar & psi_auto() const
Conformal factor at the current time step jtime.
double ang_mom_adm() const
ADM angular Momentum.
virtual const Sym_tensor & aa_comp() const
Conformal representation of the traceless part of the extrinsic curvature: Returns the value at the ...
Scalar decouple
Function used to construct from the total .
void update_aa()
Conformal representation of the traceless part of the extrinsic curvature: .
virtual const Scalar & n_auto() const
Lapse function at the current time step jtime.
virtual ~Isol_hor()
Destructor.
double radius
Radius of the horizon in LORENE's units.
virtual const Sym_tensor & aa_auto() const
Conformal representation of the traceless part of the extrinsic curvature: Returns the value at the ...
double mass_hor() const
Mass computed at the horizon.
void init_bhole_seul()
Initiates for a single black hole.
Evolution_std< Sym_tensor > aa_auto_evol
Values at successive time steps of the components of the conformal representation of the traceless p...
void set_gamt(const Metric &gam_tilde)
Sets the conformal metric to gam_tilde.
double radius_hor() const
Radius of the horizon.
Map_af & mp
Affine mapping.
virtual const Scalar & psi_comp() const
Conformal factor at the current time step jtime.
double axi_break() const
Breaking of the axial symmetry on the horizon.
Evolution_std< Scalar > n_auto_evol
Values at successive time steps of the lapse function .
Evolution_std< Scalar > n_comp_evol
Values at successive time steps of the lapse function .
Evolution_std< Vector > dn_evol
Values at successive time steps of the covariant derivative of the lapse with respect to the flat met...
Isol_hor(Map_af &mpi, int depth_in=3)
Standard constructor.
double boost_x
Boost velocity in x-direction.
void operator=(const Isol_hor &)
Assignment to another Isol_hor.
Evolution_std< Vector > beta_comp_evol
Values at successive time steps of the shift function .
double area_hor() const
Area of the horizon.
virtual const Vector & beta_auto() const
Shift function at the current time step jtime.
virtual void sauve(FILE *fich, bool partial_save) const
Total or partial saves in a binary file.
double omega_hor() const
Orbital velocity.
Evolution_std< Scalar > psi_comp_evol
Values at successive time steps of the lapse function .
Flat metric for tensor calculation.
Metric for tensor calculation.
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
virtual const Sym_tensor & cov() const
Read-only access to the covariant representation.
const double & get_double(int position=0) const
Returns the reference of a double stored in the list.
const Scalar & get_scalar(int position=0) const
Returns the reference of a Scalar stored in the list.
void add_scalar(const Scalar &ti, int position=0)
Adds the address of a new Scalar to the list.
Tensor field of valence 0 (or component of a tensorial field).
void raccord(int n)
Performs the matching of the nucleus with respect to the first shell.
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this.
void div_rsint_dzpuis(int ced_mult_r)
Division by but with the output flag dzpuis set to ced_mult_r .
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Tbl & set_domain(int l)
Read/write of the value in a given domain.
virtual void annule(int l_min, int l_max)
Sets the Scalar to zero in several domains.
Valeur & set_spectral_va()
Returns va (read/write version).
const Valeur & get_spectral_va() const
Returns va (read only version).
void mult_rsint()
Multiplication by everywhere; dzpuis is not changed.
void set_dzpuis(int)
Modifies the dzpuis flag.
void div_rsint()
Division by everywhere; dzpuis is not changed.
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
void set_outer_boundary(int l, double x)
Sets the value of the Scalar at the outer boundary of a given domain.
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values of the Scalar in the co...
void import(const Scalar &ci)
Assignment to another Scalar defined on a different mapping.
Class intended to describe valence-2 symmetric tensors.
Sym_tensor derive_lie(const Vector &v) const
Computes the Lie derivative of this with respect to some vector field v.
Evolution_std< Sym_tensor > hh_evol
Values at successive time steps of the components .
virtual const Sym_tensor & gam_dd() const
Induced metric (covariant components ) at the current time step (jtime ).
virtual void sauve(FILE *fich, bool partial_save) const
Total or partial saves in a binary file.
virtual const Scalar & nn() const
Lapse function N at the current time step (jtime ).
virtual const Sym_tensor & k_uu() const
Extrinsic curvature tensor (contravariant components ) at the current time step (jtime ).
Evolution_std< Sym_tensor > hata_evol
Values at successive time steps of the components .
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).
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator<<).
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).
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...
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 ).
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.
Scalar * p_ln_psi
Pointer on the logarithm of at the current time step (jtime).
Scalar * p_psi4
Pointer on the factor 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.
Metric * p_gamma
Pointer on the induced metric at the current time step (jtime).
Evolution_std< Sym_tensor > k_dd_evol
Values at successive time steps of the covariant components of the extrinsic curvature tensor .
Evolution_full< Tbl > adm_mass_evol
ADM mass at each time step, since the creation of the slice.
Evolution_std< Sym_tensor > gam_dd_evol
Values at successive time steps of the covariant components of the induced metric .
Evolution_std< Sym_tensor > gam_uu_evol
Values at successive time steps of the contravariant components of the induced metric .
virtual const Vector & beta() const
shift vector at the current time step (jtime )
const Metric & gam() const
Induced metric 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 .
const Base_val & get_base() const
Return the bases for spectral expansions (member base ).
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base ).
Tensor field of valence 1.
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Sym_tensor ope_killing_conf(const Metric &gam) const
Computes the conformal Killing operator associated with a given metric.
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Scalar & set(int)
Read/write access to a component.
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Cmp pow(const Cmp &, int)
Power .
#define R_CHEBU
base de Chebychev ordinaire (fin), dev. en 1/r
#define T_COSSIN_SI
sin impair-cos pair alternes, sin pour m=0
#define T_COSSIN_CI
cos impair-sin pair alternes, cos pour m=0
#define P_COSSIN
dev. standart
#define R_CHEBPIM_P
Cheb. pair-impair suivant m, pair pour m=0.
#define R_CHEB
base de Chebychev ordinaire (fin)
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
double zerosec_b(double(*f)(double, const Param &), const Param &par, double a, double b, double precis, int nitermax, int &niter)
Finding the zero a function on a bounded domain.
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
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 .
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.
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Coord phi
coordinate centered on the grid