LORENE
bin_ns_bh_glob.C
1 /*
2  * Copyright (c) 2005 Philippe Grandclement
3  *
4  * This file is part of LORENE.
5  *
6  * LORENE is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * LORENE is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with LORENE; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  *
20  */
21 
22 
23 
24 
25 /*
26  * $Id: bin_ns_bh_glob.C,v 1.9 2016/12/05 16:17:46 j_novak Exp $
27  * $Log: bin_ns_bh_glob.C,v $
28  * Revision 1.9 2016/12/05 16:17:46 j_novak
29  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
30  *
31  * Revision 1.8 2014/10/13 08:52:43 j_novak
32  * Lorene classes and functions now belong to the namespace Lorene.
33  *
34  * Revision 1.7 2014/10/06 15:13:01 j_novak
35  * Modified #include directives to use c++ syntax.
36  *
37  * Revision 1.6 2007/04/24 20:13:53 f_limousin
38  * Implementation of Dirichlet and Neumann BC for the lapse
39  *
40  * Revision 1.5 2006/06/23 07:09:24 p_grandclement
41  * Addition of spinning black hole
42  *
43  * Revision 1.4 2006/06/01 12:47:52 p_grandclement
44  * update of the Bin_ns_bh project
45  *
46  * Revision 1.3 2005/12/01 12:59:10 p_grandclement
47  * Files for bin_ns_bh project
48  *
49  * Revision 1.2 2005/11/30 11:09:06 p_grandclement
50  * Changes for the Bin_ns_bh project
51  *
52  * Revision 1.1 2005/08/29 15:10:15 p_grandclement
53  * Addition of things needed :
54  * 1) For BBH with different masses
55  * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
56  * WORKING YET !!!
57  *
58  *
59  * $Header: /cvsroot/Lorene/C++/Source/Bin_ns_bh/bin_ns_bh_glob.C,v 1.9 2016/12/05 16:17:46 j_novak Exp $
60  *
61  */
62 
63 
64 
65 //standard
66 #include <cstdlib>
67 #include <cmath>
68 
69 // Lorene
70 #include "nbr_spx.h"
71 #include "tenseur.h"
72 #include "bhole.h"
73 #include "bin_ns_bh.h"
74 #include "proto.h"
75 #include "utilitaires.h"
76 #include "graphique.h"
77 #include "unites.h"
78 #include "metrique.h"
79 
80 namespace Lorene {
81 double Bin_ns_bh::adm_systeme() const {
82  Cmp der_un (hole.psi_auto().dsdr()) ;
83 
84  Map_af auxi_mp (star.get_mp()) ;
85  Cmp der_deux (star.confpsi_auto().dsdr()) ;
86 
87  double masse = hole.mp.integrale_surface_infini(der_un) +
88  auxi_mp.integrale_surface_infini(der_deux) ;
89 
90  masse /= -2*M_PI ;
91  return masse ;
92 }
93 
94 double Bin_ns_bh::adm_systeme_volume() const {
95 
96  using namespace Unites ;
97 
98  Tenseur auxi_bh (flat_scalar_prod(hole.tkij_tot, hole.tkij_auto)) ;
99  Tenseur kk_bh (hole.mp) ;
100  kk_bh = 0 ;
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 ;
106  }
107  Cmp integ_bh (pow(hole.psi_tot(), 5.)*kk_bh()) ;
108  integ_bh.annule(0) ;
109  integ_bh.std_base_scal() ;
110 
111  Cmp integ_hor1 (hole.psi_tot()) ;
112 
113  Cmp tet (hole.mp) ;
114  tet = hole.mp.tet ;
115  Cmp phi (hole.mp) ;
116  phi = hole.mp.phi ;
117  Tenseur rad (hole.mp, 1, COV, hole.mp.get_bvect_cart()) ;
118  rad.set_etat_qcq() ;
119  rad.set(0) = cos(phi)*sin(tet) ;
120  rad.set(1) = sin(phi)*sin(tet) ;
121  rad.set(2) = cos(tet) ;
122 
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. ;
130  integ_hor2.std_base_scal() ;
131 
132  Tenseur auxi_ns (flat_scalar_prod(star.tkij_tot, star.tkij_auto)) ;
133  Tenseur kk_ns (star.mp) ;
134  kk_ns = 0 ;
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 ;
140  }
141  Cmp integ_ns (pow(star.confpsi_comp() + star.confpsi_auto(), 5.)*kk_ns()) ;
142  integ_ns.std_base_scal() ;
143 
144  Cmp integ_matter (pow(star.confpsi_comp()+star.confpsi_auto(), 5.)*star.ener_euler()) ;
145  integ_matter.std_base_scal() ;
146 
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 ;
151 
152  return masse ;
153 }
154 
155 double Bin_ns_bh::komar_systeme() const {
156  Cmp der_un (hole.n_auto().dsdr()) ;
157 
158  Map_af auxi_mp (star.get_mp()) ;
159  Cmp der_deux (star.n_auto().dsdr()) ;
160 
161  double masse = hole.mp.integrale_surface_infini(der_un) +
162  auxi_mp.integrale_surface_infini(der_deux) ;
163 
164  masse /= 4*M_PI ;
165 
166  return masse ;
167 }
168 
169 double Bin_ns_bh::viriel() const {
170  double adm = adm_systeme() ;
171  double komar = komar_systeme() ;
172 
173  return (adm-komar)/adm ;
174 }
175 
176 double Bin_ns_bh::moment_systeme_inf() const {
177 
178  if (omega == 0)
179  return 0 ;
180  else {
181 
182  // On construit une grille et un mapping auxiliaire :
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 ;
188  courant /= 2. ;
189  }
190  bornes[0] = 0 ;
191  bornes[nzones] = __infinity ;
192 
193  Map_af mapping (*hole.mp.get_mg(), bornes) ;
194 
195  delete [] bornes ;
196 
197  // On construit k_total dessus :
198  Tenseur_sym k_total (mapping, 2, CON, mapping.get_bvect_cart()) ;
199  k_total.set_etat_qcq() ;
200 
201  Tenseur shift_un (mapping, 1, CON, mapping.get_bvect_cart()) ;
202  shift_un.set_etat_qcq() ;
203 
204  Tenseur shift_deux (mapping, 1, CON, mapping.get_bvect_cart()) ;
205  shift_deux.set_etat_qcq() ;
206 
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()) ;
212 
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()) ;
218 
219  Tenseur shift_tot (shift_un+shift_deux) ;
220  shift_tot.set_std_base() ;
221  shift_tot.annule(0, nzones-2) ;
222  // On enleve les residus
223  shift_tot.inc2_dzpuis() ;
224  shift_tot.dec2_dzpuis() ;
225 
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. ;
232  }
233 
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() ;
237 
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)) ;
244 
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)) ;
251 
252  // On multiplie par y et x :
253  integrant_un.va = integrant_un.va.mult_st() ;
254  integrant_un.va = integrant_un.va.mult_sp() ;
255 
256  integrant_deux.va = integrant_deux.va.mult_st() ;
257  integrant_deux.va = integrant_deux.va.mult_cp() ;
258 
259  double moment = mapping.integrale_surface_infini (-integrant_un+integrant_deux) ;
260 
261  moment /= 8*M_PI ;
262  return moment ;
263  }
264 }
265 
266 double Bin_ns_bh::moment_systeme_hor() const {
267 
268  using namespace Unites ;
269 
270  if (omega == 0)
271  return 0 ;
272  else {
273  //Contribution du trou noir :
274  Cmp xa_bh (hole.mp) ;
275  xa_bh = hole.mp.xa ;
276  xa_bh.std_base_scal() ;
277 
278  Cmp ya_bh (hole.mp) ;
279  ya_bh = hole.mp.ya ;
280  ya_bh.std_base_scal() ;
281 
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()) ;
290 
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 ;
295 
296  // Contribution de l'étoile :
297  Cmp xa_ns (star.mp) ;
298  xa_ns = star.mp.xa ;
299  xa_ns.std_base_scal() ;
300 
301  Cmp ya_ns (star.mp) ;
302  ya_ns = star.mp.ya ;
303  ya_ns.std_base_scal() ;
304 
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() ;
308 
309  double moment_ns = integrant_ns.integrale() * ggrav ;
310  return moment_ns + moment_bh ;
311  }
312 }
313 
314 Tbl Bin_ns_bh::linear_momentum_systeme_inf() const {
315 
316  Tbl res (3) ;
317  res.set_etat_qcq() ;
318 
319  // On construit une grille et un mapping auxiliaire :
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 ;
325  courant /= 2. ;
326  }
327  bornes[0] = 0 ;
328  bornes[nzones] = __infinity ;
329 
330  Map_af mapping (*hole.mp.get_mg(), bornes) ;
331 
332  delete [] bornes ;
333 
334  // On construit k_total dessus :
335  Tenseur_sym k_total (mapping, 2, CON, mapping.get_bvect_cart()) ;
336  k_total.set_etat_qcq() ;
337 
338  Tenseur shift_un (mapping, 1, CON, mapping.get_bvect_cart()) ;
339  shift_un.set_etat_qcq() ;
340 
341  Tenseur shift_deux (mapping, 1, CON, mapping.get_bvect_cart()) ;
342  shift_deux.set_etat_qcq() ;
343 
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()) ;
349 
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()) ;
355 
356  shift_un.set_std_base() ;
357  shift_deux.set_std_base() ;
358 
359  Tenseur shift_tot (shift_un+shift_deux) ;
360  shift_tot.set_std_base() ;
361  shift_tot.annule(0, nzones-2) ;
362 
363  Cmp compy (shift_tot(1)) ;
364  compy.mult_r_zec() ;
365 
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)) ;
374 
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. ;
381  }
382 
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)) ;
390 
391  res.set(comp) = mapping.integrale_surface_infini (integrant)/8/M_PI ;
392  }
393  return res ;
394 }
395 
396 double Bin_ns_bh::distance_propre_axe_bh (const int nr) const {
397 
398  double x_bh = hole.mp.get_ori_x() + hole.rayon ;
399 
400  // Les coefficients du changement de variable :
401  double pente = -2./x_bh ;
402  double constante = - 1. ;
403 
404  double ksi ; // variable d'integration.
405  double xabs ; // x reel.
406  double air_un, theta_un, phi_un ; // coordonnee spheriques 1
407  double air_deux, theta_deux, phi_deux ; // coordonnee spheriques 2
408 
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 ;
413 
414  for (int i=0 ; i<nr ; i++) {
415  ksi = -cos (M_PI*i/(nr-1)) ;
416  xabs = (ksi+constante)/pente ;
417 
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) ;
420 
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.) ;
423  }
424 
425  // On prend les coefficients de la fonction
426  cfrcheb(deg, deg, coloc, deg, coef) ;
427 
428  // On integre
429  double* som = new double[nr] ;
430  som[0] = 2 ;
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)
434  som[i] = 0 ;
435 
436  double res = 0 ;
437  for (int i=0 ; i<nr ; i++)
438  res += som[i]*coef[i] ;
439 
440  res /= pente ;
441 
442  delete [] deg ;
443  delete [] coef ;
444  delete [] coloc ;
445  delete [] som ;
446 
447  return res ;
448 }
449 
450 double Bin_ns_bh::distance_propre_axe_ns (const int nr) const {
451 
452  double x_ns = star.mp.get_ori_x() - star.mp.val_r (star.nzet, -1, M_PI/2, M_PI) ;
453 
454  // Les coefficients du changement de variable :
455  double pente = 2./x_ns ;
456  double constante = - 1. ;
457 
458  double ksi ; // variable d'integration.
459  double xabs ; // x reel.
460  double air_un, theta_un, phi_un ; // coordonnee spheriques 1
461  double air_deux, theta_deux, phi_deux ; // coordonnee spheriques 2
462 
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 ;
467 
468  for (int i=0 ; i<nr ; i++) {
469  ksi = -cos (M_PI*i/(nr-1)) ;
470  xabs = (ksi-constante)/pente ;
471 
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) ;
474 
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.) ;
477  }
478 
479  // On prend les coefficients de la fonction
480  cfrcheb(deg, deg, coloc, deg, coef) ;
481 
482  // On integre
483  double* som = new double[nr] ;
484  som[0] = 2 ;
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)
488  som[i] = 0 ;
489 
490  double res = 0 ;
491  for (int i=0 ; i<nr ; i++)
492  res += som[i]*coef[i] ;
493 
494  res /= pente ;
495 
496  delete [] deg ;
497  delete [] coef ;
498  delete [] coloc ;
499  delete [] som ;
500 
501  return res ;
502 }
503 
504 double Bin_ns_bh::smarr() const {
505 
506  using namespace Unites ;
507 
508  // The tests
509  Tenseur psiq_t (pow(star.get_confpsi()(), 4.)) ;
510  psiq_t.set_std_base() ;
511 
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 ;
518  }
519  Metrique met (furmet, false) ;
520 
521  Tenseur_sym kij (star.get_tkij_tot()/psiq_t) ;
522  kij.change_triad(star.mp.get_bvect_cart()) ;
523  kij.dec2_dzpuis() ;
524  Tenseur_sym kij_cov (manipule (kij, met)) ;
525  Tenseur shift (star.get_shift()) ;
526  shift.change_triad(star.mp.get_bvect_cart()) ;
527 
528  Tenseur aime (star.mp, 1, CON, star.mp.get_bvect_cart()) ;
529  aime.set_etat_qcq() ;
530  aime.set(0) = -omega*star.mp.ya ;
531  aime.set(1) = omega*star.mp.xa ;
532  aime.set(2) = 0 ;
533  aime.annule(star.mp.get_mg()->get_nzone()-1) ;
534  aime.set_std_base() ;
535  shift = shift - aime ;
536 
537  // La matière :
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))) ;
542 
543  // La partie avec la matière :
544  Cmp psiq (pow(star.get_confpsi()(), 4.)) ;
545 
546  Cmp integ_matter (star.get_nnn()()*(star.get_ener_euler()() + star.get_s_euler()())
547  - 2*(star.get_ener_euler()()+star.get_press()())*psiq*flat_scalar_prod(u_euler, shift)()) ;
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 ;
552 
553  // Integrale sur horizon :
554  Cmp tet (hole.mp) ;
555  tet = hole.mp.tet ;
556  Cmp phi (hole.mp) ;
557  phi = hole.mp.phi ;
558  Tenseur rad (hole.mp, 1, COV, hole.mp.get_bvect_cart()) ;
559  rad.set_etat_qcq() ;
560  rad.set(0) = cos(phi)*sin(tet) ;
561  rad.set(1) = sin(phi)*sin(tet) ;
562  rad.set(2) = cos(tet) ;
563 
564  Cmp temp (hole.mp) ;
565  temp.annule_hard() ;
566  temp.set_dzpuis(2) ;
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() ;
572 
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()) ;
578  hor_term /= 4*M_PI ;
579 
580  double m_test = hor_term + matter_term + 2*omega*moment_systeme_inf() +
581  2*(hole.omega_local-omega)*hole.local_momentum() ;
582 
583  return m_test ;
584  }
585 }
Coord xa
Absolute x coordinate.
Definition: map.h:748
double integrale_surface(const Cmp &ci, double rayon) const
Performs the surface integration of ci on the sphere of radius rayon .
const Map & get_mp() const
Returns the mapping.
Definition: etoile.h:662
double integrale_surface_infini(const Cmp &ci) const
Performs the surface integration of ci at infinity.
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:801
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:783
const Tenseur & get_shift() const
Returns the total shift vector .
Definition: etoile.h:733
double get_ori_x() const
Returns the x coordinate of the origin.
Definition: map.h:786
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 confpsi_comp
Part of the conformal factor $$ generated principaly by the companion star.
Definition: et_bin_nsbh.h:109
const Tenseur & get_psi_tot() const
Returns the total .
Definition: bhole.h:424
int get_nzet() const
Returns the number of domains occupied by the star.
Definition: etoile.h:665
double omega_local
local angular velocity
Definition: bhole.h:276
Tenseur psi_auto
Part of generated by the hole.
Definition: bhole.h:290
const Tenseur & get_n_tot() const
Returns the total N .
Definition: bhole.h:407
Tenseur press
Fluid pressure.
Definition: etoile.h:464
void set(const Map *mp, Mtbl *(*construct)(const Map *))
Semi-constructor from a mapping and a method.
Definition: coord.C:137
Tenseur n_auto
Part of the lapse { N} generated principaly by the star.
Definition: et_bin_nsbh.h:85
Cmp cos(const Cmp &)
Cosine.
Definition: cmp_math.C:97
Coord tet
coordinate centered on the grid
Definition: map.h:737
Tenseur shift_auto
Part of generated by the hole.
Definition: bhole.h:297
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition: etoile.h:892
Coord phi
coordinate centered on the grid
Definition: map.h:738
const Tenseur & get_press() const
Returns the fluid pressure.
Definition: etoile.h:685
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:477
Et_bin_nsbh star
The neutron star.
Definition: bin_ns_bh.h:131
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition: tenseur.h:707
Map_af & mp
Affine mapping.
Definition: bhole.h:273
const Tenseur & get_s_euler() const
Returns the trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:691
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
double get_rayon() const
Returns the radius of the horizon.
Definition: bhole.h:347
const Tenseur & get_nnn() const
Returns the total lapse function N.
Definition: etoile.h:730
double omega
Angular velocity with respect to an asymptotically inertial observer.
Definition: bin_ns_bh.h:139
const Tenseur & get_confpsi() const
Returns the part of the conformal factor $$.
Definition: et_bin_nsbh.h:257
double rayon
Radius of the horizon in LORENE&#39;s units.
Definition: bhole.h:274
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Tenseur_sym tkij_auto
Auto .
Definition: bhole.h:307
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
const Tenseur & get_ener_euler() const
Returns the total energy density with respect to the Eulerian observer.
Definition: etoile.h:688
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition: cmp.C:647
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:435
Tenseur_sym tkij_tot
Total .
Definition: bhole.h:308
Coord ya
Absolute y coordinate.
Definition: map.h:749
Tenseur confpsi_auto
Part of the conformal factor $$ generated principaly by the star.
Definition: et_bin_nsbh.h:104
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:809
const Tenseur_sym & get_tkij_tot() const
Returns the total extrinsic curvature tensor $K^{ij}$ generated by { shift_auto} and { shift_comp}...
Definition: et_bin_nsbh.h:308
Tenseur psi_tot
Total .
Definition: bhole.h:292
Tenseur_sym tkij_tot
Total extrinsic curvature tensor $K^{ij}$ generated by { shift_auto} and { shift_comp}.
Definition: et_bin_nsbh.h:152
Bhole hole
The black hole.
Definition: bin_ns_bh.h:134
Cmp sin(const Cmp &)
Sine.
Definition: cmp_math.C:72
Tenseur n_auto
Part of N generated by the hole.
Definition: bhole.h:286
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor $K^{ij}$ generated by { shift_auto}.
Definition: et_bin_nsbh.h:146
const Tenseur & get_u_euler() const
Returns the fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:697
void convert_absolute(double xx, double yy, double zz, double &rr, double &theta, double &pphi) const
Determines the coordinates corresponding to given absolute Cartesian coordinates (X...
Definition: map.C:305
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:468
Tenseur manipule(const Tenseur &, const Metrique &, int idx)
Raise or lower the index idx depending on its type, using the given Metrique .