75#include "utilitaires.h"
81double Bin_ns_bh::adm_systeme()
const {
82 Cmp der_un (
hole.psi_auto().dsdr()) ;
84 Map_af auxi_mp (
star.get_mp()) ;
85 Cmp der_deux (
star.confpsi_auto().dsdr()) ;
87 double masse =
hole.mp.integrale_surface_infini(der_un) +
88 auxi_mp.integrale_surface_infini(der_deux) ;
94double Bin_ns_bh::adm_systeme_volume()
const {
96 using namespace Unites ;
99 Tenseur kk_bh (
hole.mp) ;
101 Tenseur work_bh(
hole.mp) ;
102 work_bh.set_etat_qcq() ;
103 for (
int i=0 ; i<3 ; i++) {
104 work_bh.set() = auxi_bh(i, i) ;
105 kk_bh = kk_bh + work_bh ;
107 Cmp integ_bh (
pow(
hole.psi_tot(), 5.)*kk_bh()) ;
109 integ_bh.std_base_scal() ;
111 Cmp integ_hor1 (
hole.psi_tot()) ;
117 Tenseur rad (
hole.mp, 1, COV,
hole.mp.get_bvect_cart()) ;
123 Cmp integ_hor2 (
hole.mp) ;
124 integ_hor2.annule_hard() ;
125 integ_hor2.set_dzpuis(2) ;
126 for (
int m=0 ; m<3 ; m++)
127 for (
int n=0 ; n<3 ; n++)
128 integ_hor2 += rad(m)*rad(n)*
hole.tkij_tot(m,n) ;
129 integ_hor2 *=
pow(
hole.psi_tot(),3.)/4. ;
133 Tenseur kk_ns (
star.mp) ;
135 Tenseur work_ns(
star.mp) ;
136 work_ns.set_etat_qcq() ;
137 for (
int i=0 ; i<3 ; i++) {
138 work_ns.set() = auxi_ns(i, i) ;
139 kk_ns = kk_ns + work_ns ;
141 Cmp integ_ns (
pow(
star.confpsi_comp() +
star.confpsi_auto(), 5.)*kk_ns()) ;
142 integ_ns.std_base_scal() ;
144 Cmp integ_matter (
pow(
star.confpsi_comp()+
star.confpsi_auto(), 5.)*
star.ener_euler()) ;
145 integ_matter.std_base_scal() ;
147 double masse = (integ_bh.integrale()+integ_ns.integrale())/16/M_PI +
148 hole.mp.integrale_surface(integ_hor1,
hole.rayon)/
hole.rayon/4/M_PI +
149 hole.mp.integrale_surface(integ_hor2,
hole.rayon)/2/M_PI
150 + integ_matter.integrale()*ggrav ;
155double Bin_ns_bh::komar_systeme()
const {
156 Cmp der_un (
hole.n_auto().dsdr()) ;
158 Map_af auxi_mp (
star.get_mp()) ;
159 Cmp der_deux (
star.n_auto().dsdr()) ;
161 double masse =
hole.mp.integrale_surface_infini(der_un) +
162 auxi_mp.integrale_surface_infini(der_deux) ;
169double Bin_ns_bh::viriel()
const {
170 double adm = adm_systeme() ;
171 double komar = komar_systeme() ;
173 return (adm-komar)/adm ;
176double Bin_ns_bh::moment_systeme_inf()
const {
183 int nzones =
hole.mp.get_mg()->get_nzone() ;
184 double* bornes =
new double [nzones+1] ;
185 double courant = fabs(
hole.mp.get_ori_x()-
star.mp.get_ori_x())+1 ;
186 for (
int i=nzones-1 ; i>0 ; i--) {
187 bornes[i] = courant ;
191 bornes[nzones] = __infinity ;
193 Map_af mapping (*
hole.mp.get_mg(), bornes) ;
198 Tenseur_sym k_total (mapping, 2, CON, mapping.get_bvect_cart()) ;
199 k_total.set_etat_qcq() ;
201 Tenseur shift_un (mapping, 1, CON, mapping.get_bvect_cart()) ;
202 shift_un.set_etat_qcq() ;
204 Tenseur shift_deux (mapping, 1, CON, mapping.get_bvect_cart()) ;
205 shift_deux.set_etat_qcq() ;
207 shift_un.set_triad (*
hole.shift_auto.get_triad()) ;
208 shift_un.set(0).import(
hole.shift_auto(0)) ;
209 shift_un.set(1).import(
hole.shift_auto(1)) ;
210 shift_un.set(2).import(
hole.shift_auto(2)) ;
211 shift_un.change_triad (mapping.get_bvect_cart()) ;
213 shift_deux.set_triad (*
star.shift_auto.get_triad()) ;
214 shift_deux.set(0).import(
star.shift_auto(0)) ;
215 shift_deux.set(1).import(
star.shift_auto(1)) ;
216 shift_deux.set(2).import(
star.shift_auto(2)) ;
217 shift_deux.change_triad(mapping.get_bvect_cart()) ;
219 Tenseur shift_tot (shift_un+shift_deux) ;
220 shift_tot.set_std_base() ;
221 shift_tot.annule(0, nzones-2) ;
223 shift_tot.inc2_dzpuis() ;
224 shift_tot.dec2_dzpuis() ;
226 Tenseur grad (shift_tot.gradient()) ;
227 Tenseur trace (grad(0, 0)+grad(1, 1)+grad(2, 2)) ;
228 for (
int i=0 ; i<3 ; i++) {
229 k_total.set(i, i) = grad(i, i)-trace()/3. ;
230 for (
int j=i+1 ; j<3 ; j++)
231 k_total.set(i, j) = (grad(i, j)+grad(j, i))/2. ;
234 for (
int lig=0 ; lig<3 ; lig++)
235 for (
int col=lig ; col<3 ; col++)
236 k_total.set(lig, col).mult_r_zec() ;
238 Tenseur vecteur_un (mapping, 1, CON, mapping.get_bvect_cart()) ;
239 vecteur_un.set_etat_qcq() ;
240 for (
int i=0 ; i<3 ; i++)
241 vecteur_un.set(i) = k_total(0, i) ;
242 vecteur_un.change_triad (mapping.get_bvect_spher()) ;
243 Cmp integrant_un (vecteur_un(0)) ;
245 Tenseur vecteur_deux (mapping, 1, CON, mapping.get_bvect_cart()) ;
246 vecteur_deux.set_etat_qcq() ;
247 for (
int i=0 ; i<3 ; i++)
248 vecteur_deux.set(i) = k_total(1, i) ;
249 vecteur_deux.change_triad (mapping.get_bvect_spher()) ;
250 Cmp integrant_deux (vecteur_deux(0)) ;
253 integrant_un.va = integrant_un.va.mult_st() ;
254 integrant_un.va = integrant_un.va.mult_sp() ;
256 integrant_deux.va = integrant_deux.va.mult_st() ;
257 integrant_deux.va = integrant_deux.va.mult_cp() ;
259 double moment = mapping.integrale_surface_infini (-integrant_un+integrant_deux) ;
266double Bin_ns_bh::moment_systeme_hor()
const {
268 using namespace Unites ;
274 Cmp xa_bh (
hole.mp) ;
276 xa_bh.std_base_scal() ;
278 Cmp ya_bh (
hole.mp) ;
280 ya_bh.std_base_scal() ;
282 Tenseur vecteur_bh (
hole.mp, 1, CON,
hole.mp.get_bvect_cart()) ;
283 vecteur_bh.set_etat_qcq() ;
284 for (
int i=0 ; i<3 ; i++)
285 vecteur_bh.set(i) = (-ya_bh*
hole.tkij_tot(0, i)+
286 xa_bh *
hole.tkij_tot(1, i)) ;
287 vecteur_bh.set_std_base() ;
288 vecteur_bh.annule(
hole.mp.get_mg()->get_nzone()-1) ;
289 vecteur_bh.change_triad (
hole.mp.get_bvect_spher()) ;
291 Cmp integrant_bh (
pow(
hole.psi_tot(), 6)*vecteur_bh(0)) ;
292 integrant_bh.std_base_scal() ;
293 double moment_bh =
hole.mp.integrale_surface
294 (integrant_bh,
hole.rayon)/8/M_PI ;
297 Cmp xa_ns (
star.mp) ;
299 xa_ns.std_base_scal() ;
301 Cmp ya_ns (
star.mp) ;
303 ya_ns.std_base_scal() ;
305 Cmp integrant_ns (
pow(
star.confpsi_auto()+
star.confpsi_comp(), 10)*(
star.ener_euler()+
star.press())*
306 (xa_ns*
star.u_euler(1) - ya_ns*
star.u_euler(0))) ;
307 integrant_ns.std_base_scal() ;
309 double moment_ns = integrant_ns.integrale() * ggrav ;
310 return moment_ns + moment_bh ;
314Tbl Bin_ns_bh::linear_momentum_systeme_inf()
const {
320 int nzones =
hole.mp.get_mg()->get_nzone() ;
321 double* bornes =
new double [nzones+1] ;
322 double courant = fabs(
hole.mp.get_ori_x()-
star.mp.get_ori_x())+1 ;
323 for (
int i=nzones-1 ; i>0 ; i--) {
324 bornes[i] = courant ;
328 bornes[nzones] = __infinity ;
330 Map_af mapping (*
hole.mp.get_mg(), bornes) ;
335 Tenseur_sym k_total (mapping, 2, CON, mapping.get_bvect_cart()) ;
336 k_total.set_etat_qcq() ;
338 Tenseur shift_un (mapping, 1, CON, mapping.get_bvect_cart()) ;
339 shift_un.set_etat_qcq() ;
341 Tenseur shift_deux (mapping, 1, CON, mapping.get_bvect_cart()) ;
342 shift_deux.set_etat_qcq() ;
344 shift_un.set_triad (*
hole.shift_auto.get_triad()) ;
345 shift_un.set(0).import(
hole.shift_auto(0)) ;
346 shift_un.set(1).import(
hole.shift_auto(1)) ;
347 shift_un.set(2).import(
hole.shift_auto(2)) ;
348 shift_un.change_triad (mapping.get_bvect_cart()) ;
350 shift_deux.set_triad (*
star.shift_auto.get_triad()) ;
351 shift_deux.set(0).import(
star.shift_auto(0)) ;
352 shift_deux.set(1).import(
star.shift_auto(1)) ;
353 shift_deux.set(2).import(
star.shift_auto(2)) ;
354 shift_deux.change_triad(mapping.get_bvect_cart()) ;
356 shift_un.set_std_base() ;
357 shift_deux.set_std_base() ;
359 Tenseur shift_tot (shift_un+shift_deux) ;
360 shift_tot.set_std_base() ;
361 shift_tot.annule(0, nzones-2) ;
363 Cmp compy (shift_tot(1)) ;
366 int nr = mapping.get_mg()->get_nr(nzones-1) ;
367 int nt = mapping.get_mg()->get_nt(nzones-1) ;
368 int np = mapping.get_mg()->get_np(nzones-1) ;
369 Tbl val_inf (nt*np) ;
370 val_inf.set_etat_qcq() ;
371 for (
int k=0 ; k<np ; k++)
372 for (
int j=0 ; j<nt ; j++)
373 val_inf.set(k*nt + j) = fabs(compy (nzones-1, k, j, nr-1)) ;
375 Tenseur grad (shift_tot.gradient()) ;
376 Tenseur trace (grad(0, 0)+grad(1, 1)+grad(2, 2)) ;
377 for (
int i=0 ; i<3 ; i++) {
378 k_total.set(i, i) = grad(i, i)-trace()/3. ;
379 for (
int j=i+1 ; j<3 ; j++)
380 k_total.set(i, j) = (grad(i, j)+grad(j, i))/2. ;
383 for (
int comp=0 ; comp<3 ; comp++) {
384 Tenseur vecteur (mapping, 1, CON, mapping.get_bvect_cart()) ;
385 vecteur.set_etat_qcq() ;
386 for (
int i=0 ; i<3 ; i++)
387 vecteur.set(i) = k_total(i, comp) ;
388 vecteur.change_triad (mapping.get_bvect_spher()) ;
389 Cmp integrant (vecteur(0)) ;
391 res.set(comp) = mapping.integrale_surface_infini (integrant)/8/M_PI ;
396double Bin_ns_bh::distance_propre_axe_bh (
const int nr)
const {
398 double x_bh =
hole.mp.get_ori_x() +
hole.rayon ;
401 double pente = -2./x_bh ;
402 double constante = - 1. ;
406 double air_un, theta_un, phi_un ;
407 double air_deux, theta_deux, phi_deux ;
409 double* coloc =
new double[nr] ;
410 double* coef =
new double[nr] ;
411 int* deg =
new int[3] ;
412 deg[0] = 1 ; deg[1] = 1 ; deg[2] = nr ;
414 for (
int i=0 ; i<nr ; i++) {
415 ksi = -
cos (M_PI*i/(nr-1)) ;
416 xabs = (ksi+constante)/pente ;
418 hole.mp.convert_absolute (xabs, 0, 0, air_un, theta_un, phi_un) ;
419 star.mp.convert_absolute (xabs, 0, 0, air_deux, theta_deux, phi_deux) ;
421 coloc[i] =
pow (
hole.psi_auto().val_point (air_un, theta_un, phi_un) +
422 star.confpsi_auto().val_point (air_deux, theta_deux, phi_deux), 2.) ;
426 cfrcheb(deg, deg, coloc, deg, coef) ;
429 double* som =
new double[nr] ;
431 for (
int i=2 ; i<nr ; i+=2)
432 som[i] = 1./(i+1)-1./(i-1) ;
433 for (
int i=1 ; i<nr ; i+=2)
437 for (
int i=0 ; i<nr ; i++)
438 res += som[i]*coef[i] ;
450double Bin_ns_bh::distance_propre_axe_ns (
const int nr)
const {
452 double x_ns =
star.mp.get_ori_x() -
star.mp.val_r (
star.nzet, -1, M_PI/2, M_PI) ;
455 double pente = 2./x_ns ;
456 double constante = - 1. ;
460 double air_un, theta_un, phi_un ;
461 double air_deux, theta_deux, phi_deux ;
463 double* coloc =
new double[nr] ;
464 double* coef =
new double[nr] ;
465 int* deg =
new int[3] ;
466 deg[0] = 1 ; deg[1] = 1 ; deg[2] = nr ;
468 for (
int i=0 ; i<nr ; i++) {
469 ksi = -
cos (M_PI*i/(nr-1)) ;
470 xabs = (ksi-constante)/pente ;
472 hole.mp.convert_absolute (xabs, 0, 0, air_un, theta_un, phi_un) ;
473 star.mp.convert_absolute (xabs, 0, 0, air_deux, theta_deux, phi_deux) ;
475 coloc[i] =
pow (
hole.psi_auto().val_point (air_un, theta_un, phi_un) +
476 star.confpsi_auto().val_point (air_deux, theta_deux, phi_deux), 2.) ;
480 cfrcheb(deg, deg, coloc, deg, coef) ;
483 double* som =
new double[nr] ;
485 for (
int i=2 ; i<nr ; i+=2)
486 som[i] = 1./(i+1)-1./(i-1) ;
487 for (
int i=1 ; i<nr ; i+=2)
491 for (
int i=0 ; i<nr ; i++)
492 res += som[i]*coef[i] ;
504double Bin_ns_bh::smarr()
const {
506 using namespace Unites ;
509 Tenseur psiq_t (
pow(
star.get_confpsi()(), 4.)) ;
510 psiq_t.set_std_base() ;
512 Tenseur_sym furmet (
star.mp, 2, CON,
star.mp.get_bvect_cart()) ;
513 furmet.set_etat_qcq() ;
514 for (
int i=0 ; i< 3 ; i++) {
515 furmet.set(i,i) = 1/psiq_t() ;
516 for(
int j=i+1 ; j<3 ; j++)
517 furmet.set(i,j) = 0 ;
519 Metrique met (furmet,
false) ;
521 Tenseur_sym kij (
star.get_tkij_tot()/psiq_t) ;
522 kij.change_triad(
star.mp.get_bvect_cart()) ;
524 Tenseur_sym kij_cov (
manipule (kij, met)) ;
525 Tenseur shift (
star.get_shift()) ;
526 shift.change_triad(
star.mp.get_bvect_cart()) ;
528 Tenseur aime (
star.mp, 1, CON,
star.mp.get_bvect_cart()) ;
529 aime.set_etat_qcq() ;
533 aime.annule(
star.mp.get_mg()->get_nzone()-1) ;
534 aime.set_std_base() ;
535 shift = shift - aime ;
538 Tenseur u_euler (
star.get_u_euler()) ;
539 u_euler.change_triad(
star.mp.get_bvect_cart()) ;
540 Tenseur u_i_bas (
manipule(u_euler, met)) ;
541 Tenseur mat (qpig*(
star.get_nnn()*(
star.get_ener_euler() +
star.get_s_euler()) - 2*(
star.get_ener_euler()+
star.get_press())*
contract(u_i_bas, 0, shift, 0))) ;
544 Cmp psiq (
pow(
star.get_confpsi()(), 4.)) ;
546 Cmp integ_matter (
star.get_nnn()()*(
star.get_ener_euler()() +
star.get_s_euler()())
548 integ_matter = integ_matter *
pow(
star.get_confpsi()(),6.) ;
549 integ_matter.std_base_scal() ;
550 integ_matter.annule(
star.get_nzet(),
star.get_mp().get_mg()->get_nzone()-1) ;
551 double matter_term = integ_matter.integrale()*qpig/4/M_PI ;
558 Tenseur rad (
hole.mp, 1, COV,
hole.mp.get_bvect_cart()) ;
567 for (
int m=0 ; m<3 ; m++)
568 for (
int n=0 ; n<3 ; n++)
569 temp += rad(m)*rad(n)*
hole.tkij_tot(m,n) ;
570 temp *=
pow(
hole.psi_tot(),2.) ;
571 temp.std_base_scal() ;
573 Cmp integ_hor ((
hole.get_n_tot()().dsdr()-
hole.get_n_tot()()*temp)
574 *
pow(
hole.get_psi_tot()(), 2)) ;
575 integ_hor.std_base_scal() ;
576 integ_hor.raccord(1) ;
577 double hor_term =
hole.mp.integrale_surface(integ_hor,
hole.get_rayon()) ;
580 double m_test = hor_term + matter_term + 2*
omega*moment_systeme_inf() +
Et_bin_nsbh star
The neutron star.
Bhole hole
The black hole.
double omega
Angular velocity with respect to an asymptotically inertial observer.
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Cmp sin(const Cmp &)
Sine.
Cmp pow(const Cmp &, int)
Power .
Cmp cos(const Cmp &)
Cosine.
Tenseur flat_scalar_prod(const Tenseur &t1, const Tenseur &t2)
Scalar product of two Tenseur when the metric is : performs the contraction of the last index of t1 w...
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Tenseur manipule(const Tenseur &, const Metrique &, int idx)
Raise or lower the index idx depending on its type, using the given Metrique .
Coord phi
coordinate centered on the grid
Coord tet
coordinate centered on the grid