LORENE
et_bin_bhns_extr_velpot.C
1/*
2 * Method of class Et_bin_bhns_extr to compute the velocity scalar
3 * potential $\psi$ by solving the continuity equation
4 * in the Kerr-Schild background metric or in the conformally flat one
5 * with extreme mass ratio
6 *
7 * (see file et_bin_bhns_extr.h for documentation).
8 *
9 */
10
11/*
12 * Copyright (c) 2004-2005 Keisuke Taniguchi
13 *
14 * This file is part of LORENE.
15 *
16 * LORENE is free software; you can redistribute it and/or modify
17 * it under the terms of the GNU General Public License version 2
18 * as published by the Free Software Foundation.
19 *
20 * LORENE is distributed in the hope that it will be useful,
21 * but WITHOUT ANY WARRANTY; without even the implied warranty of
22 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23 * GNU General Public License for more details.
24 *
25 * You should have received a copy of the GNU General Public License
26 * along with LORENE; if not, write to the Free Software
27 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28 *
29 */
30
31char et_bin_bhns_extr_velpot_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_velpot.C,v 1.3 2014/10/13 08:52:55 j_novak Exp $" ;
32
33/*
34 * $Id: et_bin_bhns_extr_velpot.C,v 1.3 2014/10/13 08:52:55 j_novak Exp $
35 * $Log: et_bin_bhns_extr_velpot.C,v $
36 * Revision 1.3 2014/10/13 08:52:55 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.2 2005/02/28 23:17:18 k_taniguchi
40 * Modification to include the case of the conformally flat background metric
41 *
42 * Revision 1.1 2004/11/30 20:51:56 k_taniguchi
43 * *** empty log message ***
44 *
45 *
46 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_velpot.C,v 1.3 2014/10/13 08:52:55 j_novak Exp $
47 *
48 */
49
50// C headers
51//#include <math.h>
52
53// Lorene headers
54#include "et_bin_bhns_extr.h"
55#include "scalar.h"
56#include "metrique.h"
57#include "etoile.h"
58#include "eos.h"
59#include "param.h"
60#include "coord.h"
61#include "unites.h"
62
63// Local prototype
64namespace Lorene {
65Cmp raccord_c1(const Cmp& uu, int l1) ;
66
67double Et_bin_bhns_extr::velocity_pot_extr(const double& mass,
68 const double& sepa,
69 int mermax, double precis,
70 double relax) {
71
72 using namespace Unites ;
73
74 if (kerrschild) {
75
76 if (eos.identify() == 5 || eos.identify() == 4 ||
77 eos.identify() == 3) {
78
79 cout << "Etoile_bin::velocity_pot_extr" << endl ;
80 cout << "The code has not prepared for this kind of EOS!!!"
81 << endl;
82 abort() ;
83
84 } // End of strange stars case
85 else {
86
87 int nzm1 = mp.get_mg()->get_nzone() - 1 ;
88
89 //--------------------------------
90 // Specific relativistic enthalpy ---> hhh
91 //--------------------------------
92
93 Tenseur hhh = exp(unsurc2 * ent) ; // = 1 at the Newtonian limit
94 hhh.set_std_base() ;
95
96 //--------------------------------------------
97 // Computation of W^i = - A^2 h Gamma_n B^i/N
98 //--------------------------------------------
99
100 Tenseur www = - a_car * hhh * gam_euler * bsn ;
101
103 // components on the mapping Cartesian basis
104
105 //-------------------------------------------------
106 // Constant value of W^i at the center of the star
107 //-------------------------------------------------
108
109 Tenseur v_orb(mp, 1, CON, mp.get_bvect_cart()) ;
110
111 v_orb.set_etat_qcq() ;
112 for (int i=0; i<3; i++) {
113 v_orb.set(i) = www(i)(0, 0, 0, 0) ;
114 }
115
116 v_orb.set_triad( *(www.get_triad()) ) ;
117 v_orb.set_std_base() ;
118
119 // Some auxiliary terms
120 // --------------------
121
122 const Coord& xx = mp.x ;
123 const Coord& yy = mp.y ;
124 const Coord& zz = mp.z ;
125
126 Tenseur r_bh(mp) ;
127 r_bh.set_etat_qcq() ;
128 r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
129 r_bh.set_std_base() ;
130
131 Tenseur msr(mp) ;
132 msr = ggrav * mass / r_bh ;
133 msr.set_std_base() ;
134
135 Tenseur lapse_bh2(mp) ; // lapse_bh * lapse_bh
136 lapse_bh2 = 1. / (1.+2.*msr) ;
137 lapse_bh2.set_std_base() ;
138
139 Tenseur lapse_bh3(mp) ; // lapse_bh * lapse_bh * lapse_bh
140 lapse_bh3 = 1./pow(1.+2.*msr, 1.5) ;
141 lapse_bh3.set_std_base() ;
142
143 Tenseur xx_con(mp, 1, CON, mp.get_bvect_cart()) ;
144 xx_con.set_etat_qcq() ;
145 xx_con.set(0) = xx + sepa ;
146 xx_con.set(1) = yy ;
147 xx_con.set(2) = zz ;
148 xx_con.set_std_base() ;
149
150 Tenseur xsr_con(mp, 1, CON, mp.get_bvect_cart()) ;
151 xsr_con = xx_con / r_bh ;
152 xsr_con.set_std_base() ;
153
154 Tenseur xx_cov(mp, 1, COV, mp.get_bvect_cart()) ;
155 xx_cov.set_etat_qcq() ;
156 xx_cov.set(0) = xx + sepa ;
157 xx_cov.set(1) = yy ;
158 xx_cov.set(2) = zz ;
159 xx_cov.set_std_base() ;
160
161 Tenseur xsr_cov(mp, 1, COV, mp.get_bvect_cart()) ;
162 xsr_cov = xx_cov / r_bh ;
163 xsr_cov.set_std_base() ;
164
165 // X_i/r_bh in the spherical coordinate
166 const Coord& rr = mp.r ;
167 const Coord& st = mp.sint ;
168 const Coord& ct = mp.cost ;
169 const Coord& sp = mp.sinp ;
170 const Coord& cp = mp.cosp ;
171
172 Tenseur rr_spher_cov(mp, 1, COV, mp.get_bvect_spher()) ;
173 rr_spher_cov.set_etat_qcq() ;
174 rr_spher_cov.set(0) = rr + sepa*st*cp ;
175 rr_spher_cov.set(1) = sepa*ct*cp ;
176 rr_spher_cov.set(2) = - sepa*sp ;
177 rr_spher_cov.set_std_base() ;
178
179 Tenseur xsr_spher_cov(mp, 1, COV, mp.get_bvect_spher()) ;
180 xsr_spher_cov = rr_spher_cov / r_bh ;
181 xsr_spher_cov.set_std_base() ;
182
183 // l^j D_j H_ent
184 Tenseur ldent = flat_scalar_prod(xsr_con, ent.gradient()) ;
185
186 // l^j D_j beta_auto
187 Tenseur ldbeta = flat_scalar_prod(xsr_con, beta_auto.gradient()) ;
188
189 // l^j D_j psi0
190 Tenseur ldpsi0 = flat_scalar_prod(xsr_con, psi0.gradient()) ;
191
192 // l^i D_i (l^j D_j psi0)
193 Tenseur ldldpsi0 = flat_scalar_prod(xsr_con, ldpsi0.gradient()) ;
194
195 //-------------------------------------------------
196 // Source and coefficients a,b for poisson_compact
197 // (idenpendent from psi0)
198 //-------------------------------------------------
199
200 Cmp dndh_log = eos.der_nbar_ent(ent(), nzet) ;
201
202 // In order to avoid any division by zero in the computation of
203 // zeta_h the value of dndh_log is set to 1
204 // in the external domains:
205 for (int l=nzet; l <= nzm1; l++) {
206 dndh_log.set(l) = 1 ;
207 }
208
209 double erreur ;
210
211 Tenseur zeta_h( ent() / dndh_log ) ;
212 zeta_h.set_std_base() ;
213
214 Tenseur tmp_zeta = 1 - unsurc2 * zeta_h ;
215 tmp_zeta.set_std_base() ;
216
217 Tenseur bb = tmp_zeta * (ent.gradient_spher()
218 - 2.*lapse_bh2 * msr * ldent
219 * xsr_spher_cov)
220 + unsurc2 * zeta_h * (beta_auto.gradient_spher()
221 - 2.*lapse_bh2 * msr * ldbeta
222 * xsr_spher_cov)
223 - unsurc2 * 2. * zeta_h * lapse_bh2 * lapse_bh2 * msr / r_bh
224 * (1.+4.*msr) * xsr_spher_cov ;
225
226 Tenseur entmb = ent - beta_auto ;
227
228 Tenseur source = flat_scalar_prod(www - v_orb, ent.gradient())
229 + unsurc2*zeta_h*( flat_scalar_prod(v_orb, entmb.gradient())
231 / gam_euler )
232 + 2.*lapse_bh2 * msr * flat_scalar_prod(xsr_cov, v_orb)
233 * flat_scalar_prod(xsr_con, ent.gradient())
234 + unsurc2 * 2. * zeta_h
235 * (lapse_bh2*msr*(ldldpsi0
236 - flat_scalar_prod(xsr_cov, v_orb)
237 * (flat_scalar_prod(xsr_con, entmb.gradient())
238 - lapse_bh2 * (1.+4.*msr) / r_bh))
239 + a_car * hhh * gam_euler * lapse_bh3 * msr * (1.+3.*msr)
240 / r_bh) ;
241
242 source.annule(nzet, nzm1) ;
243
244 //----------------------------------------------------
245 // Resolution by means of Map_radial::poisson_compact
246 //----------------------------------------------------
247
248 Param par ;
249 int niter ;
250 par.add_int(mermax) ;
251 par.add_double(precis, 0) ;
252 par.add_double(relax, 1) ;
253 par.add_int_mod(niter) ;
254
255 if (psi0.get_etat() == ETATZERO) {
257 psi0.set() = 0 ;
258 }
259
260 source.set().va.ylm() ;
261
262 mp.poisson_compact(source(), zeta_h(), bb, par, psi0.set() ) ;
263
264 //-----------------------
265 // Check of the solution
266 //-----------------------
267
268 Tenseur bb_dpsi0 = flat_scalar_prod( bb, psi0.gradient_spher() ) ;
269
270 Cmp oper = zeta_h() * psi0().laplacien() + bb_dpsi0() ;
271
272 source.set().va.ylm_i() ;
273
274 erreur = diffrel(oper, source())(0) ;
275
276 cout << "Check of the resolution of the continuity equation : "
277 << endl ;
278 cout << "norme(source) : " << norme(source())(0)
279 << " diff oper/source : " << erreur << endl ;
280
281 //--------------------------
282 // Computation of grad(psi)
283 //--------------------------
284
285 // The computation is done component by component because
286 // psi0.gradient() is a covariant vector, whereas v_orb is a
287 // contravariant one.
288
290
291 for (int i=0; i<3; i++) {
292 d_psi.set(i) = (psi0.gradient())(i) + v_orb(i) ;
293 }
294
295 d_psi.set_triad( *(v_orb.get_triad()) ) ;
296
297 // C^1 continuation of d_psi outside the star
298 // (to ensure a smooth enthalpy field accross the stellar surface)
299 // ----------------------------------------------------------------
300
301 d_psi.annule(nzet, nzm1) ;
302 for (int i=0; i<3; i++) {
303 d_psi.set(i) = raccord_c1(d_psi(i), nzet) ;
304 }
305
306 assert( d_psi.get_triad() == &(mp.get_bvect_cart()) ) ;
307
309
310 return erreur ;
311
312 }
313
314 }
315 else {
316
317 if (eos.identify() == 5 || eos.identify() == 4 ||
318 eos.identify() == 3) {
319
320 cout << "Etoile_bin::velocity_pot_extr" << endl ;
321 cout << "The code has not prepared for this kind of EOS!!!"
322 << endl;
323 abort() ;
324
325 } // End of strange stars case
326 else {
327
328 int nzm1 = mp.get_mg()->get_nzone() - 1 ;
329
330 //--------------------------------
331 // Specific relativistic enthalpy ---> hhh
332 //--------------------------------
333
334 Tenseur hhh = exp(unsurc2 * ent) ; // = 1 at the Newtonian limit
335 hhh.set_std_base() ;
336
337 //--------------------------------------------
338 // Computation of W^i = - A^2 h Gamma_n B^i/N
339 //--------------------------------------------
340
341 Tenseur www = - a_car * hhh * gam_euler * bsn ;
342
344 // components on the mapping Cartesian basis
345
346 //-------------------------------------------------
347 // Constant value of W^i at the center of the star
348 //-------------------------------------------------
349
350 Tenseur v_orb(mp, 1, CON, mp.get_bvect_cart()) ;
351
352 v_orb.set_etat_qcq() ;
353 for (int i=0; i<3; i++) {
354 v_orb.set(i) = www(i)(0, 0, 0, 0) ;
355 }
356
357 v_orb.set_triad( *(www.get_triad()) ) ;
358 v_orb.set_std_base() ;
359
360 // Some auxiliary terms
361 // --------------------
362
363 const Coord& xx = mp.x ;
364 const Coord& yy = mp.y ;
365 const Coord& zz = mp.z ;
366
367 Tenseur r_bh(mp) ;
368 r_bh.set_etat_qcq() ;
369 r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
370 r_bh.set_std_base() ;
371
372 Tenseur msr(mp) ;
373 msr = ggrav * mass / r_bh ;
374 msr.set_std_base() ;
375
376 Tenseur xx_cov(mp, 1, COV, mp.get_bvect_cart()) ;
377 xx_cov.set_etat_qcq() ;
378 xx_cov.set(0) = xx + sepa ;
379 xx_cov.set(1) = yy ;
380 xx_cov.set(2) = zz ;
381 xx_cov.set_std_base() ;
382
383 // X_i in the spherical coordinate
384 const Coord& rr = mp.r ;
385 const Coord& st = mp.sint ;
386 const Coord& ct = mp.cost ;
387 const Coord& sp = mp.sinp ;
388 const Coord& cp = mp.cosp ;
389
390 Tenseur rr_spher_cov(mp, 1, COV, mp.get_bvect_spher()) ;
391 rr_spher_cov.set_etat_qcq() ;
392 rr_spher_cov.set(0) = rr + sepa*st*cp ;
393 rr_spher_cov.set(1) = sepa*ct*cp ;
394 rr_spher_cov.set(2) = - sepa*sp ;
395 rr_spher_cov.set_std_base() ;
396
397 Tenseur tmp_bh(mp) ;
398 tmp_bh = 0.5 * msr * msr / (1.-0.25*msr*msr) / r_bh / r_bh ;
399 tmp_bh.set_std_base() ;
400
401 //-------------------------------------------------
402 // Source and coefficients a,b for poisson_compact
403 // (idenpendent from psi0)
404 //-------------------------------------------------
405
406 Cmp dndh_log = eos.der_nbar_ent(ent(), nzet) ;
407
408 // In order to avoid any division by zero in the computation of
409 // zeta_h the value of dndh_log is set to 1
410 // in the external domains:
411 for (int l=nzet; l <= nzm1; l++) {
412 dndh_log.set(l) = 1 ;
413 }
414
415 double erreur ;
416
417 Tenseur zeta_h( ent() / dndh_log ) ;
418 zeta_h.set_std_base() ;
419
420 Tenseur tmp_zeta = 1 - unsurc2 * zeta_h ;
421 tmp_zeta.set_std_base() ;
422
423 Tenseur bb = tmp_zeta * ent.gradient_spher()
424 + unsurc2 * zeta_h * (beta_auto.gradient_spher()
425 + tmp_bh * rr_spher_cov) ;
426
427 Tenseur entmb = ent - beta_auto ;
428
429 Tenseur source = flat_scalar_prod(www - v_orb, ent.gradient())
430 + unsurc2*zeta_h*( flat_scalar_prod(v_orb, entmb.gradient())
431 - tmp_bh * flat_scalar_prod(v_orb, xx_cov)
433 / gam_euler ) ;
434
435 source.annule(nzet, nzm1) ;
436
437 //----------------------------------------------------
438 // Resolution by means of Map_radial::poisson_compact
439 //----------------------------------------------------
440
441 Param par ;
442 int niter ;
443 par.add_int(mermax) ;
444 par.add_double(precis, 0) ;
445 par.add_double(relax, 1) ;
446 par.add_int_mod(niter) ;
447
448 if (psi0.get_etat() == ETATZERO) {
450 psi0.set() = 0 ;
451 }
452
453 source.set().va.ylm() ;
454
455 mp.poisson_compact(source(), zeta_h(), bb, par, psi0.set() ) ;
456
457 //-----------------------
458 // Check of the solution
459 //-----------------------
460
461 Tenseur bb_dpsi0 = flat_scalar_prod( bb, psi0.gradient_spher() ) ;
462
463 Cmp oper = zeta_h() * psi0().laplacien() + bb_dpsi0() ;
464
465 source.set().va.ylm_i() ;
466
467 erreur = diffrel(oper, source())(0) ;
468
469 cout << "Check of the resolution of the continuity equation : "
470 << endl ;
471 cout << "norme(source) : " << norme(source())(0)
472 << " diff oper/source : " << erreur << endl ;
473
474 //--------------------------
475 // Computation of grad(psi)
476 //--------------------------
477
478 // The computation is done component by component because
479 // psi0.gradient() is a covariant vector, whereas v_orb is a
480 // contravariant one.
481
483
484 for (int i=0; i<3; i++) {
485 d_psi.set(i) = (psi0.gradient())(i) + v_orb(i) ;
486 }
487
488 d_psi.set_triad( *(v_orb.get_triad()) ) ;
489
490 // C^1 continuation of d_psi outside the star
491 // (to ensure a smooth enthalpy field accross the stellar surface)
492 // ----------------------------------------------------------------
493
494 d_psi.annule(nzet, nzm1) ;
495 for (int i=0; i<3; i++) {
496 d_psi.set(i) = raccord_c1(d_psi(i), nzet) ;
497 }
498
499 assert( d_psi.get_triad() == &(mp.get_bvect_cart()) ) ;
500
502
503 return erreur ;
504
505 }
506
507 }
508
509}
510}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Valeur va
The numerical value of the Cmp
Definition cmp.h:464
Tbl & set(int l)
Read/write of the value in a given domain.
Definition cmp.h:724
Active physical coordinates and mapping derivatives.
Definition coord.h:90
virtual int identify() const =0
Returns a number to identify the sub-classe of Eos the object belongs to.
Cmp der_nbar_ent(const Cmp &ent, int nzet, int l_min=0, const Param *par=0x0) const
Computes the logarithmic derivative from the log-enthalpy and extra parameters.
Definition eos.C:407
bool kerrschild
Indicator of the background metric: true for the Kerr-Shild metric, false for the conformally flat on...
double velocity_pot_extr(const double &mass, const double &sepa, int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
Tenseur psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case)
Definition etoile.h:833
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition etoile.h:950
Tenseur d_psi
Gradient of (in the irrotational case) (Cartesian components with respect to ref_triad )
Definition etoile.h:838
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur 's...
Definition etoile.h:828
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:432
const Eos & eos
Equation of state of the stellar matter.
Definition etoile.h:451
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition etoile.h:471
Map & mp
Mapping associated with the star.
Definition etoile.h:429
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition etoile.h:457
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition etoile.h:506
Tenseur a_car
Total conformal factor .
Definition etoile.h:515
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition etoile.h:442
Coord cosp
Definition map.h:724
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping,...
Definition map.h:791
Coord y
y coordinate centered on the grid
Definition map.h:727
Coord sint
Definition map.h:721
virtual void poisson_compact(const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const =0
Resolution of the elliptic equation in the case where the stellar interior is covered by a single do...
Coord r
r coordinate centered on the grid
Definition map.h:718
Coord sinp
Definition map.h:723
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition map.h:783
Coord z
z coordinate centered on the grid
Definition map.h:728
Coord x
x coordinate centered on the grid
Definition map.h:726
Coord cost
Definition map.h:722
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
int get_nzone() const
Returns the number of domains.
Definition grilles.h:448
Parameter storage.
Definition param.h:125
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition param.C:315
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition param.C:385
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition param.C:246
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:301
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition tenseur.h:704
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:824
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tenseur.C:636
void annule(int l)
Sets the Tenseur to zero in a given domain.
Definition tenseur.C:900
const Tenseur & gradient() const
Returns the gradient of *this (Cartesian coordinates)
Definition tenseur.C:1542
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1170
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition tenseur.C:1548
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition tenseur.C:668
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition tenseur.C:674
int get_etat() const
Returns the logical state.
Definition tenseur.h:707
void ylm()
Computes the coefficients of *this.
Definition valeur_ylm.C:138
void ylm_i()
Inverse of ylm()
Cmp exp(const Cmp &)
Exponential.
Definition cmp_math.C:270
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:504
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition cmp_math.C:481
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:348
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...
Lorene prototypes.
Definition app_hor.h:64
Standard units of space, time and mass.