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 char binary_global_C[] = "$Header: /cvsroot/Lorene/C++/Source/Binary/binary_global.C,v 1.15 2006/08/01 14:26:50 f_limousin Exp $" ;
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 #include "math.h"
00088
00089
00090 #include "nbr_spx.h"
00091 #include "binary.h"
00092 #include "unites.h"
00093 #include "metric.h"
00094
00095
00096
00097
00098
00099 double Binary::mass_adm() const {
00100
00101 using namespace Unites ;
00102 if (p_mass_adm == 0x0) {
00103
00104 p_mass_adm = new double ;
00105
00106 *p_mass_adm = 0 ;
00107
00108 const Map_af map0 (et[0]->get_mp()) ;
00109 const Metric& flat = (et[0]->get_flat()) ;
00110
00111 Vector dpsi(0.5*(et[0]->get_lnq() -
00112 et[0]->get_logn()).derive_cov(flat)) ;
00113
00114 Vector ww (0.125*(contract(et[0]->get_hij().derive_cov(flat), 1, 2)
00115 - (et[0]->get_hij().trace(flat)).derive_con(flat))) ;
00116
00117 dpsi.change_triad(map0.get_bvect_spher()) ;
00118 ww.change_triad(map0.get_bvect_spher()) ;
00119
00120
00121 Scalar integrand (dpsi(1) + 0*ww(1)) ;
00122
00123 *p_mass_adm = map0.integrale_surface_infini (integrand) / (-qpig/2.) ;
00124
00125 }
00126
00127 return *p_mass_adm ;
00128
00129 }
00130
00131 double Binary::mass_adm_vol() const {
00132
00133 using namespace Unites ;
00134
00135 double massadm ;
00136 massadm = 0. ;
00137
00138 for (int i=0; i<=1; i++) {
00139
00140
00141 const Scalar& psi4 = et[i]->get_psi4() ;
00142 Scalar psi (pow(psi4, 0.25)) ;
00143 psi.std_spectral_base() ;
00144 const Scalar& ener_euler = et[i]->get_ener_euler() ;
00145 const Scalar& kcar_auto = et[i]->get_kcar_auto() ;
00146 const Scalar& kcar_comp = et[i]->get_kcar_comp() ;
00147 const Metric& gtilde = et[i]->get_gtilde() ;
00148 const Metric& flat = et[i]->get_flat() ;
00149 const Sym_tensor& hij = et[i]->get_hij() ;
00150 const Sym_tensor& hij_auto = et[i]->get_hij_auto() ;
00151 const Vector& dcov_logn = et[i]->get_dcov_logn() ;
00152 const Vector& dcov_phi = et[i]->get_dcov_phi() ;
00153 const Vector& dcov_lnq = 2*dcov_phi + dcov_logn ;
00154 const Scalar& lnq_auto = et[i]->get_lnq_auto() ;
00155 const Scalar& logn_auto = et[i]->get_logn_auto() ;
00156 const Scalar& phi_auto = 0.5 * (lnq_auto - logn_auto) ;
00157
00158 const Tensor& dcov_hij_auto = hij_auto.derive_cov(flat) ;
00159 const Tensor& dcov_gtilde = gtilde.cov().derive_cov(flat) ;
00160 const Tensor& dcov_phi_auto = phi_auto.derive_cov(flat) ;
00161 const Tensor& dcov_logn_auto = logn_auto.derive_cov(flat) ;
00162 const Tensor& dcov_lnq_auto = lnq_auto.derive_cov(flat) ;
00163 Tensor dcovdcov_lnq_auto = lnq_auto.derive_cov(flat).derive_cov(flat) ;
00164 dcovdcov_lnq_auto.inc_dzpuis() ;
00165 Tensor dcovdcov_logn_auto = logn_auto.derive_cov(flat).derive_cov(flat) ;
00166 dcovdcov_logn_auto.inc_dzpuis() ;
00167
00168
00169 Scalar source = - psi4 % (qpig*ener_euler + (kcar_auto + kcar_comp)/4.)
00170 - 0*2*contract(contract(gtilde.con(), 0, dcov_phi, 0),
00171 0, dcov_phi_auto, 0, true) ;
00172
00173
00174 source += 4*contract(hij, 0, 1, dcov_logn * dcov_phi_auto, 0, 1) +
00175 2*contract(hij, 0, 1, dcov_phi * dcov_phi_auto, 0, 1) +
00176 0.0625 * contract(gtilde.con(), 0, 1, contract(
00177 dcov_hij_auto, 0, 1, dcov_gtilde, 0, 1), 0, 1) -
00178 0.125 * contract(gtilde.con(), 0, 1, contract(dcov_hij_auto,
00179 0, 1, dcov_gtilde, 0, 2), 0, 1) -
00180 contract(hij,0,1,dcovdcov_lnq_auto + dcov_lnq_auto*dcov_lnq,0,1) +
00181 contract(hij,0,1,dcovdcov_logn_auto + dcov_logn_auto*dcov_logn,0,1) ;
00182
00183 source = source * psi ;
00184
00185 source.std_spectral_base() ;
00186
00187 massadm += - source.integrale()/qpig ;
00188 }
00189
00190 return massadm ;
00191 }
00192
00193
00194
00195
00196
00197 double Binary::mass_kom() const {
00198
00199 using namespace Unites ;
00200
00201 if (p_mass_kom == 0x0) {
00202
00203 p_mass_kom = new double ;
00204
00205 *p_mass_kom = 0 ;
00206
00207 const Tensor& logn = et[0]->get_logn() ;
00208 const Metric& flat = (et[0]->get_flat()) ;
00209 const Sym_tensor& hij = (et[0]->get_hij()) ;
00210 Map_af map0 (et[0]->get_mp()) ;
00211
00212 Vector vect = logn.derive_con(flat) +
00213 contract(hij, 1, logn.derive_cov(flat), 0) ;
00214 vect.change_triad(map0.get_bvect_spher()) ;
00215 Scalar integrant (vect(1)) ;
00216
00217 *p_mass_kom = map0.integrale_surface_infini (integrant) / qpig ;
00218
00219 }
00220
00221 return *p_mass_kom ;
00222
00223 }
00224
00225 double Binary::mass_kom_vol() const {
00226
00227 using namespace Unites ;
00228
00229 double masskom ;
00230 masskom = 0. ;
00231
00232 for (int i=0; i<=1; i++) {
00233
00234
00235 const Scalar& psi4 = et[i]->get_psi4() ;
00236 const Scalar& ener_euler = et[i]->get_ener_euler() ;
00237 const Scalar& s_euler = et[i]->get_s_euler() ;
00238 const Scalar& kcar_auto = et[i]->get_kcar_auto() ;
00239 const Scalar& kcar_comp = et[i]->get_kcar_comp() ;
00240 const Metric& gtilde = et[i]->get_gtilde() ;
00241 const Metric& flat = et[i]->get_flat() ;
00242 const Sym_tensor& hij = et[i]->get_hij() ;
00243 const Scalar& logn = et[i]->get_logn_auto() + et[i]->get_logn_comp() ;
00244 const Scalar& logn_auto = et[i]->get_logn_auto() ;
00245 Scalar nn = exp(logn) ;
00246 nn.std_spectral_base() ;
00247
00248 const Tensor& dcov_logn_auto = logn_auto.derive_cov(flat) ;
00249 const Vector& dcov_logn = et[i]->get_dcov_logn() ;
00250 const Vector& dcon_logn = et[i]->get_dcon_logn() ;
00251 const Vector& dcov_phi = et[i]->get_dcov_phi() ;
00252 Tensor dcovdcov_logn_auto = (logn_auto.derive_cov(flat))
00253 .derive_cov(flat) ;
00254 dcovdcov_logn_auto.inc_dzpuis() ;
00255
00256 Scalar source = qpig * psi4 % (ener_euler + s_euler) ;
00257 source += psi4 % (kcar_auto + kcar_comp) ;
00258 source += - 0*contract(dcov_logn_auto, 0, dcon_logn, 0, true)
00259 - 2. * contract(contract(gtilde.con(), 0, dcov_phi, 0), 0,
00260 dcov_logn_auto, 0, true) ;
00261 source += - contract(hij, 0, 1, dcovdcov_logn_auto +
00262 dcov_logn_auto*dcov_logn, 0, 1) ;
00263
00264 source = source / qpig * nn ;
00265
00266 source.std_spectral_base() ;
00267
00268 masskom += source.integrale() ;
00269
00270 }
00271
00272 return masskom ;
00273
00274 }
00275
00276
00277
00278
00279
00280
00281 const Tbl& Binary::angu_mom() const {
00282
00283 using namespace Unites ;
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448
00449
00450
00451
00452
00453
00454
00455
00456
00457
00458
00459
00460
00461
00462
00463
00464
00465
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475
00476
00477
00478
00479
00480
00481
00482
00483
00484
00485
00486 if (p_angu_mom == 0x0) {
00487
00488 p_angu_mom = new Tbl(3) ;
00489 p_angu_mom->annule_hard() ;
00490
00491
00492 Base_vect_cart bvect_ref(0.) ;
00493
00494 for (int i=0; i<=1; i++) {
00495
00496 const Map& mp = et[i]->get_mp() ;
00497 int nzm1 = mp.get_mg()->get_nzone() - 1 ;
00498
00499
00500
00501
00502 double r0 = mp.val_r(nzm1-1, 1, 0, 0) ;
00503 double sigma = 1.*r0 ;
00504
00505 Scalar rr (mp) ;
00506 rr = mp.r ;
00507
00508 Scalar ff (mp) ;
00509 ff = exp( -(rr - r0)*(rr - r0)/sigma/sigma ) ;
00510 for (int ii=0; ii<nzm1; ii++)
00511 ff.set_domain(ii) = 1. ;
00512 ff.set_outer_boundary(nzm1, 0) ;
00513 ff.std_spectral_base() ;
00514
00515
00516 Vector vphi(mp, CON, bvect_ref) ;
00517 Scalar yya (mp) ;
00518 yya = mp.ya ;
00519 Scalar xxa (mp) ;
00520 xxa = mp.xa ;
00521 vphi.set(1) = - yya * ff ;
00522 vphi.set(2) = xxa * ff ;
00523 vphi.set(3) = 0 ;
00524
00525 vphi.set(1).set_outer_boundary(nzm1, 0) ;
00526 vphi.set(2).set_outer_boundary(nzm1, 0) ;
00527
00528 vphi.std_spectral_base() ;
00529 vphi.change_triad(mp.get_bvect_cart()) ;
00530
00531
00532
00533 const Scalar& ee = et[i]->get_ener_euler() ;
00534 const Scalar& pp = et[i]->get_press() ;
00535 const Scalar& psi4 = et[i]->get_psi4() ;
00536 Scalar rho = pow(psi4, double(2.5)) * (ee+pp) ;
00537 rho.std_spectral_base() ;
00538
00539 Vector jmom = rho * (et[i]->get_u_euler()) ;
00540
00541 const Metric& gtilde = et[i]->get_gtilde() ;
00542 const Metric_flat flat (mp.flat_met_cart()) ;
00543
00544 Vector vphi_cov = vphi.up_down(gtilde) ;
00545
00546 Scalar integrand = contract(jmom, 0, vphi_cov, 0) ;
00547
00548 p_angu_mom->set(2) += integrand.integrale() ;
00549
00550
00551
00552
00553 const Sym_tensor& aij = et[i]->get_tkij_auto() ;
00554 rho = pow(psi4, double(1.5)) ;
00555 rho.std_spectral_base() ;
00556
00557
00558 Itbl type (2) ;
00559 type.set(0) = CON ;
00560 type.set(1) = COV ;
00561
00562 Tensor dcov_vphi (mp, 2, type, mp.get_bvect_cart()) ;
00563 dcov_vphi.set(1,1) = 0. ;
00564 dcov_vphi.set(2,1) = ff ;
00565 dcov_vphi.set(3,1) = 0. ;
00566 dcov_vphi.set(2,2) = 0. ;
00567 dcov_vphi.set(3,2) = 0. ;
00568 dcov_vphi.set(3,3) = 0. ;
00569 dcov_vphi.set(1,2) = -ff ;
00570 dcov_vphi.set(1,3) = 0. ;
00571 dcov_vphi.set(2,3) = 0. ;
00572 dcov_vphi.inc_dzpuis(2) ;
00573
00574 Connection gamijk (gtilde, flat) ;
00575 const Tensor& deltaijk = gamijk.get_delta() ;
00576
00577
00578 Sym_tensor kill_phi (mp, COV, mp.get_bvect_cart()) ;
00579 kill_phi = contract(gtilde.cov(), 1, dcov_vphi +
00580 contract(deltaijk, 2, vphi, 0), 0) +
00581 contract(dcov_vphi + contract(deltaijk, 2, vphi, 0), 0,
00582 gtilde.cov(), 1) ;
00583
00584 integrand = rho * contract(aij, 0, 1, kill_phi, 0, 1) ;
00585
00586 p_angu_mom->set(2) += integrand.integrale() / (4*qpig) ;
00587
00588
00589 }
00590
00591 }
00592
00593 return *p_angu_mom ;
00594
00595 }
00596
00597
00598
00599
00600
00601
00602
00603 double Binary::total_ener() const {
00604
00605
00606
00607
00608
00609
00610
00611
00612
00613
00614 return *p_total_ener ;
00615
00616 }
00617
00618
00619
00620
00621
00622
00623 double Binary::virial() const {
00624
00625 if (p_virial == 0x0) {
00626
00627 p_virial = new double ;
00628
00629 *p_virial = 1. - mass_kom() / mass_adm() ;
00630
00631 }
00632
00633 return *p_virial ;
00634
00635 }