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 char et_bin_bhns_extr_eq_ylm_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_eq_ylm.C,v 1.5 2005/02/28 23:12:28 k_taniguchi Exp $" ;
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 #include <math.h>
00060
00061
00062 #include "et_bin_bhns_extr.h"
00063 #include "etoile.h"
00064 #include "map.h"
00065 #include "coord.h"
00066 #include "param.h"
00067 #include "eos.h"
00068 #include "graphique.h"
00069 #include "utilitaires.h"
00070 #include "unites.h"
00071
00072 void Et_bin_bhns_extr::equil_bhns_extr_ylm_ks(double ent_c,
00073 const double& mass,
00074 const double& sepa,
00075 double* nu_int,
00076 double* beta_int,
00077 double* shift_int,
00078 int mermax,
00079 int mermax_poisson,
00080 double relax_poisson,
00081 double relax_ylm,
00082 int mermax_potvit,
00083 double relax_potvit,
00084 int np_filter,
00085 double thres_adapt,
00086 Tbl& diff) {
00087
00088
00089
00090 using namespace Unites ;
00091
00092 assert( kerrschild == true ) ;
00093
00094
00095
00096
00097 const Mg3d* mg = mp.get_mg() ;
00098 int nz = mg->get_nzone() ;
00099
00100
00101 Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
00102
00103
00104 int l_b = nzet - 1 ;
00105 int i_b = mg->get_nr(l_b) - 1 ;
00106 int k_b ;
00107 int j_b ;
00108
00109
00110 double ent_b = 0 ;
00111
00112
00113
00114
00115 double& diff_ent = diff.set(0) ;
00116 double& diff_vel_pot = diff.set(1) ;
00117 double& diff_logn = diff.set(2) ;
00118 double& diff_beta = diff.set(3) ;
00119 double& diff_shift_x = diff.set(4) ;
00120 double& diff_shift_y = diff.set(5) ;
00121 double& diff_shift_z = diff.set(6) ;
00122
00123
00124
00125
00126 Param par_adapt ;
00127 int nitermax = 100 ;
00128 int niter ;
00129 int adapt_flag = 1 ;
00130
00131
00132 int nz_search = nzet ;
00133
00134 double precis_secant = 1.e-14 ;
00135 double alpha_r ;
00136 double reg_map = 1. ;
00137
00138 Tbl ent_limit(nz) ;
00139
00140 par_adapt.add_int(nitermax, 0) ;
00141
00142 par_adapt.add_int(nzet, 1) ;
00143
00144
00145 par_adapt.add_int(nz_search, 2) ;
00146
00147 par_adapt.add_int(adapt_flag, 3) ;
00148
00149
00150 par_adapt.add_int(j_b, 4) ;
00151
00152 par_adapt.add_int(k_b, 5) ;
00153
00154 par_adapt.add_int_mod(niter, 0) ;
00155
00156 par_adapt.add_double(precis_secant, 0) ;
00157
00158
00159 par_adapt.add_double(reg_map, 1) ;
00160
00161 par_adapt.add_double(alpha_r, 2) ;
00162
00163 par_adapt.add_tbl(ent_limit, 0) ;
00164
00165
00166
00167
00168
00169 double precis_poisson = 1.e-16 ;
00170
00171 Param par_poisson1 ;
00172
00173 par_poisson1.add_int(mermax_poisson, 0) ;
00174 par_poisson1.add_double(relax_poisson, 0) ;
00175 par_poisson1.add_double(precis_poisson, 1) ;
00176 par_poisson1.add_int_mod(niter, 0) ;
00177
00178 par_poisson1.add_cmp_mod( ssjm1_logn ) ;
00179
00180
00181
00182
00183 Param par_poisson2 ;
00184
00185 par_poisson2.add_int(mermax_poisson, 0) ;
00186 par_poisson2.add_double(relax_poisson, 0) ;
00187 par_poisson2.add_double(precis_poisson, 1) ;
00188 par_poisson2.add_int_mod(niter, 0) ;
00189
00190 par_poisson2.add_cmp_mod( ssjm1_beta ) ;
00191
00192
00193
00194
00195 Param par_poisson_vect ;
00196
00197 par_poisson_vect.add_int(mermax_poisson, 0) ;
00198
00199 par_poisson_vect.add_double(relax_poisson, 0) ;
00200 par_poisson_vect.add_double(precis_poisson, 1) ;
00201 par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
00202 par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
00203 par_poisson_vect.add_int_mod(niter, 0) ;
00204
00205
00206
00207
00208
00209 Tenseur pot_ext = logn_comp + pot_centri + loggam ;
00210
00211 Tenseur ent_jm1 = ent ;
00212
00213 Tenseur source(mp) ;
00214
00215
00216 Tenseur source_shift(mp, 1, CON, ref_triad) ;
00217
00218
00219
00220 int l_max = 4 ;
00221
00222 if (l_max > 6) {
00223 cout << "We don't have those ylm programmed yet!!!!" << endl ;
00224 abort() ;
00225 }
00226
00227 int nylm = (l_max+1) * (l_max+1) ;
00228
00229 Cmp** ylmvec = new Cmp* [nylm] ;
00230
00231
00232
00233
00234
00235 for(int mer=0 ; mer<mermax ; mer++ ) {
00236
00237 cout << "-----------------------------------------------" << endl ;
00238 cout << "step: " << mer << endl ;
00239 cout << "diff_ent = " << diff_ent << endl ;
00240
00241
00242
00243
00244
00245
00246 if (irrotational) {
00247 diff_vel_pot = velocity_pot_extr(mass, sepa, mermax_potvit,
00248 precis_poisson, relax_potvit) ;
00249 }
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259 double logn_auto_c = logn_auto()(0, 0, 0, 0) ;
00260 double pot_ext_c = pot_ext()(0, 0, 0, 0) ;
00261
00262
00263
00264
00265
00266 double alpha_r2 = 0 ;
00267 for (int k=0; k<mg->get_np(l_b); k++) {
00268 for (int j=0; j<mg->get_nt(l_b); j++) {
00269
00270 double pot_ext_b = pot_ext()(l_b, k, j, i_b) ;
00271 double logn_auto_b = logn_auto()(l_b, k, j, i_b) ;
00272
00273
00274 double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b)
00275 / ( logn_auto_b - logn_auto_c ) ;
00276
00277 if (alpha_r2_jk > alpha_r2) {
00278 alpha_r2 = alpha_r2_jk ;
00279 k_b = k ;
00280 j_b = j ;
00281 }
00282
00283 }
00284 }
00285
00286 alpha_r = sqrt(alpha_r2) ;
00287
00288 cout << "k_b, j_b, alpha_r: " << k_b << " " << j_b << " "
00289 << alpha_r << endl ;
00290
00291
00292
00293
00294 logn_auto = alpha_r2 * logn_auto ;
00295 logn_auto_regu = alpha_r2 * logn_auto_regu ;
00296 logn_auto_c = logn_auto()(0, 0, 0, 0) ;
00297
00298
00299
00300
00301
00302
00303 (logn_auto().va).smooth(nzet, (logn_auto.set()).va) ;
00304
00305
00306
00307
00308
00309
00310 ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
00311
00312
00313
00314
00315
00316
00317 double dentdx = ent().dsdx()(0, 0, 0, 0) ;
00318 double dentdy = ent().dsdy()(0, 0, 0, 0) ;
00319
00320 cout << "dH/dx|_center = " << dentdx << endl ;
00321 cout << "dH/dy|_center = " << dentdy << endl ;
00322
00323 double dec_fact = 1. ;
00324
00325 Tenseur func_in(mp) ;
00326 func_in.set_etat_qcq() ;
00327 func_in.set() = 1. - dec_fact * (dentdx/ent_c) * mp.x
00328 - dec_fact * (dentdy/ent_c) * mp.y ;
00329 func_in.set().annule(nzet, nz-1) ;
00330 func_in.set_std_base() ;
00331
00332 Tenseur func_ex(mp) ;
00333 func_ex.set_etat_qcq() ;
00334 func_ex.set() = 1. ;
00335 func_ex.set().annule(0, nzet-1) ;
00336 func_ex.set_std_base() ;
00337
00338
00339
00340 ent.set() = ent() * (func_in() + func_ex()) ;
00341
00342 (ent().va).smooth(nzet, (ent.set()).va) ;
00343
00344 double dentdx_new = ent().dsdx()(0, 0, 0, 0) ;
00345 double dentdy_new = ent().dsdy()(0, 0, 0, 0) ;
00346 cout << "dH/dx|_new = " << dentdx_new << endl ;
00347 cout << "dH/dy|_new = " << dentdy_new << endl ;
00348
00349
00350
00351
00352
00353
00354
00355
00356 double dent_eq = ent().dsdr().val_point(ray_eq_pi(),M_PI/2.,M_PI) ;
00357 double dent_pole = ent().dsdr().val_point(ray_pole(),0.,0.) ;
00358 double rap_dent = fabs( dent_eq / dent_pole ) ;
00359 cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
00360
00361 if ( rap_dent < thres_adapt ) {
00362 adapt_flag = 0 ;
00363 cout << "******* FROZEN MAPPING *********" << endl ;
00364 }
00365 else{
00366 adapt_flag = 1 ;
00367
00368 }
00369
00370 ent_limit.set_etat_qcq() ;
00371 for (int l=0; l<nzet; l++) {
00372 ent_limit.set(l) = ent()(l, k_b, j_b, i_b) ;
00373 }
00374
00375 ent_limit.set(nzet-1) = ent_b ;
00376 Map_et mp_prev = mp_et ;
00377
00378 mp.adapt(ent(), par_adapt) ;
00379
00380
00381
00382
00383
00384 mp_prev.homothetie(alpha_r) ;
00385
00386 mp.reevaluate(&mp_prev, nzet+1, ent.set()) ;
00387
00388 double fact_resize = 1. / alpha_r ;
00389 for (int l=nzet; l<nz-1; l++) {
00390 mp_et.resize(l, fact_resize) ;
00391 }
00392 mp_et.resize_extr(fact_resize) ;
00393
00394 double ent_s_max = -1 ;
00395 int k_s_max = -1 ;
00396 int j_s_max = -1 ;
00397 for (int k=0; k<mg->get_np(l_b); k++) {
00398 for (int j=0; j<mg->get_nt(l_b); j++) {
00399 double xx = fabs( ent()(l_b, k, j, i_b) ) ;
00400 if (xx > ent_s_max) {
00401 ent_s_max = xx ;
00402 k_s_max = k ;
00403 j_s_max = j ;
00404 }
00405 }
00406 }
00407 cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
00408 << " and nzet : " << endl ;
00409 cout << " " << ent_s_max << " reached for k = " << k_s_max
00410 << " and j = " << j_s_max << endl ;
00411
00412
00413
00414
00415
00416
00417
00418 int oldindex = -1 ;
00419 for (int l=0; l<=l_max; l++) {
00420 for (int m=0; m<= l; m++) {
00421
00422 int index = l*l + 2*m ;
00423 if(m > 0) index -= 1 ;
00424 if(index != oldindex+1) {
00425 cout << "error!! " << l << " " << m << " " << index
00426 << " " << oldindex << endl ;
00427 abort() ;
00428 }
00429
00430
00431 if ((l+m) % 2 == 0) {
00432 ylmvec[index] = new Cmp (logn_auto()) ;
00433 } else {
00434 ylmvec[index] = new Cmp (shift_auto(2)) ;
00435 }
00436 oldindex = index ;
00437 if (m > 0) {
00438 index += 1 ;
00439 if((l+m) % 2 == 0) {
00440 ylmvec[index] = new Cmp (logn_auto()) ;
00441 } else {
00442 ylmvec[index] = new Cmp (shift_auto(2)) ;
00443 }
00444 oldindex = index ;
00445 }
00446 }
00447 }
00448 if(oldindex+1 != nylm) {
00449 cout << "ERROR!!! " << oldindex << " "<< nylm << endl ;
00450 abort() ;
00451 }
00452
00453
00454
00455
00456 get_ylm(nylm, ylmvec);
00457
00458
00459
00460
00461
00462 equation_of_state() ;
00463
00464
00465
00466
00467
00468
00469 hydro_euler_extr(mass, sepa) ;
00470
00471
00472
00473
00474
00475
00476
00477 const Coord& xx = mp.x ;
00478 const Coord& yy = mp.y ;
00479 const Coord& zz = mp.z ;
00480
00481 Tenseur r_bh(mp) ;
00482 r_bh.set_etat_qcq() ;
00483 r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
00484 r_bh.set_std_base() ;
00485
00486 Tenseur xx_cov(mp, 1, COV, ref_triad) ;
00487 xx_cov.set_etat_qcq() ;
00488 xx_cov.set(0) = xx + sepa ;
00489 xx_cov.set(1) = yy ;
00490 xx_cov.set(2) = zz ;
00491 xx_cov.set_std_base() ;
00492
00493 Tenseur xsr_cov(mp, 1, COV, ref_triad) ;
00494 xsr_cov = xx_cov / r_bh ;
00495 xsr_cov.set_std_base() ;
00496
00497 Tenseur msr(mp) ;
00498 msr = ggrav * mass / r_bh ;
00499 msr.set_std_base() ;
00500
00501 Tenseur lapse_bh(mp) ;
00502 lapse_bh = 1. / sqrt( 1.+2.*msr ) ;
00503 lapse_bh.set_std_base() ;
00504
00505 Tenseur lapse_bh2(mp) ;
00506 lapse_bh2 = 1. / (1.+2.*msr) ;
00507 lapse_bh2.set_std_base() ;
00508
00509 Tenseur ldnu(mp) ;
00510 ldnu = flat_scalar_prod_desal(xsr_cov, d_logn_auto) ;
00511 ldnu.set_std_base() ;
00512
00513 Tenseur ldbeta(mp) ;
00514 ldbeta = flat_scalar_prod_desal(xsr_cov, d_beta_auto) ;
00515 ldbeta.set_std_base() ;
00516
00517 Tenseur lltkij(mp) ;
00518 lltkij.set_etat_qcq() ;
00519 lltkij.set() = 0 ;
00520 lltkij.set_std_base() ;
00521
00522 for (int i=0; i<3; i++)
00523 for (int j=0; j<3; j++)
00524 lltkij.set() += xsr_cov(i) * xsr_cov(j) * tkij_auto(i, j) ;
00525
00526 Tenseur lshift(mp) ;
00527 lshift = flat_scalar_prod_desal(xsr_cov, shift_auto) ;
00528 lshift.set_std_base() ;
00529
00530 Tenseur d_ldnu(mp, 1, COV, ref_triad) ;
00531 d_ldnu = ldnu.gradient() ;
00532 d_ldnu.change_triad(ref_triad) ;
00533
00534 Tenseur ldldnu(mp) ;
00535 ldldnu = flat_scalar_prod_desal(xsr_cov, d_ldnu) ;
00536 ldldnu.set_std_base() ;
00537
00538 Tenseur d_ldbeta(mp, 1, COV, ref_triad) ;
00539 d_ldbeta = ldbeta.gradient() ;
00540 d_ldbeta.change_triad(ref_triad) ;
00541
00542 Tenseur ldldbeta(mp) ;
00543 ldldbeta = flat_scalar_prod_desal(xsr_cov, d_ldbeta) ;
00544 ldldbeta.set_std_base() ;
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554 if (relativistic) {
00555
00556 source = qpig * a_car % (ener_euler + s_euler) + akcar_auto
00557 - flat_scalar_prod_desal(d_logn_auto, d_beta_auto)
00558 + 2.*lapse_bh2 % msr % (ldnu % ldbeta + ldldnu)
00559 + lapse_bh2 % lapse_bh2 % msr % (2.*(ldnu + 4.*msr % ldnu)
00560 - ldbeta) / r_bh
00561 - (4.*a_car % lapse_bh2 % lapse_bh2 % msr / 3. / nnn / r_bh)
00562 * (2.+3.*msr) * (3.+4.*msr) * lltkij
00563 + (2.*a_car % lapse_bh2 % lapse_bh2 % lapse_bh % msr
00564 / nnn / r_bh / r_bh) * (2.+10.*msr+9.*msr%msr) * lshift
00565 + (4.*pow(lapse_bh2, 3.) % msr % msr / 3. / r_bh / r_bh)
00566 % (2.*(a_car%lapse_bh2/nnn/nnn - 1.) * pow(2.+3.*msr, 2.)
00567 + (a_car - 1.) % pow(1.+3.*msr, 2.)
00568 - 3.*(a_car%lapse_bh/nnn - 1.)*(2.+10.*msr+9.*msr%msr)) ;
00569
00570 }
00571 else {
00572 cout << "The computation of BH-NS binary systems"
00573 << " should be relativistic !!!" << endl ;
00574 abort() ;
00575 }
00576
00577 source.set_std_base() ;
00578
00579
00580
00581
00582 double* intvec_logn = new double[nylm] ;
00583
00584
00585 get_integrals(nylm, intvec_logn, ylmvec, source()) ;
00586
00587
00588 if(nu_int[0] != -1.0) {
00589 for (int i=0; i<nylm; i++) {
00590 intvec_logn[i] *= relax_ylm ;
00591 intvec_logn[i] += (1.0 - relax_ylm) * nu_int[i] ;
00592 }
00593 }
00594
00595 source().poisson_ylm(par_poisson1, logn_auto.set(), nylm,
00596 intvec_logn) ;
00597
00598 for (int i=0; i<nylm; i++) {
00599 nu_int[i] = intvec_logn[i] ;
00600 }
00601
00602 delete [] intvec_logn ;
00603
00604
00605
00606
00607 logn_auto_regu = logn_auto ;
00608
00609
00610
00611
00612 Tbl tdiff_logn = diffrel(logn_auto().laplacien(), source()) ;
00613
00614 cout <<
00615 "Relative error in the resolution of the equation for logn_auto : "
00616 << endl ;
00617
00618 for (int l=0; l<nz; l++) {
00619 cout << tdiff_logn(l) << " " ;
00620 }
00621 cout << endl ;
00622 diff_logn = max(abs(tdiff_logn)) ;
00623
00624 if (relativistic) {
00625
00626
00627
00628
00629
00630
00631
00632
00633 source = qpig * a_car % s_euler + 0.75 * akcar_auto
00634 - 0.5 * flat_scalar_prod_desal(d_logn_auto, d_logn_auto)
00635 - 0.5 * flat_scalar_prod_desal(d_beta_auto, d_beta_auto)
00636 + lapse_bh2 % msr % (ldnu%ldnu + ldbeta%ldbeta + 2.*ldldbeta)
00637 + lapse_bh2 % lapse_bh2 % msr % (2.*(1.+4.*msr) * ldbeta
00638 - ldnu) / r_bh
00639 - (a_car % lapse_bh2 % lapse_bh2 % msr / nnn / r_bh)
00640 * (2.+3.*msr) * (3.+4.*msr) * lltkij
00641 + (2.*a_car % lapse_bh2 % lapse_bh2 % lapse_bh % msr
00642 / nnn / r_bh / r_bh) * (2.+10.*msr+9.*msr%msr) * lshift
00643 + (2.*pow(lapse_bh2, 3.) % msr % msr / r_bh / r_bh)
00644 % ((a_car%lapse_bh2/nnn/nnn - 1.) * pow(2.+3.*msr, 2.)
00645 + (a_car - 1.) * pow(1.+3.*msr, 2.)
00646 - 2.*(a_car%lapse_bh/nnn - 1.)*(2.+10.*msr+9.*msr%msr)) ;
00647
00648 source.set_std_base() ;
00649
00650
00651
00652
00653 double* intvec_beta = new double[nylm] ;
00654
00655
00656 get_integrals(nylm, intvec_beta, ylmvec, source()) ;
00657
00658 if(beta_int[0] != -1.0) {
00659 for (int i=0; i<nylm; i++) {
00660 intvec_beta[i] *= relax_ylm ;
00661 intvec_beta[i] += (1.0 - relax_ylm) * beta_int[i] ;
00662 }
00663 }
00664
00665 source().poisson_ylm(par_poisson2, beta_auto.set(),
00666 nylm, intvec_beta) ;
00667
00668 for (int i=0; i<nylm; i++) {
00669 beta_int[i] = intvec_beta[i] ;
00670 }
00671
00672 delete [] intvec_beta ;
00673
00674
00675
00676
00677 Tbl tdiff_beta = diffrel(beta_auto().laplacien(), source()) ;
00678
00679 cout << "Relative error in the resolution of the equation for "
00680 << "beta_auto : " << endl ;
00681 for (int l=0; l<nz; l++) {
00682 cout << tdiff_beta(l) << " " ;
00683 }
00684 cout << endl ;
00685 diff_beta = max(abs(tdiff_beta)) ;
00686
00687
00688
00689
00690
00691
00692
00693
00694 Tenseur xx_con(mp, 1, CON, ref_triad) ;
00695 xx_con.set_etat_qcq() ;
00696 xx_con.set(0) = xx + sepa ;
00697 xx_con.set(1) = yy ;
00698 xx_con.set(2) = zz ;
00699 xx_con.set_std_base() ;
00700
00701 Tenseur xsr_con(mp, 1, CON, ref_triad) ;
00702 xsr_con = xx_con / r_bh ;
00703 xsr_con.set_std_base() ;
00704
00705
00706
00707
00708 Tenseur shift_auto_local = shift_auto ;
00709 shift_auto_local.change_triad( mp.get_bvect_cart() ) ;
00710
00711
00712
00713
00714
00715 Tenseur dn = shift_auto_local.gradient() ;
00716
00717
00718 dn.change_triad(ref_triad) ;
00719
00720
00721 Tenseur divn = contract(dn, 0, 1) ;
00722
00723
00724 Tenseur ldn_con = contract(xsr_con, 0, dn, 0) ;
00725
00726
00727 Tenseur ldn_local = ldn_con ;
00728 ldn_local.change_triad( mp.get_bvect_cart() ) ;
00729 Tenseur dldn = ldn_local.gradient() ;
00730 dldn.change_triad(ref_triad) ;
00731
00732
00733 Tenseur ldldn = contract(xsr_con, 0, dldn, 0) ;
00734
00735
00736 Tenseur ldn_cov = contract(xsr_cov, 0, dn, 1) ;
00737
00738
00739 Tenseur lldn_cov = contract(xsr_con, 0, ldn_cov, 0) ;
00740
00741
00742 Tenseur eldn(mp, 1, CON, ref_triad) ;
00743 eldn.set_etat_qcq() ;
00744 eldn.set(0) = ldn_cov(0) ;
00745 eldn.set(1) = ldn_cov(1) ;
00746 eldn.set(2) = ldn_cov(2) ;
00747 eldn.set_std_base() ;
00748
00749
00750 Tenseur ldivn = xsr_con % divn ;
00751
00752
00753 Tenseur ldivn_local = ldivn ;
00754 ldivn_local.change_triad( mp.get_bvect_cart() ) ;
00755 Tenseur dldivn = ldivn_local.gradient() ;
00756 dldivn.change_triad(ref_triad) ;
00757
00758
00759 Tenseur ldldivn = contract(xsr_con, 0, dldivn, 0) ;
00760
00761
00762 Tenseur ln = contract(xsr_cov, 0, shift_auto, 0) ;
00763
00764 Tenseur vtmp = 6. * d_beta_auto - 8. * d_logn_auto ;
00765
00766 Tenseur lvtmp = contract(xsr_con, 0, vtmp, 0) ;
00767
00768
00769 Tenseur evtmp(mp, 1, CON, ref_triad) ;
00770 evtmp.set_etat_qcq() ;
00771 evtmp.set(0) = vtmp(0) ;
00772 evtmp.set(1) = vtmp(1) ;
00773 evtmp.set(2) = vtmp(2) ;
00774 evtmp.set_std_base() ;
00775
00776
00777 Tenseur lapse_ns(mp) ;
00778 lapse_ns = exp(logn_auto) ;
00779 lapse_ns.set_std_base() ;
00780
00781
00782
00783
00784 source_shift = (-4.*qpig) * nnn % a_car % (ener_euler + press)
00785 % u_euler
00786 + nnn % flat_scalar_prod_desal(tkij_auto, vtmp)
00787 - 2.*nnn % lapse_bh2 * msr / r_bh
00788 % flat_scalar_prod_desal(tkij_auto, xsr_cov)
00789 + 2.*lapse_bh2 * msr * (3.*ldldn + ldldivn) / 3.
00790 - lapse_bh2 * msr / r_bh
00791 * (4.*ldivn - lapse_bh2 % (3.*ldn_con + 8.*msr * ldn_con)
00792 - (eldn + 2.*lapse_bh2*(9.+11.*msr)*lldn_cov%xsr_con) / 3.)
00793 - 2.*lapse_bh2 % lapse_bh2 * msr / r_bh / r_bh
00794 * ( (4.+11.*msr) * shift_auto
00795 - lapse_bh2 * (12.+51.*msr+46.*msr*msr) * ln % xsr_con )
00796 / 3.
00797 + 8.*pow(lapse_bh2, 4.) * msr / r_bh / r_bh
00798 % (lapse_ns - 1.) * (2.+10.*msr+9.*msr*msr) * xsr_con / 3.
00799 + 2.*pow(lapse_bh2, 3.) * msr / r_bh * (2.+3.*msr)
00800 * ( (1.+2.*msr) * evtmp - (3.+2.*msr) * lvtmp * xsr_con) / 3. ;
00801
00802 source_shift.set_std_base() ;
00803
00804
00805
00806
00807
00808
00809 for (int i=0; i<3; i++) {
00810 for (int l=0; l<nz; l++) {
00811 if (source_shift(i).get_etat() != ETATZERO)
00812 source_shift.set(i).filtre_phi(np_filter, l) ;
00813 }
00814 }
00815
00816
00817
00818
00819 source_shift.change_triad( mp.get_bvect_cart() ) ;
00820
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832 double lambda_shift = double(1) / double(3) ;
00833
00834 double* intvec_shift = new double[4*nylm] ;
00835 for (int i=0; i<3; i++) {
00836 double* intvec = new double[nylm];
00837 get_integrals(nylm, intvec, ylmvec, source_shift(i));
00838
00839 for (int j=0; j<nylm; j++) {
00840 intvec_shift[i*nylm+j] = intvec[j] ;
00841 }
00842 if(shift_int[i*nylm] != -1.0) {
00843 for (int j=0; j<nylm; j++) {
00844 intvec_shift[i*nylm+j] *= relax_ylm ;
00845 intvec_shift[i*nylm+j] += (1.0 - relax_ylm)
00846 * shift_int[i*nylm+j] ;
00847 }
00848 }
00849 delete [] intvec ;
00850 }
00851 Tenseur source_scal (-skxk(source_shift)) ;
00852 Cmp soscal (source_scal()) ;
00853 double* intvec2 = new double[nylm] ;
00854 get_integrals(nylm, intvec2, ylmvec, soscal) ;
00855
00856 for (int j=0; j<nylm; j++) {
00857 intvec_shift[3*nylm+j] = intvec2[j] ;
00858 }
00859 if(shift_int[3*nylm] != -1.0) {
00860 for (int i=0; i<nylm; i++) {
00861 intvec_shift[3*nylm+i] *= relax_ylm ;
00862 intvec_shift[3*nylm+i] += (1.0 - relax_ylm)
00863 *shift_int[3*nylm+i] ;
00864 }
00865 }
00866
00867 delete [] intvec2 ;
00868
00869 source_shift.poisson_vect_ylm(lambda_shift, par_poisson_vect,
00870 shift_auto, w_shift,
00871 khi_shift, nylm, intvec_shift) ;
00872
00873 for (int i=0; i<4*nylm; i++) {
00874 shift_int[i] = intvec_shift[i] ;
00875 }
00876
00877 delete[] intvec_shift ;
00878
00879
00880
00881
00882
00883 Tenseur divna = contract(shift_auto.gradient(), 0, 1) ;
00884
00885
00886
00887 Tenseur graddivn = divna.gradient() ;
00888
00889
00890
00891 Tenseur lap_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
00892 lap_shift.set_etat_qcq() ;
00893 for (int i=0; i<3; i++) {
00894 lap_shift.set(i) = shift_auto(i).laplacien()
00895 + lambda_shift * graddivn(i) ;
00896 }
00897
00898 Tbl tdiff_shift_x = diffrel(lap_shift(0), source_shift(0)) ;
00899 Tbl tdiff_shift_y = diffrel(lap_shift(1), source_shift(1)) ;
00900 Tbl tdiff_shift_z = diffrel(lap_shift(2), source_shift(2)) ;
00901
00902 cout << "Relative error in the resolution of the equation for "
00903 << "shift_auto : " << endl ;
00904 cout << "x component : " ;
00905 for (int l=0; l<nz; l++) {
00906 cout << tdiff_shift_x(l) << " " ;
00907 }
00908 cout << endl ;
00909 cout << "y component : " ;
00910 for (int l=0; l<nz; l++) {
00911 cout << tdiff_shift_y(l) << " " ;
00912 }
00913 cout << endl ;
00914 cout << "z component : " ;
00915 for (int l=0; l<nz; l++) {
00916 cout << tdiff_shift_z(l) << " " ;
00917 }
00918 cout << endl ;
00919
00920 diff_shift_x = max(abs(tdiff_shift_x)) ;
00921 diff_shift_y = max(abs(tdiff_shift_y)) ;
00922 diff_shift_z = max(abs(tdiff_shift_z)) ;
00923
00924
00925
00926
00927
00928
00929
00930 shift_auto.change_triad( ref_triad ) ;
00931
00932 }
00933
00934
00935
00936
00937
00938 Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
00939 diff_ent = diff_ent_tbl(0) ;
00940 for (int l=1; l<nzet; l++) {
00941 diff_ent += diff_ent_tbl(l) ;
00942 }
00943 diff_ent /= nzet ;
00944
00945 ent_jm1 = ent ;
00946
00947 for (int i=0; i<nylm; i++) {
00948 delete ylmvec[i];
00949 }
00950
00951 }
00952
00953
00954
00955
00956
00957 delete[] ylmvec ;
00958
00959 }
00960
00961
00962 void Et_bin_bhns_extr::equil_bhns_extr_ylm_cf(double ent_c,
00963 const double& mass,
00964 const double& sepa,
00965 double* nu_int,
00966 double* beta_int,
00967 double* shift_int,
00968 int mermax,
00969 int mermax_poisson,
00970 double relax_poisson,
00971 double relax_ylm,
00972 int mermax_potvit,
00973 double relax_potvit,
00974 int np_filter,
00975 double thres_adapt,
00976 Tbl& diff) {
00977
00978
00979
00980 using namespace Unites ;
00981
00982 assert( kerrschild == false ) ;
00983
00984
00985
00986
00987 const Mg3d* mg = mp.get_mg() ;
00988 int nz = mg->get_nzone() ;
00989
00990
00991 Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
00992
00993
00994 int l_b = nzet - 1 ;
00995 int i_b = mg->get_nr(l_b) - 1 ;
00996 int k_b ;
00997 int j_b ;
00998
00999
01000 double ent_b = 0 ;
01001
01002
01003
01004
01005 double& diff_ent = diff.set(0) ;
01006 double& diff_vel_pot = diff.set(1) ;
01007 double& diff_logn = diff.set(2) ;
01008 double& diff_beta = diff.set(3) ;
01009 double& diff_shift_x = diff.set(4) ;
01010 double& diff_shift_y = diff.set(5) ;
01011 double& diff_shift_z = diff.set(6) ;
01012
01013
01014
01015
01016 Param par_adapt ;
01017 int nitermax = 100 ;
01018 int niter ;
01019 int adapt_flag = 1 ;
01020
01021
01022 int nz_search = nzet ;
01023
01024 double precis_secant = 1.e-14 ;
01025 double alpha_r ;
01026 double reg_map = 1. ;
01027
01028 Tbl ent_limit(nz) ;
01029
01030 par_adapt.add_int(nitermax, 0) ;
01031
01032 par_adapt.add_int(nzet, 1) ;
01033
01034
01035 par_adapt.add_int(nz_search, 2) ;
01036
01037 par_adapt.add_int(adapt_flag, 3) ;
01038
01039
01040 par_adapt.add_int(j_b, 4) ;
01041
01042 par_adapt.add_int(k_b, 5) ;
01043
01044 par_adapt.add_int_mod(niter, 0) ;
01045
01046 par_adapt.add_double(precis_secant, 0) ;
01047
01048
01049 par_adapt.add_double(reg_map, 1) ;
01050
01051 par_adapt.add_double(alpha_r, 2) ;
01052
01053 par_adapt.add_tbl(ent_limit, 0) ;
01054
01055
01056
01057
01058
01059 double precis_poisson = 1.e-16 ;
01060
01061 Param par_poisson1 ;
01062
01063 par_poisson1.add_int(mermax_poisson, 0) ;
01064 par_poisson1.add_double(relax_poisson, 0) ;
01065 par_poisson1.add_double(precis_poisson, 1) ;
01066 par_poisson1.add_int_mod(niter, 0) ;
01067
01068 par_poisson1.add_cmp_mod( ssjm1_logn ) ;
01069
01070
01071
01072
01073 Param par_poisson2 ;
01074
01075 par_poisson2.add_int(mermax_poisson, 0) ;
01076 par_poisson2.add_double(relax_poisson, 0) ;
01077 par_poisson2.add_double(precis_poisson, 1) ;
01078 par_poisson2.add_int_mod(niter, 0) ;
01079
01080 par_poisson2.add_cmp_mod( ssjm1_beta ) ;
01081
01082
01083
01084
01085 Param par_poisson_vect ;
01086
01087 par_poisson_vect.add_int(mermax_poisson, 0) ;
01088
01089 par_poisson_vect.add_double(relax_poisson, 0) ;
01090 par_poisson_vect.add_double(precis_poisson, 1) ;
01091 par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
01092 par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
01093 par_poisson_vect.add_int_mod(niter, 0) ;
01094
01095
01096
01097
01098
01099 Tenseur pot_ext = logn_comp + pot_centri + loggam ;
01100
01101 Tenseur ent_jm1 = ent ;
01102
01103 Tenseur source(mp) ;
01104
01105
01106 Tenseur source_shift(mp, 1, CON, ref_triad) ;
01107
01108
01109
01110 int l_max = 4 ;
01111
01112 if (l_max > 6) {
01113 cout << "We don't have those ylm programmed yet!!!!" << endl ;
01114 abort() ;
01115 }
01116
01117 int nylm = (l_max+1) * (l_max+1) ;
01118
01119 Cmp** ylmvec = new Cmp* [nylm] ;
01120
01121
01122
01123
01124
01125 for(int mer=0 ; mer<mermax ; mer++ ) {
01126
01127 cout << "-----------------------------------------------" << endl ;
01128 cout << "step: " << mer << endl ;
01129 cout << "diff_ent = " << diff_ent << endl ;
01130
01131
01132
01133
01134
01135
01136 if (irrotational) {
01137 diff_vel_pot = velocity_pot_extr(mass, sepa, mermax_potvit,
01138 precis_poisson, relax_potvit) ;
01139 }
01140
01141
01142
01143
01144
01145
01146
01147
01148
01149 double logn_auto_c = logn_auto()(0, 0, 0, 0) ;
01150 double pot_ext_c = pot_ext()(0, 0, 0, 0) ;
01151
01152
01153
01154
01155
01156 double alpha_r2 = 0 ;
01157 for (int k=0; k<mg->get_np(l_b); k++) {
01158 for (int j=0; j<mg->get_nt(l_b); j++) {
01159
01160 double pot_ext_b = pot_ext()(l_b, k, j, i_b) ;
01161 double logn_auto_b = logn_auto()(l_b, k, j, i_b) ;
01162
01163
01164 double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b)
01165 / ( logn_auto_b - logn_auto_c ) ;
01166
01167 if (alpha_r2_jk > alpha_r2) {
01168 alpha_r2 = alpha_r2_jk ;
01169 k_b = k ;
01170 j_b = j ;
01171 }
01172
01173 }
01174 }
01175
01176 alpha_r = sqrt(alpha_r2) ;
01177
01178 cout << "k_b, j_b, alpha_r: " << k_b << " " << j_b << " "
01179 << alpha_r << endl ;
01180
01181
01182
01183
01184 logn_auto = alpha_r2 * logn_auto ;
01185 logn_auto_regu = alpha_r2 * logn_auto_regu ;
01186 logn_auto_c = logn_auto()(0, 0, 0, 0) ;
01187
01188
01189
01190
01191
01192
01193 (logn_auto().va).smooth(nzet, (logn_auto.set()).va) ;
01194
01195
01196
01197
01198
01199
01200 ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
01201
01202
01203
01204
01205
01206 double dentdx = ent().dsdx()(0, 0, 0, 0) ;
01207 double dentdy = ent().dsdy()(0, 0, 0, 0) ;
01208
01209 cout << "dH/dx|_center = " << dentdx << endl ;
01210 cout << "dH/dy|_center = " << dentdy << endl ;
01211
01212 double dec_fact = 1. ;
01213
01214 Tenseur func_in(mp) ;
01215 func_in.set_etat_qcq() ;
01216 func_in.set() = 1. - dec_fact * (dentdx/ent_c) * mp.x ;
01217 func_in.set().annule(nzet, nz-1) ;
01218 func_in.set_std_base() ;
01219
01220 Tenseur func_ex(mp) ;
01221 func_ex.set_etat_qcq() ;
01222 func_ex.set() = 1. ;
01223 func_ex.set().annule(0, nzet-1) ;
01224 func_ex.set_std_base() ;
01225
01226
01227
01228 ent.set() = ent() * (func_in() + func_ex()) ;
01229
01230 (ent().va).smooth(nzet, (ent.set()).va) ;
01231
01232 double dentdx_new = ent().dsdx()(0, 0, 0, 0) ;
01233
01234 cout << "dH/dx|_new = " << dentdx_new << endl ;
01235
01236
01237
01238
01239
01240
01241
01242
01243 double dent_eq = ent().dsdr().val_point(ray_eq_pi(),M_PI/2.,M_PI) ;
01244 double dent_pole = ent().dsdr().val_point(ray_pole(),0.,0.) ;
01245 double rap_dent = fabs( dent_eq / dent_pole ) ;
01246 cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
01247
01248 if ( rap_dent < thres_adapt ) {
01249 adapt_flag = 0 ;
01250 cout << "******* FROZEN MAPPING *********" << endl ;
01251 }
01252 else{
01253 adapt_flag = 1 ;
01254
01255 }
01256
01257 ent_limit.set_etat_qcq() ;
01258 for (int l=0; l<nzet; l++) {
01259 ent_limit.set(l) = ent()(l, k_b, j_b, i_b) ;
01260 }
01261
01262 ent_limit.set(nzet-1) = ent_b ;
01263 Map_et mp_prev = mp_et ;
01264
01265 mp.adapt(ent(), par_adapt) ;
01266
01267
01268
01269
01270
01271 mp_prev.homothetie(alpha_r) ;
01272
01273 mp.reevaluate_symy(&mp_prev, nzet+1, ent.set()) ;
01274
01275 double fact_resize = 1. / alpha_r ;
01276 for (int l=nzet; l<nz-1; l++) {
01277 mp_et.resize(l, fact_resize) ;
01278 }
01279 mp_et.resize_extr(fact_resize) ;
01280
01281 double ent_s_max = -1 ;
01282 int k_s_max = -1 ;
01283 int j_s_max = -1 ;
01284 for (int k=0; k<mg->get_np(l_b); k++) {
01285 for (int j=0; j<mg->get_nt(l_b); j++) {
01286 double xx = fabs( ent()(l_b, k, j, i_b) ) ;
01287 if (xx > ent_s_max) {
01288 ent_s_max = xx ;
01289 k_s_max = k ;
01290 j_s_max = j ;
01291 }
01292 }
01293 }
01294 cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
01295 << " and nzet : " << endl ;
01296 cout << " " << ent_s_max << " reached for k = " << k_s_max
01297 << " and j = " << j_s_max << endl ;
01298
01299
01300
01301
01302
01303
01304
01305 int oldindex = -1 ;
01306 for (int l=0; l<=l_max; l++) {
01307 for (int m=0; m<= l; m++) {
01308
01309 int index = l*l + 2*m ;
01310 if(m > 0) index -= 1 ;
01311 if(index != oldindex+1) {
01312 cout << "error!! " << l << " " << m << " " << index
01313 << " " << oldindex << endl ;
01314 abort() ;
01315 }
01316
01317
01318 if ((l+m) % 2 == 0) {
01319 ylmvec[index] = new Cmp (logn_auto()) ;
01320 } else {
01321 ylmvec[index] = new Cmp (shift_auto(2)) ;
01322 }
01323 oldindex = index ;
01324 if (m > 0) {
01325 index += 1 ;
01326 if((l+m) % 2 == 0) {
01327 ylmvec[index] = new Cmp (logn_auto()) ;
01328 } else {
01329 ylmvec[index] = new Cmp (shift_auto(2)) ;
01330 }
01331 oldindex = index ;
01332 }
01333 }
01334 }
01335 if(oldindex+1 != nylm) {
01336 cout << "ERROR!!! " << oldindex << " "<< nylm << endl ;
01337 abort() ;
01338 }
01339
01340
01341
01342
01343 get_ylm(nylm, ylmvec);
01344
01345
01346
01347
01348
01349 equation_of_state() ;
01350
01351
01352
01353
01354
01355
01356 hydro_euler_extr(mass, sepa) ;
01357
01358
01359
01360
01361
01362
01363
01364 const Coord& xx = mp.x ;
01365 const Coord& yy = mp.y ;
01366 const Coord& zz = mp.z ;
01367
01368 Tenseur r_bh(mp) ;
01369 r_bh.set_etat_qcq() ;
01370 r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
01371 r_bh.set_std_base() ;
01372
01373 Tenseur xx_cov(mp, 1, COV, ref_triad) ;
01374 xx_cov.set_etat_qcq() ;
01375 xx_cov.set(0) = xx + sepa ;
01376 xx_cov.set(1) = yy ;
01377 xx_cov.set(2) = zz ;
01378 xx_cov.set_std_base() ;
01379
01380 Tenseur msr(mp) ;
01381 msr = ggrav * mass / r_bh ;
01382 msr.set_std_base() ;
01383
01384 Tenseur tmp(mp) ;
01385 tmp = 1. / ( 1. - 0.25*msr*msr ) ;
01386 tmp.set_std_base() ;
01387
01388 Tenseur xdnu(mp) ;
01389 xdnu = flat_scalar_prod_desal(xx_cov, d_logn_auto) ;
01390 xdnu.set_std_base() ;
01391
01392 Tenseur xdbeta(mp) ;
01393 xdbeta = flat_scalar_prod_desal(xx_cov, d_beta_auto) ;
01394 xdbeta.set_std_base() ;
01395
01396
01397
01398
01399
01400
01401
01402
01403 if (relativistic) {
01404
01405 source = qpig * a_car % (ener_euler + s_euler) + akcar_auto
01406 - flat_scalar_prod_desal(d_logn_auto, d_beta_auto)
01407 - 0.5 * tmp % msr % msr % xdnu / r_bh / r_bh
01408 - tmp % msr % xdbeta / r_bh / r_bh ;
01409
01410 }
01411 else {
01412 cout << "The computation of BH-NS binary systems"
01413 << " should be relativistic !!!" << endl ;
01414 abort() ;
01415 }
01416
01417 source.set_std_base() ;
01418
01419
01420
01421
01422 double* intvec_logn = new double[nylm] ;
01423
01424
01425 get_integrals(nylm, intvec_logn, ylmvec, source()) ;
01426
01427
01428 if(nu_int[0] != -1.0) {
01429 for (int i=0; i<nylm; i++) {
01430 intvec_logn[i] *= relax_ylm ;
01431 intvec_logn[i] += (1.0 - relax_ylm) * nu_int[i] ;
01432 }
01433 }
01434
01435 source().poisson_ylm(par_poisson1, logn_auto.set(), nylm,
01436 intvec_logn) ;
01437
01438 for (int i=0; i<nylm; i++) {
01439 nu_int[i] = intvec_logn[i] ;
01440 }
01441
01442 delete [] intvec_logn ;
01443
01444
01445
01446
01447 logn_auto_regu = logn_auto ;
01448
01449
01450
01451
01452 Tbl tdiff_logn = diffrel(logn_auto().laplacien(), source()) ;
01453
01454 cout <<
01455 "Relative error in the resolution of the equation for logn_auto : "
01456 << endl ;
01457
01458 for (int l=0; l<nz; l++) {
01459 cout << tdiff_logn(l) << " " ;
01460 }
01461 cout << endl ;
01462 diff_logn = max(abs(tdiff_logn)) ;
01463
01464 if (relativistic) {
01465
01466
01467
01468
01469
01470
01471
01472
01473 source = qpig * a_car % s_euler + 0.75 * akcar_auto
01474 - 0.5 * flat_scalar_prod_desal(d_logn_auto, d_logn_auto)
01475 - 0.5 * flat_scalar_prod_desal(d_beta_auto, d_beta_auto)
01476 - tmp % msr % xdnu / r_bh / r_bh
01477 - 0.5 * tmp % msr %msr % xdbeta / r_bh / r_bh ;
01478
01479 source.set_std_base() ;
01480
01481
01482
01483
01484 double* intvec_beta = new double[nylm] ;
01485
01486
01487 get_integrals(nylm, intvec_beta, ylmvec, source()) ;
01488
01489 if(beta_int[0] != -1.0) {
01490 for (int i=0; i<nylm; i++) {
01491 intvec_beta[i] *= relax_ylm ;
01492 intvec_beta[i] += (1.0 - relax_ylm) * beta_int[i] ;
01493 }
01494 }
01495
01496 source().poisson_ylm(par_poisson2, beta_auto.set(),
01497 nylm, intvec_beta) ;
01498
01499 for (int i=0; i<nylm; i++) {
01500 beta_int[i] = intvec_beta[i] ;
01501 }
01502
01503 delete [] intvec_beta ;
01504
01505
01506
01507
01508 Tbl tdiff_beta = diffrel(beta_auto().laplacien(), source()) ;
01509
01510 cout << "Relative error in the resolution of the equation for "
01511 << "beta_auto : " << endl ;
01512 for (int l=0; l<nz; l++) {
01513 cout << tdiff_beta(l) << " " ;
01514 }
01515 cout << endl ;
01516 diff_beta = max(abs(tdiff_beta)) ;
01517
01518
01519
01520
01521
01522
01523
01524
01525 Tenseur bhtmp(mp, 1, COV, ref_triad) ;
01526 bhtmp.set_etat_qcq() ;
01527 bhtmp = tmp % msr % (3.*msr-8.) % xx_cov / r_bh / r_bh ;
01528 bhtmp.set_std_base() ;
01529
01530 Tenseur vtmp = 6. * d_beta_auto - 8. * d_logn_auto ;
01531
01532
01533
01534
01535 source_shift = (-4.*qpig) * nnn % a_car % (ener_euler + press)
01536 % u_euler
01537 + nnn % flat_scalar_prod_desal(tkij_auto, vtmp+bhtmp) ;
01538
01539 source_shift.set_std_base() ;
01540
01541
01542
01543
01544
01545
01546 for (int i=0; i<3; i++) {
01547 for (int l=0; l<nz; l++) {
01548 if (source_shift(i).get_etat() != ETATZERO)
01549 source_shift.set(i).filtre_phi(np_filter, l) ;
01550 }
01551 }
01552
01553
01554
01555
01556 source_shift.change_triad( mp.get_bvect_cart() ) ;
01557
01558
01559
01560
01561
01562
01563
01564
01565
01566
01567
01568
01569 double lambda_shift = double(1) / double(3) ;
01570
01571 double* intvec_shift = new double[4*nylm] ;
01572 for (int i=0; i<3; i++) {
01573 double* intvec = new double[nylm];
01574 get_integrals(nylm, intvec, ylmvec, source_shift(i));
01575
01576 for (int j=0; j<nylm; j++) {
01577 intvec_shift[i*nylm+j] = intvec[j] ;
01578 }
01579 if(shift_int[i*nylm] != -1.0) {
01580 for (int j=0; j<nylm; j++) {
01581 intvec_shift[i*nylm+j] *= relax_ylm ;
01582 intvec_shift[i*nylm+j] += (1.0 - relax_ylm)
01583 * shift_int[i*nylm+j] ;
01584 }
01585 }
01586 delete [] intvec ;
01587 }
01588 Tenseur source_scal (-skxk(source_shift)) ;
01589 Cmp soscal (source_scal()) ;
01590 double* intvec2 = new double[nylm] ;
01591 get_integrals(nylm, intvec2, ylmvec, soscal) ;
01592
01593 for (int j=0; j<nylm; j++) {
01594 intvec_shift[3*nylm+j] = intvec2[j] ;
01595 }
01596 if(shift_int[3*nylm] != -1.0) {
01597 for (int i=0; i<nylm; i++) {
01598 intvec_shift[3*nylm+i] *= relax_ylm ;
01599 intvec_shift[3*nylm+i] += (1.0 - relax_ylm)
01600 *shift_int[3*nylm+i] ;
01601 }
01602 }
01603
01604 delete [] intvec2 ;
01605
01606 source_shift.poisson_vect_ylm(lambda_shift, par_poisson_vect,
01607 shift_auto, w_shift,
01608 khi_shift, nylm, intvec_shift) ;
01609
01610 for (int i=0; i<4*nylm; i++) {
01611 shift_int[i] = intvec_shift[i] ;
01612 }
01613
01614 delete[] intvec_shift ;
01615
01616
01617
01618
01619
01620 Tenseur divna = contract(shift_auto.gradient(), 0, 1) ;
01621
01622
01623
01624 Tenseur graddivn = divna.gradient() ;
01625
01626
01627
01628 Tenseur lap_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
01629 lap_shift.set_etat_qcq() ;
01630 for (int i=0; i<3; i++) {
01631 lap_shift.set(i) = shift_auto(i).laplacien()
01632 + lambda_shift * graddivn(i) ;
01633 }
01634
01635 Tbl tdiff_shift_x = diffrel(lap_shift(0), source_shift(0)) ;
01636 Tbl tdiff_shift_y = diffrel(lap_shift(1), source_shift(1)) ;
01637 Tbl tdiff_shift_z = diffrel(lap_shift(2), source_shift(2)) ;
01638
01639 cout << "Relative error in the resolution of the equation for "
01640 << "shift_auto : " << endl ;
01641 cout << "x component : " ;
01642 for (int l=0; l<nz; l++) {
01643 cout << tdiff_shift_x(l) << " " ;
01644 }
01645 cout << endl ;
01646 cout << "y component : " ;
01647 for (int l=0; l<nz; l++) {
01648 cout << tdiff_shift_y(l) << " " ;
01649 }
01650 cout << endl ;
01651 cout << "z component : " ;
01652 for (int l=0; l<nz; l++) {
01653 cout << tdiff_shift_z(l) << " " ;
01654 }
01655 cout << endl ;
01656
01657 diff_shift_x = max(abs(tdiff_shift_x)) ;
01658 diff_shift_y = max(abs(tdiff_shift_y)) ;
01659 diff_shift_z = max(abs(tdiff_shift_z)) ;
01660
01661
01662
01663
01664
01665
01666
01667 shift_auto.change_triad( ref_triad ) ;
01668
01669 }
01670
01671
01672
01673
01674
01675 Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
01676 diff_ent = diff_ent_tbl(0) ;
01677 for (int l=1; l<nzet; l++) {
01678 diff_ent += diff_ent_tbl(l) ;
01679 }
01680 diff_ent /= nzet ;
01681
01682 ent_jm1 = ent ;
01683
01684 for (int i=0; i<nylm; i++) {
01685 delete ylmvec[i];
01686 }
01687
01688 }
01689
01690
01691
01692
01693
01694 delete[] ylmvec ;
01695
01696 }
01697