00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 char time_slice_C[] = "$Header: /cvsroot/Lorene/C++/Source/Time_slice/time_slice.C,v 1.14 2008/12/02 15:02:21 j_novak Exp $" ;
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090 #include <assert.h>
00091
00092
00093 #include "time_slice.h"
00094 #include "utilitaires.h"
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105 Time_slice::Time_slice(const Scalar& lapse_in, const Vector& shift_in,
00106 const Sym_tensor& gamma_in, const Sym_tensor& kk_in,
00107 int depth_in)
00108 : depth(depth_in),
00109 scheme_order(depth_in-1),
00110 jtime(0),
00111 the_time(0., depth_in),
00112 gam_dd_evol(depth_in),
00113 gam_uu_evol(depth_in),
00114 k_dd_evol(depth_in),
00115 k_uu_evol(depth_in),
00116 n_evol(lapse_in, depth_in),
00117 beta_evol(shift_in, depth_in),
00118 trk_evol(depth_in),
00119 adm_mass_evol() {
00120
00121 set_der_0x0() ;
00122
00123 double time_init = the_time[jtime] ;
00124
00125 p_gamma = new Metric(gamma_in) ;
00126
00127 if (gamma_in.get_index_type(0) == COV) {
00128 gam_dd_evol.update(gamma_in, jtime, time_init) ;
00129 }
00130 else {
00131 gam_uu_evol.update(gamma_in, jtime, time_init) ;
00132 }
00133
00134 if (kk_in.get_index_type(0) == COV) {
00135 k_dd_evol.update(kk_in, jtime, time_init) ;
00136 }
00137 else {
00138 k_uu_evol.update(kk_in, jtime, time_init) ;
00139 }
00140
00141 trk_evol.update( kk_in.trace(*p_gamma), jtime, the_time[jtime] ) ;
00142
00143 }
00144
00145
00146
00147
00148 Time_slice::Time_slice(const Scalar& lapse_in, const Vector& shift_in,
00149 const Evolution_std<Sym_tensor>& gamma_in)
00150 : depth(gamma_in.get_size()),
00151 scheme_order(gamma_in.get_size()-1),
00152 jtime(0),
00153 the_time(0., gamma_in.get_size()),
00154 gam_dd_evol( gamma_in.get_size() ),
00155 gam_uu_evol( gamma_in.get_size() ),
00156 k_dd_evol( gamma_in.get_size() ),
00157 k_uu_evol( gamma_in.get_size() ),
00158 n_evol(lapse_in, gamma_in.get_size() ),
00159 beta_evol(shift_in, gamma_in.get_size() ),
00160 trk_evol(gamma_in.get_size() ),
00161 adm_mass_evol() {
00162
00163 cerr <<
00164 "Time_slice constuctor from evolution of gamma not implemented yet !\n" ;
00165 abort() ;
00166
00167 set_der_0x0() ;
00168
00169 }
00170
00171
00172
00173 Time_slice::Time_slice(const Map& mp, const Base_vect& triad, int depth_in)
00174 : depth(depth_in),
00175 scheme_order(depth_in-1),
00176 jtime(0),
00177 the_time(0., depth_in),
00178 gam_dd_evol(depth_in),
00179 gam_uu_evol(depth_in),
00180 k_dd_evol(depth_in),
00181 k_uu_evol(depth_in),
00182 n_evol(depth_in),
00183 beta_evol(depth_in),
00184 trk_evol(depth_in),
00185 adm_mass_evol() {
00186
00187 double time_init = the_time[jtime] ;
00188
00189 const Base_vect_spher* ptriad_s =
00190 dynamic_cast<const Base_vect_spher*>(&triad) ;
00191 bool spher = (ptriad_s != 0x0) ;
00192
00193 if (spher) {
00194 gam_dd_evol.update( mp.flat_met_spher().cov(), jtime, time_init) ;
00195 }
00196 else {
00197 assert( dynamic_cast<const Base_vect_cart*>(&triad) != 0x0) ;
00198 gam_dd_evol.update( mp.flat_met_cart().cov(), jtime, time_init) ;
00199 }
00200
00201
00202
00203 Sym_tensor ktmp(mp, COV, triad) ;
00204 ktmp.set_etat_zero() ;
00205 k_dd_evol.update(ktmp, jtime, time_init) ;
00206
00207
00208 Scalar tmp(mp) ;
00209 tmp.set_etat_one() ;
00210 tmp.std_spectral_base() ;
00211 n_evol.update(tmp, jtime, time_init) ;
00212
00213
00214 Vector btmp(mp, CON, triad) ;
00215 btmp.set_etat_zero() ;
00216 beta_evol.update(btmp, jtime, time_init) ;
00217
00218
00219 tmp.set_etat_zero() ;
00220 trk_evol.update(tmp, jtime, time_init) ;
00221
00222 set_der_0x0() ;
00223 }
00224
00225
00226
00227
00228 Time_slice::Time_slice(const Map& mp, const Base_vect& triad, FILE* fich,
00229 bool partial_read, int depth_in)
00230 : depth(depth_in),
00231 the_time(depth_in),
00232 gam_dd_evol(depth_in),
00233 gam_uu_evol(depth_in),
00234 k_dd_evol(depth_in),
00235 k_uu_evol(depth_in),
00236 n_evol(depth_in),
00237 beta_evol(depth_in),
00238 trk_evol(depth_in),
00239 adm_mass_evol() {
00240
00241
00242
00243
00244 int depth_file ;
00245 fread_be(&depth_file, sizeof(int), 1, fich) ;
00246 if (depth_file != depth_in) {
00247 cout <<
00248 "Time_slice constructor from file: the depth read in file \n"
00249 << " is different from that given in the argument list : \n"
00250 << " depth_file = " << depth_file
00251 << " <-> depth_in " << depth_in << " !" << endl ;
00252 abort() ;
00253 }
00254 fread_be(&scheme_order, sizeof(int), 1, fich) ;
00255 fread_be(&jtime, sizeof(int), 1, fich) ;
00256
00257
00258
00259 int jmin = jtime - depth + 1 ;
00260 int indicator ;
00261 for (int j=jmin; j<=jtime; j++) {
00262 fread_be(&indicator, sizeof(int), 1, fich) ;
00263 if (indicator == 1) {
00264 double xx ;
00265 fread_be(&xx, sizeof(double), 1, fich) ;
00266 the_time.update(xx, j, xx) ;
00267 }
00268 }
00269
00270
00271
00272
00273
00274 for (int j=jmin; j<=jtime; j++) {
00275 fread_be(&indicator, sizeof(int), 1, fich) ;
00276 if (indicator == 1) {
00277 Scalar nn_file(mp, *(mp.get_mg()), fich) ;
00278 n_evol.update(nn_file, j, the_time[j]) ;
00279 }
00280 }
00281
00282
00283 for (int j=jmin; j<=jtime; j++) {
00284 fread_be(&indicator, sizeof(int), 1, fich) ;
00285 if (indicator == 1) {
00286 Vector beta_file(mp, triad, fich) ;
00287 beta_evol.update(beta_file, j, the_time[j]) ;
00288 }
00289 }
00290
00291
00292
00293 if (!partial_read) {
00294 cout <<
00295 "Time_slice constructor from file: the case of full reading\n"
00296 << " is not ready yet !" << endl ;
00297 abort() ;
00298 }
00299
00300 set_der_0x0() ;
00301
00302 }
00303
00304
00305
00306
00307 Time_slice::Time_slice(const Time_slice& tin)
00308 : depth(tin.depth),
00309 scheme_order(tin.scheme_order),
00310 jtime(tin.jtime),
00311 the_time(tin.the_time),
00312 gam_dd_evol(tin.gam_dd_evol),
00313 gam_uu_evol(tin.gam_uu_evol),
00314 k_dd_evol(tin.k_dd_evol),
00315 k_uu_evol(tin.k_uu_evol),
00316 n_evol(tin.n_evol),
00317 beta_evol(tin.beta_evol),
00318 trk_evol(tin.trk_evol),
00319 adm_mass_evol(tin.adm_mass_evol) {
00320
00321 set_der_0x0() ;
00322 }
00323
00324
00325
00326 Time_slice::Time_slice(int depth_in)
00327 : depth(depth_in),
00328 scheme_order(depth_in-1),
00329 jtime(0),
00330 the_time(0., depth_in),
00331 gam_dd_evol(depth_in),
00332 gam_uu_evol(depth_in),
00333 k_dd_evol(depth_in),
00334 k_uu_evol(depth_in),
00335 n_evol(depth_in),
00336 beta_evol(depth_in),
00337 trk_evol(depth_in),
00338 adm_mass_evol() {
00339
00340 set_der_0x0() ;
00341 }
00342
00343
00344
00345
00346
00347
00348
00349
00350 Time_slice::~Time_slice(){
00351
00352 Time_slice::del_deriv() ;
00353
00354 }
00355
00356
00357
00358
00359
00360 void Time_slice::del_deriv() const {
00361
00362 if (p_gamma != 0x0) delete p_gamma ;
00363
00364 set_der_0x0() ;
00365
00366 adm_mass_evol.downdate(jtime) ;
00367 }
00368
00369
00370 void Time_slice::set_der_0x0() const {
00371
00372 p_gamma = 0x0 ;
00373
00374 }
00375
00376
00377
00378
00379
00380
00381 void Time_slice::operator=(const Time_slice& tin) {
00382
00383 depth = tin.depth;
00384 scheme_order = tin.scheme_order ;
00385 jtime = tin.jtime;
00386 the_time = tin.the_time;
00387 gam_dd_evol = tin.gam_dd_evol;
00388 gam_uu_evol = tin.gam_uu_evol;
00389 k_dd_evol = tin.k_dd_evol;
00390 k_uu_evol = tin.k_uu_evol;
00391 n_evol = tin.n_evol;
00392 beta_evol = tin.beta_evol;
00393 trk_evol = tin.trk_evol ;
00394
00395 del_deriv() ;
00396
00397 }
00398
00399
00400
00401
00402
00403
00404 ostream& Time_slice::operator>>(ostream& flux) const {
00405
00406 flux << "\n------------------------------------------------------------\n"
00407 << "Lorene class : " << typeid(*this).name() << '\n' ;
00408 flux << "Number of stored slices : " << depth
00409 << " order of time scheme : " << scheme_order << '\n'
00410 << "Time label t = " << the_time[jtime]
00411 << " index of time step j = " << jtime << '\n' << '\n' ;
00412 if (adm_mass_evol.is_known(jtime)) {
00413 flux << "ADM mass : " << adm_mass() << endl ;
00414 }
00415
00416 flux << "Max. of absolute values of the various fields in each domain: \n" ;
00417 if (gam_dd_evol.is_known(jtime)) {
00418 maxabs( gam_dd_evol[jtime], "gam_{ij}", flux) ;
00419 }
00420 if (gam_uu_evol.is_known(jtime)) {
00421 maxabs( gam_uu_evol[jtime], "gam^{ij}", flux) ;
00422 }
00423 if (k_dd_evol.is_known(jtime)) {
00424 maxabs( k_dd_evol[jtime], "K_{ij}", flux) ;
00425 }
00426 if (k_uu_evol.is_known(jtime)) {
00427 maxabs( k_uu_evol[jtime], "K^{ij}", flux) ;
00428 }
00429 if (n_evol.is_known(jtime)) {
00430 maxabs( n_evol[jtime], "N", flux) ;
00431 }
00432 if (beta_evol.is_known(jtime)) {
00433 maxabs( beta_evol[jtime], "beta^i", flux) ;
00434 }
00435 if (trk_evol.is_known(jtime)) {
00436 maxabs( trk_evol[jtime], "tr K", flux) ;
00437 }
00438
00439 if (p_gamma != 0x0) flux << "Metric gamma is up to date" << endl ;
00440
00441 return flux ;
00442
00443 }
00444
00445
00446 ostream& operator<<(ostream& flux, const Time_slice& sigma) {
00447
00448 sigma >> flux ;
00449 return flux ;
00450
00451 }
00452
00453
00454 void Time_slice::save(const char* rootname) const {
00455
00456
00457
00458 char* filename = new char[ strlen(rootname)+10 ] ;
00459 strcpy(filename, rootname) ;
00460 char nomj[7] ;
00461 sprintf(nomj, "%06d", jtime) ;
00462 strcat(filename, nomj) ;
00463 strcat(filename, ".d") ;
00464
00465 FILE* fich = fopen(filename, "w") ;
00466 if (fich == 0x0) {
00467 cout << "Problem in opening file " << filename << " ! " << endl ;
00468 perror(" reason") ;
00469 abort() ;
00470 }
00471
00472
00473
00474 const Map& map = nn().get_mp() ;
00475 const Mg3d& mgrid = *(map.get_mg()) ;
00476 const Base_vect& triad = *(beta().get_triad()) ;
00477
00478 mgrid.sauve(fich) ;
00479 map.sauve(fich) ;
00480 triad.sauve(fich) ;
00481
00482 fwrite_be(&depth, sizeof(int), 1, fich) ;
00483
00484
00485
00486 bool partial_save = false ;
00487 sauve(fich, partial_save) ;
00488
00489
00490
00491
00492 fclose(fich) ;
00493
00494 delete [] filename ;
00495
00496 }
00497
00498
00499
00500 void Time_slice::sauve(FILE* fich, bool partial_save) const {
00501
00502
00503
00504
00505 fwrite_be(&depth, sizeof(int), 1, fich) ;
00506 fwrite_be(&scheme_order, sizeof(int), 1, fich) ;
00507 fwrite_be(&jtime, sizeof(int), 1, fich) ;
00508
00509
00510
00511 int jmin = jtime - depth + 1 ;
00512 for (int j=jmin; j<=jtime; j++) {
00513 int indicator = (the_time.is_known(j)) ? 1 : 0 ;
00514 fwrite_be(&indicator, sizeof(int), 1, fich) ;
00515 if (indicator == 1) {
00516 double xx = the_time[j] ;
00517 fwrite_be(&xx, sizeof(double), 1, fich) ;
00518 }
00519 }
00520
00521
00522
00523
00524
00525 nn() ;
00526 for (int j=jmin; j<=jtime; j++) {
00527 int indicator = (n_evol.is_known(j)) ? 1 : 0 ;
00528 fwrite_be(&indicator, sizeof(int), 1, fich) ;
00529 if (indicator == 1) n_evol[j].sauve(fich) ;
00530 }
00531
00532
00533 beta() ;
00534 for (int j=jmin; j<=jtime; j++) {
00535 int indicator = (beta_evol.is_known(j)) ? 1 : 0 ;
00536 fwrite_be(&indicator, sizeof(int), 1, fich) ;
00537 if (indicator == 1) beta_evol[j].sauve(fich) ;
00538 }
00539
00540
00541
00542 if (!partial_save) {
00543
00544 cout << "Time_slice::sauve: the full writing is not ready yet !"
00545 << endl ;
00546 abort() ;
00547 }
00548
00549
00550 }
00551
00552
00553
00554
00555
00556
00557
00558
00559
00560
00561
00562
00563
00564