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
00029
00030
00031 char CTS_gen[] = "$Header: /cvsroot/Lorene/C++/Source/Isol_hor/CTS_gen.C,v 1.5 2008/08/19 06:42:00 j_novak Exp $" ;
00032
00033
00034
00035
00036
00037 #include "headcpp.h"
00038
00039
00040 #include <stdlib.h>
00041 #include <assert.h>
00042
00043
00044 #include "time_slice.h"
00045 #include "isol_hor.h"
00046 #include "metric.h"
00047 #include "evolution.h"
00048 #include "unites.h"
00049 #include "graphique.h"
00050 #include "utilitaires.h"
00051 #include "param.h"
00052 #include "vector.h"
00053
00054 void Isol_hor::init_data_CTS_gen(int, double, int, int,
00055 int solve_lapse, int solve_psi, int solve_shift,
00056 double precis, double relax_nn ,
00057 double relax_psi, double relax_beta ,
00058 int niter, double a, double ) {
00059
00060 using namespace Unites ;
00061
00062
00063
00064
00065 double lambda = 0. ;
00066
00067 double ttime = the_time[jtime] ;
00068 const Base_vect& triad = *(ff.get_triad()) ;
00069
00070
00071 ofstream conv("resconv.d") ;
00072 ofstream kss_out("kss.d") ;
00073 conv << " # diff_nn diff_psi diff_beta " << endl ;
00074
00075
00076 Scalar nntilde_j = nn() * pow(psi(), a) ;
00077 nntilde_j.std_spectral_base() ;
00078 Scalar psi_j = psi() ;
00079 Vector beta_j = beta() ;
00080
00081
00082
00083
00084
00085 for (int mer=0; mer<niter; mer++) {
00086
00087
00088
00089
00090 Scalar temp_scal_0 (mp) ;
00091 Scalar temp_scal_2 (mp) ;
00092 Scalar temp_scal_3 (mp) ;
00093 Scalar temp_scal_4 (mp) ;
00094
00095 Vector temp_vect_2 (mp, CON, triad) ;
00096 Vector temp_vect_3 (mp, CON, triad) ;
00097 Vector temp_vect_4 (mp, CON, triad) ;
00098
00099 Scalar psi_j_a = pow(psi_j, a) ;
00100 psi_j_a.std_spectral_base() ;
00101 Scalar psi_j_ma = pow(psi_j, -a) ;
00102 psi_j_ma.std_spectral_base() ;
00103 Vector beta_j_d = beta_j.down(0, met_gamt) ;
00104 Scalar nnpsia_j = nntilde_j * psi_j_a /( psi_j*psi_j*psi_j*psi_j*psi_j*psi_j) ;
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115 cout << " relax_nn = " << relax_nn << endl ;
00116 cout << " relax_psi = " << relax_psi << endl ;
00117 cout << " relax_beta = " << relax_beta << endl ;
00118
00119
00120
00121
00122
00123
00124
00125
00126 Scalar sour_psi (mp) ;
00127
00128
00129
00130 temp_scal_3 = 1./8. * met_gamt.ricci_scal() * psi_j ;
00131 temp_scal_3.inc_dzpuis() ;
00132 temp_scal_4 = 1./12. * trK*trK * psi_j*psi_j*psi_j*psi_j*psi_j
00133 - 1./32. * psi_j*psi_j*psi_j*psi_j*psi_j * psi_j_ma*psi_j_ma / (nntilde_j*nntilde_j) *
00134 contract(beta_j_d.ope_killing_conf(met_gamt), 0, 1, beta_j.ope_killing_conf(met_gamt), 0, 1) ;
00135
00136 sour_psi = temp_scal_3 + temp_scal_4 ;
00137
00138
00139
00140
00141 temp_scal_3 = - contract (hh(), 0, 1, psi_j.derive_cov(ff).derive_cov(ff), 0, 1) ;
00142 temp_scal_3.inc_dzpuis() ;
00143 temp_scal_4 = - contract (hdirac(), 0, psi_j.derive_cov(ff), 0) ;
00144
00145 sour_psi += temp_scal_3 + temp_scal_4 ;
00146
00147 sour_psi.annule_domain(0) ;
00148
00149
00150
00151
00152 Scalar sour_nntilde (mp) ;
00153
00154
00155
00156 temp_scal_2 = - psi_j*psi_j*psi_j*psi_j * psi_j_ma * trK_point ;
00157 temp_scal_2.inc_dzpuis(2) ;
00158
00159 temp_scal_3 = psi_j*psi_j*psi_j*psi_j * psi_j_ma * contract (beta_j, 0, trK.derive_cov(ff), 0)
00160 - a/8. * nntilde_j * met_gamt.ricci_scal() ;
00161 temp_scal_3.inc_dzpuis() ;
00162
00163 temp_scal_4 = nntilde_j * (4-a)/12. * psi_j*psi_j*psi_j*psi_j * trK*trK
00164 - 2 * (a+1) * contract(psi_j.derive_cov(ff), 0, nntilde_j.derive_con(met_gamt), 0) / psi_j
00165 - a*(a+1) * nntilde_j * contract(psi_j.derive_cov(met_gamt), 0, psi_j.derive_con(met_gamt), 0) / (psi_j * psi_j)
00166 + (a+8)/32. * psi_j*psi_j*psi_j*psi_j * psi_j_ma*psi_j_ma/nntilde_j
00167 * contract(beta_j_d.ope_killing_conf(met_gamt), 0, 1, beta_j.ope_killing_conf(met_gamt), 0, 1) ;
00168
00169 sour_nntilde = temp_scal_2 + temp_scal_3 + temp_scal_4 ;
00170
00171
00172
00173 temp_scal_3 = - contract (hh(), 0, 1, nntilde_j.derive_cov(ff).derive_cov(ff), 0, 1) ;
00174 temp_scal_3.inc_dzpuis() ;
00175
00176 temp_scal_4 = - contract (hdirac(), 0, nntilde_j.derive_cov(ff), 0) ;
00177
00178 sour_nntilde += temp_scal_3 + temp_scal_4 ;
00179 sour_nntilde.annule_domain(0) ;
00180
00181
00182
00183
00184 Vector sour_beta(mp, CON, triad) ;
00185
00186
00187
00188 temp_vect_3 = 4./3. * psi_j_a * nntilde_j * trK.derive_con(met_gamt)
00189 - contract(met_gamt.ricci(), 1, beta_j, 0).up(0, met_gamt) ;
00190 temp_vect_3.inc_dzpuis() ;
00191
00192 temp_vect_4 = contract(beta_j.ope_killing_conf(met_gamt), 1, nnpsia_j.derive_cov(ff), 0) / nnpsia_j ;
00193
00194 sour_beta = temp_vect_3 + temp_vect_4 ;
00195
00196
00197
00198
00199 temp_vect_3 = (lambda - 1./3.) * contract (beta_j.derive_cov(ff).derive_con(ff), 0, 1)
00200 - contract(contract(met_gamt.connect().get_delta(), 1, beta_j, 0).derive_con(ff), 1, 2)
00201 - contract(hh(), 0, 1, beta_j.derive_cov(met_gamt).derive_cov(ff), 1, 2)
00202 - 1./3. * contract (hh(), 1, beta_j.divergence(ff).derive_cov(ff), 0) ;
00203 temp_vect_3.inc_dzpuis() ;
00204
00205 temp_vect_4 = - contract(hdirac(), 0, beta_j.derive_cov(met_gamt), 1)
00206 - contract(met_gamt.connect().get_delta(), 1, 2, beta_j.derive_con(met_gamt), 0, 1) ;
00207
00208 sour_beta += temp_vect_3 + temp_vect_4 ;
00209
00210
00211
00212
00213
00214
00215
00216
00217 Scalar tmp = psi_j * psi_j * psi_j * trK
00218 - psi_j * met_gamt.radial_vect().divergence(met_gamt)
00219 - psi_j*psi_j*psi_j*psi_j_ma/(2. * nntilde_j)
00220 * contract( beta_j_d.ope_killing_conf(met_gamt), 0, 1, met_gamt.radial_vect()
00221 * met_gamt.radial_vect(), 0, 1)
00222 - 1./3. * psi_j*psi_j*psi_j * trK * contract( met_gamt.cov(), 0, 1, met_gamt.radial_vect()
00223 * met_gamt.radial_vect(), 0, 1)
00224 - 4 * ( met_gamt.radial_vect()(2) * psi_j.derive_cov(ff)(2)
00225 + met_gamt.radial_vect()(3) * psi_j.derive_cov(ff)(3) ) ;
00226
00227 tmp = tmp / (4 * met_gamt.radial_vect()(1)) ;
00228
00229
00230 Valeur psi_bound (mp.get_mg()->get_angu() ) ;
00231
00232 int nnp = mp.get_mg()->get_np(1) ;
00233 int nnt = mp.get_mg()->get_nt(1) ;
00234
00235 psi_bound = 1 ;
00236
00237 for (int k=0 ; k<nnp ; k++)
00238 for (int j=0 ; j<nnt ; j++)
00239 psi_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00240
00241 psi_bound.std_base_scal() ;
00242
00243
00244
00245
00246
00247 const Coord& y = mp.y ;
00248 const Coord& x = mp.x ;
00249
00250 Scalar xx(mp) ;
00251 xx = x ;
00252 xx.std_spectral_base() ;
00253
00254 Scalar yy(mp) ;
00255 yy = y ;
00256 yy.std_spectral_base() ;
00257
00258 Vector tmp_vect = nntilde_j * psi_j_a/(psi_j * psi_j) * met_gamt.radial_vect() ;
00259 tmp_vect.change_triad(mp.get_bvect_cart() ) ;
00260
00261
00262
00263 Valeur lim_x (mp.get_mg()->get_angu()) ;
00264
00265 lim_x = 1 ;
00266
00267 for (int k=0 ; k<nnp ; k++)
00268 for (int j=0 ; j<nnt ; j++)
00269 lim_x.set(0, k, j, 0) = omega * yy.val_grid_point(1, k, j, 0)
00270 + tmp_vect(1).val_grid_point(1, k, j, 0) ;
00271
00272 lim_x.set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
00273
00274
00275
00276
00277 Valeur lim_y (mp.get_mg()->get_angu()) ;
00278
00279 lim_y = 1 ;
00280
00281 for (int k=0 ; k<nnp ; k++)
00282 for (int j=0 ; j<nnt ; j++)
00283 lim_y.set(0, k, j, 0) = - omega * xx.val_grid_point(1, k, j, 0)
00284 + tmp_vect(2).val_grid_point(1, k, j, 0) ;
00285
00286 lim_y.set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
00287
00288
00289
00290 Valeur lim_z (mp.get_mg()->get_angu()) ;
00291
00292 lim_z = 1 ;
00293
00294 for (int k=0 ; k<nnp ; k++)
00295 for (int j=0 ; j<nnt ; j++)
00296 lim_z.set(0, k, j, 0) = tmp_vect(3).val_grid_point(1, k, j, 0) ;
00297
00298 lim_z.set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
00299
00300
00301
00302
00303
00304
00305
00306
00307 Vector stilde_j = met_gamt.radial_vect() ;
00308 Scalar btilde_j = contract (beta_j_d, 0, stilde_j, 0) ;
00309
00310
00311
00312 Scalar hh_tilde = contract(stilde_j.derive_cov(met_gamt), 0, 1) ;
00313 hh_tilde.dec_dzpuis(2) ;
00314
00315
00316 tmp_vect = btilde_j * stilde_j ;
00317 Vector VV_j = btilde_j * stilde_j - beta_j ;
00318
00319
00320 Scalar accel = 2 * contract( VV_j, 0, contract( stilde_j, 0, stilde_j.down(0,
00321 met_gamt).derive_cov(met_gamt), 1), 0) ;
00322
00323
00324
00325 Sym_tensor qq_spher = met_gamt.con() - stilde_j * stilde_j ;
00326 Scalar div_VV = contract( qq_spher.down(0, met_gamt), 0, 1, VV_j.derive_cov(met_gamt), 0, 1) ;
00327
00328 Vector angular (mp, CON, mp.get_bvect_cart()) ;
00329 Scalar yya (mp) ;
00330 yya = mp.ya ;
00331 Scalar xxa (mp) ;
00332 xxa = mp.xa ;
00333
00334 angular.set(1) = - omega * yya ;
00335 angular.set(2) = omega * xxa ;
00336 angular.set(3).annule_hard() ;
00337
00338
00339 angular.set(1).set_spectral_va()
00340 .set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
00341 angular.set(2).set_spectral_va()
00342 .set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
00343 angular.set(3).set_spectral_va()
00344 .set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
00345
00346 angular.change_triad(mp.get_bvect_spher()) ;
00347
00348
00349
00350 Scalar corr_nn_kappa (mp) ;
00351 corr_nn_kappa = nntilde_j * psi_j_a - psi_j*psi_j*contract (beta_j_d, 0, met_gamt.radial_vect(), 0) ;
00352
00353 corr_nn_kappa = sqrt(sqrt (corr_nn_kappa*corr_nn_kappa)) ;
00354 corr_nn_kappa.std_spectral_base() ;
00355
00356
00357 Scalar kss = - kappa_hor() * psi_j_ma / nntilde_j ;
00358 kss.inc_dzpuis(2) ;
00359 kss += contract(stilde_j, 0, nntilde_j.derive_cov(ff), 0)/ (psi_j*psi_j*nntilde_j)
00360 + a * contract(stilde_j, 0, psi_j.derive_cov(ff), 0) / (psi_j*psi_j*psi_j)
00361 + 0.1 * contract (nntilde_j.derive_cov(ff), 0, met_gamt.radial_vect(), 0) * corr_nn_kappa ;
00362
00363
00364
00365
00366
00367 double rho = 5. ;
00368 Scalar beta_r_corr = (rho - 1) * btilde_j * hh_tilde;
00369 beta_r_corr.inc_dzpuis(2) ;
00370 Scalar nn_KK = nntilde_j * psi_j_a * trK ;
00371 nn_KK.inc_dzpuis(2) ;
00372
00373 beta_r_corr.set_dzpuis(2) ;
00374 nn_KK.set_dzpuis(2) ;
00375 accel.set_dzpuis(2) ;
00376 div_VV.set_dzpuis(2) ;
00377
00378 Scalar b_tilde_new (mp) ;
00379 b_tilde_new = 2 * contract(stilde_j, 0, btilde_j.derive_cov(ff), 0)
00380 + beta_r_corr
00381 - 3 * nntilde_j * psi_j_a * kss + nn_KK + accel + div_VV ;
00382
00383 b_tilde_new = b_tilde_new / (hh_tilde * rho) ;
00384
00385 b_tilde_new.dec_dzpuis(2) ;
00386
00387 tmp_vect.set(1) = b_tilde_new * stilde_j(1) + angular(1) ;
00388 tmp_vect.set(2) = b_tilde_new * stilde_j(2) + angular(2) ;
00389 tmp_vect.set(3) = b_tilde_new * stilde_j(3) + angular(3) ;
00390
00391 tmp_vect.change_triad(mp.get_bvect_cart() ) ;
00392
00393
00394 Valeur lim_x_AEI (mp.get_mg()->get_angu()) ;
00395
00396 lim_x_AEI = 1 ;
00397
00398 for (int k=0 ; k<nnp ; k++)
00399 for (int j=0 ; j<nnt ; j++)
00400 lim_x_AEI.set(0, k, j, 0) = tmp_vect(1).val_grid_point(1, k, j, 0) ;
00401
00402 lim_x_AEI.set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
00403
00404
00405
00406 Valeur lim_y_AEI (mp.get_mg()->get_angu()) ;
00407
00408 lim_y_AEI = 1 ;
00409
00410 for (int k=0 ; k<nnp ; k++)
00411 for (int j=0 ; j<nnt ; j++)
00412 lim_y_AEI.set(0, k, j, 0) = tmp_vect(2).val_grid_point(1, k, j, 0) ;
00413
00414 lim_y_AEI.set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
00415
00416
00417 Valeur lim_z_AEI (mp.get_mg()->get_angu()) ;
00418
00419 lim_z_AEI = 1 ;
00420
00421 for (int k=0 ; k<nnp ; k++)
00422 for (int j=0 ; j<nnt ; j++)
00423 lim_z_AEI.set(0, k, j, 0) = tmp_vect(3).val_grid_point(1, k, j, 0) ;
00424
00425 lim_z_AEI.set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436
00437
00438
00439
00440 Scalar darea = psi_j* psi_j* psi_j* psi_j
00441 * sqrt( met_gamt.cov()(2,2) * met_gamt.cov()(3,3) -
00442 met_gamt.cov()(2,3) * met_gamt.cov()(2,3) ) ;
00443
00444 darea.std_spectral_base() ;
00445
00446 Scalar area_int (darea) ;
00447 area_int.raccord(1) ;
00448 double area = mp.integrale_surface(area_int, radius + 1e-15) ;
00449
00450
00451 double radius = area / (4. * M_PI);
00452 radius = pow(radius, 1./2.) ;
00453
00454
00455 Vector phi (mp, CON, *(ff.get_triad()) ) ;
00456 tmp = 1 ;
00457 tmp.std_spectral_base() ;
00458 tmp.mult_rsint() ;
00459 phi.set(1) = 0. ;
00460 phi.set(2) = 0. ;
00461 phi.set(3) = tmp ;
00462
00463 Scalar kksphi = psi_j*psi_j*psi_j_ma/(2.*nntilde_j) *
00464 contract( beta_j_d.ope_killing_conf(met_gamt), 0, 1, met_gamt.radial_vect()* phi, 0, 1) +
00465 1./3. * trK * psi_j*psi_j *
00466 contract( met_gamt.cov(), 0, 1, met_gamt.radial_vect()* phi, 0, 1) ;
00467
00468 kksphi = kksphi / (8. * M_PI) ;
00469 Scalar kkspsi_int = kksphi * darea ;
00470
00471 double ang_mom = mp.integrale_surface(kkspsi_int, radius + 1e-15) ;
00472
00473
00474 double kappa_kerr = (pow( radius, 4) - 4 * pow( ang_mom, 2)) / ( 2 * pow( radius, 3)
00475 * sqrt( pow( radius, 4) + 4 * pow( ang_mom, 2) ) ) ;
00476
00477
00478
00479 rho = 5. ;
00480
00481 Scalar kappa (mp) ;
00482 if (mer < 0)
00483 kappa = kappa_kerr ;
00484 else
00485 kappa = 0.15 ;
00486 kappa.std_spectral_base() ;
00487 kappa.inc_dzpuis(2) ;
00488
00489
00490
00491 tmp = - a / psi_j * nntilde_j
00492 * contract(met_gamt.radial_vect(), 0, psi_j.derive_cov(met_gamt), 0)
00493 + 1./2. * psi_j * psi_j * psi_j_ma
00494 * contract( beta_j_d.ope_killing_conf(met_gamt), 0, 1, met_gamt.radial_vect()
00495 * met_gamt.radial_vect(), 0, 1)
00496 + 1./3. * nntilde_j * psi_j * psi_j * trK
00497 * contract( met_gamt.cov(), 0, 1, met_gamt.radial_vect()* met_gamt.radial_vect(), 0, 1)
00498 + psi_j * psi_j * psi_j_ma * kappa
00499 - met_gamt.radial_vect()(2) * nntilde_j.derive_cov(ff)(2)
00500 - met_gamt.radial_vect()(3) * nntilde_j.derive_cov(ff)(3)
00501 + ( rho - 1 ) * met_gamt.radial_vect()(1) * nntilde_j.derive_cov(ff)(1)
00502 + 0. * contract (nntilde_j.derive_cov(ff), 0, met_gamt.radial_vect(), 0) * corr_nn_kappa ;
00503
00504 tmp = tmp / (rho * met_gamt.radial_vect()(1)) ;
00505
00506
00507 Valeur nn_bound_kappa (mp.get_mg()->get_angu() ) ;
00508
00509 nn_bound_kappa = 1 ;
00510
00511 for (int k=0 ; k<nnp ; k++)
00512 for (int j=0 ; j<nnt ; j++)
00513 nn_bound_kappa.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00514
00515 nn_bound_kappa.std_base_scal() ;
00516
00517
00518
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549 rho = 0. ;
00550 tmp = 0.2 * psi_j_ma ;
00551 tmp = (tmp + rho * nntilde_j)/(1 + rho) ;
00552 tmp = tmp - 1 ;
00553
00554
00555 Valeur nn_bound_const (mp.get_mg()->get_angu() ) ;
00556
00557 nn_bound_const = 1 ;
00558
00559 for (int k=0 ; k<nnp ; k++)
00560 for (int j=0 ; j<nnt ; j++)
00561 nn_bound_const.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00562
00563 nn_bound_const.std_base_scal() ;
00564
00565 boundary_nn_Neu_kk(0) ;
00566
00567
00568
00569 rho = 5. ;
00570 tmp = btilde_j * psi_j*psi_j*psi_j_ma ;
00571 tmp = (tmp + rho * nntilde_j)/(1 + rho) ;
00572 tmp = tmp - 1 ;
00573
00574
00575 Valeur nn_bound_AEI (mp.get_mg()->get_angu() ) ;
00576
00577 nn_bound_AEI = 1 ;
00578
00579 for (int k=0 ; k<nnp ; k++)
00580 for (int j=0 ; j<nnt ; j++)
00581 nn_bound_AEI.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00582
00583 nn_bound_AEI.std_base_scal() ;
00584
00585
00586
00587
00588
00589
00590
00591 Scalar psi_jp1(mp) ;
00592 Scalar nntilde_jp1(mp) ;
00593 Vector beta_jp1(beta_j) ;
00594
00595
00596 if (solve_psi == 1) {
00597 psi_jp1 = sour_psi.poisson_neumann(psi_bound, 0) + 1. ;
00598
00599
00600 maxabs(psi_jp1.laplacian() - sour_psi,
00601 "Absolute error in the resolution of the equation for Psi") ;
00602
00603 psi_jp1 = relax_psi * psi_jp1 + (1 - relax_psi) * psi_j ;
00604 }
00605
00606
00607 if (solve_lapse == 1) {
00608
00609
00610
00611 nntilde_jp1 = sour_nntilde.poisson_dirichlet(nn_bound_AEI, 0) + 1. ;
00612
00613
00614 maxabs(nntilde_jp1.laplacian() - sour_nntilde,
00615 "Absolute error in the resolution of the equation for nntilde") ;
00616
00617 nntilde_jp1 = relax_nn * nntilde_jp1 + (1 - relax_nn) * nntilde_j ;
00618 }
00619
00620 if (solve_shift == 1) {
00621 double precision = 1e-8 ;
00622
00623
00624 poisson_vect_boundary(lambda, sour_beta, beta_jp1, lim_x_AEI,
00625 lim_y_AEI, lim_z_AEI, 0, precision, 20) ;
00626
00627
00628
00629 sour_beta.dec_dzpuis() ;
00630 maxabs(beta_jp1.derive_con(ff).divergence(ff)
00631 + lambda * beta_jp1.divergence(ff)
00632 .derive_con(ff) - sour_beta,
00633 "Absolute error in the resolution of the equation for beta") ;
00634
00635 cout << endl ;
00636
00637
00638 beta_jp1 = relax_beta * beta_jp1 + (1 - relax_beta) * beta_j ;
00639 }
00640
00641
00642
00643
00644
00645
00646
00647 double diff_nn, diff_psi, diff_beta ;
00648 diff_nn = 1.e-16 ;
00649 diff_psi = 1.e-16 ;
00650 diff_beta = 1.e-16 ;
00651 if (solve_lapse == 1)
00652 diff_nn = max( diffrel(nntilde_j, nntilde_jp1) ) ;
00653 if (solve_psi == 1)
00654 diff_psi = max( diffrel(psi_j, psi_jp1) ) ;
00655 if (solve_shift == 1)
00656 diff_beta = max( maxabs(beta_jp1 - beta_j) ) ;
00657
00658 cout << "step = " << mer << " : diff_psi = " << diff_psi
00659 << " diff_nn = " << diff_nn
00660 << " diff_beta = " << diff_beta << endl ;
00661 cout << "----------------------------------------------" << endl ;
00662 if ((diff_psi<precis) && (diff_nn<precis) && (diff_beta<precis))
00663 break ;
00664
00665 if (mer>0) {conv << mer << " " << log10(diff_nn) << " " << log10(diff_psi)
00666 << " " << log10(diff_beta) << endl ; } ;
00667
00668
00669
00670
00671
00672
00673
00674 psi_j = psi_jp1 ;
00675 nntilde_j = nntilde_jp1 ;
00676 beta_j = beta_jp1 ;
00677
00678
00679
00680
00681
00682
00683
00684 Scalar kkss ( psi_j_ma/(2. * nntilde_j )
00685 * contract( beta_j_d.ope_killing_conf(met_gamt), 0, 1, met_gamt.radial_vect()
00686 * met_gamt.radial_vect(), 0, 1)
00687 + 1./3. * trK * contract( met_gamt.cov(), 0, 1, met_gamt.radial_vect()
00688 * met_gamt.radial_vect(), 0, 1) ) ;
00689
00690 double max_kss = kkss.val_grid_point(1, 0, 0, 0) ;
00691 double min_kss = kkss.val_grid_point(1, 0, 0, 0) ;
00692
00693 hh_tilde = contract(met_gamt.radial_vect().derive_cov(met_gamt), 0, 1) ;
00694 double max_hh_tilde = hh_tilde.val_grid_point(1, 0, 0, 0) ;
00695 double min_hh_tilde = hh_tilde.val_grid_point(1, 0, 0, 0) ;
00696
00697
00698 for (int k=0 ; k<nnp ; k++)
00699 for (int j=0 ; j<nnt ; j++){
00700 if (kkss.val_grid_point(1, k, j, 0) > max_kss)
00701 max_kss = kkss.val_grid_point(1, k, j, 0) ;
00702 if (kkss.val_grid_point(1, k, j, 0) < min_kss)
00703 min_kss = kkss.val_grid_point(1, k, j, 0) ;
00704
00705 if (hh_tilde.val_grid_point(1, k, j, 0) > max_hh_tilde)
00706 max_hh_tilde = hh_tilde.val_grid_point(1, k, j, 0) ;
00707 if (hh_tilde.val_grid_point(1, k, j, 0) < min_hh_tilde)
00708 min_hh_tilde = hh_tilde.val_grid_point(1, k, j, 0) ;
00709 }
00710 kss_out << mer << " " << max_kss << " " << min_kss
00711 << " " << -1 * max_hh_tilde << " " << -1 * min_hh_tilde << endl ;
00712
00713
00714 if (solve_psi == 1)
00715 set_psi(psi_j) ;
00716 if (solve_lapse == 1)
00717 n_evol.update(nntilde_j * psi_j_a, jtime, ttime) ;
00718 if (solve_shift == 1)
00719 beta_evol.update(beta_j, jtime, ttime) ;
00720
00721 if (solve_shift == 1)
00722 update_aa() ;
00723
00724 }
00725
00726
00727
00728 Scalar psi_j_a = pow(psi_j, a) ;
00729 psi_j_a.std_spectral_base() ;
00730
00731
00732 if (solve_psi == 1)
00733 set_psi(psi_j) ;
00734 if (solve_lapse == 1)
00735 n_evol.update(nntilde_j * psi_j_a, jtime, ttime) ;
00736 if (solve_shift == 1)
00737 beta_evol.update(beta_j, jtime, ttime) ;
00738
00739 if (solve_shift == 1)
00740 update_aa() ;
00741
00742
00743
00744
00745
00746
00747
00748
00749
00750
00751
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766 conv.close() ;
00767 kss_out.close() ;
00768
00769
00770
00771
00772 }
00773
00774