LORENE
hole_bhns_killing.C
1 /*
2  * Methods of class Hole_bhns to compute Killing vectors
3  *
4  * (see file hole_bhns.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2006-2007 Keisuke Taniguchi
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License version 2
15  * as published by the Free Software Foundation.
16  *
17  * LORENE is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with LORENE; if not, write to the Free Software
24  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25  *
26  */
27 
28 
29 
30 /*
31  * $Id: hole_bhns_killing.C,v 1.5 2016/12/05 16:17:55 j_novak Exp $
32  * $Log: hole_bhns_killing.C,v $
33  * Revision 1.5 2016/12/05 16:17:55 j_novak
34  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35  *
36  * Revision 1.4 2014/10/13 08:53:00 j_novak
37  * Lorene classes and functions now belong to the namespace Lorene.
38  *
39  * Revision 1.3 2014/10/06 15:13:10 j_novak
40  * Modified #include directives to use c++ syntax.
41  *
42  * Revision 1.2 2008/07/02 21:01:11 k_taniguchi
43  * A bug removed.
44  *
45  * Revision 1.1 2008/05/15 19:09:53 k_taniguchi
46  * *** empty log message ***
47  *
48  *
49  * $Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_killing.C,v 1.5 2016/12/05 16:17:55 j_novak Exp $
50  *
51  */
52 
53 // C++ headers
54 //#include <>
55 
56 // C headers
57 #include <cmath>
58 
59 // Lorene headers
60 #include "hole_bhns.h"
61 #include "unites.h"
62 #include "utilitaires.h"
63 
64  //---------------------------------------------//
65  // Killing vectors on the AH //
66  //---------------------------------------------//
67 
68 namespace Lorene {
69 Vector Hole_bhns::killing_vect(const Tbl& xi_i, const double& phi_i,
70  const double& theta_i, const int& nrk_phi,
71  const int& nrk_theta) const {
72 
73  using namespace Unites ;
74 
75  assert(xi_i.get_ndim() == 1) ;
76  assert(xi_i.get_dim(0) == 3) ;
77 
78  const Mg3d* mg = mp.get_mg() ;
79  int nr = mg->get_nr(1) ;
80  int nt = mg->get_nt(1) ;
81  int np = mg->get_np(1) ;
82 
83  // Vector which is returned to the main code
84  // Spherical basis, covariant
85  Vector killing(mp, COV, mp.get_bvect_spher()) ;
86 
87  if (kerrschild) {
88 
89  cout << "Not yet prepared!!!" << endl ;
90  abort() ;
91 
92  }
93  else { // Isotropic coordinates
94 
95  // Solution of the Killing vector on the equator
96  // ---------------------------------------------
97 
98  double dp = 2. * M_PI / double(np) ; // Angular step
99 
100  // Killing vector on the equator
101  // np+1 is for the check of xi(phi=0)= xi(phi=2pi)
102  Tbl xi_t(np+1) ;
103  xi_t.set_etat_qcq() ;
104  Tbl xi_p(np+1) ;
105  xi_p.set_etat_qcq() ;
106  Tbl xi_l(np+1) ;
107  xi_l.set_etat_qcq() ;
108 
109  xi_t.set(0) = xi_i(0) ;
110  xi_p.set(0) = xi_i(1) ;
111  xi_l.set(0) = xi_i(2) ;
112 
113  Tbl xi(3) ;
114  xi.set_etat_qcq() ;
115 
116  Tbl xi_ini(3) ;
117  xi_ini.set_etat_qcq() ;
118  xi_ini.set(0) = xi_i(0) ;
119  xi_ini.set(1) = xi_i(1) ;
120  xi_ini.set(2) = xi_i(2) ;
121 
122  double pp_0 = phi_i ; // azimuthal angle phi
123 
124  for (int i=1; i<np+1; i++) {
125 
126  xi = runge_kutta_phi(xi_ini, pp_0, nrk_phi) ;
127 
128  xi_t.set(i) = xi(0) ;
129  xi_p.set(i) = xi(1) ;
130  xi_l.set(i) = xi(2) ;
131 
132  // New data for the next step
133  // -------------------------
134  pp_0 += dp ; // New angle
135  xi_ini = xi ;
136 
137  }
138 
139  cout << "Hole_bhns::killing_vect :" << endl ;
140  cout.precision(16) ;
141  cout << " xi_p(0) : " << xi_p(0) << endl ;
142  cout << " xi_p(np) : " << xi_p(np) << endl ;
143 
144  // Normalization of the Killing vector
145  // -----------------------------------
146 
147  // Initialization of the Killing vector to the phi direction
148  Scalar xi_phi(mp) ;
149  xi_phi = 0.5 ; // If we use "1." for the initialization value,
150  // the state flag becomes "ETATUN" which does not
151  // work for "set_grid_point".
152 
153  for (int k=0; k<np; k++) {
154  xi_phi.set_grid_point(0, k, nt-1, nr-1) = xi_p(k) ;
155  xi_phi.set_grid_point(1, k, nt-1, 0) = xi_p(k) ;
156  }
157  xi_phi.std_spectral_base() ;
158 
159  // Normalization of the Killing vector
160  Scalar rr(mp) ;
161  rr = mp.r ;
162  rr.std_spectral_base() ;
163 
164  Scalar st(mp) ;
165  st = mp.sint ;
166  st.std_spectral_base() ;
167 
168  Scalar source_phi(mp) ;
169  source_phi = pow(confo_tot, 2.) * rr * st / xi_phi ;
170  source_phi.std_spectral_base() ;
171 
172  double rah = rad_ah() ;
173 
174  int nn = 1000 ;
175  double hh = 2. * M_PI / double(nn) ;
176  double integ = 0. ;
177 
178  int mm ;
179  double t1, t2, t3, t4, t5 ;
180 
181  // Boole's Rule (Newton-Cotes Integral) for integration
182  // ----------------------------------------------------
183 
184  assert(nn%4 == 0) ;
185  mm = nn/4 ;
186 
187  for (int i=0; i<mm; i++) {
188 
189  t1 = hh * double(4*i) ;
190  t2 = hh * double(4*i+1) ;
191  t3 = hh * double(4*i+2) ;
192  t4 = hh * double(4*i+3) ;
193  t5 = hh * double(4*i+4) ;
194 
195  integ += (hh/45.) * (14.*source_phi.val_point(rah,M_PI/2.,t1)
196  + 64.*source_phi.val_point(rah,M_PI/2.,t2)
197  + 24.*source_phi.val_point(rah,M_PI/2.,t3)
198  + 64.*source_phi.val_point(rah,M_PI/2.,t4)
199  + 14.*source_phi.val_point(rah,M_PI/2.,t5)
200  ) ;
201 
202  }
203 
204  cout << "Hole_bhns:: t_f = " << integ << endl ;
205  double ratio = 0.5 * integ / M_PI ;
206 
207  cout << "Hole_bhns:: t_f / 2M_PI = " << ratio << endl ;
208 
209  for (int k=0; k<np; k++) {
210  xi_p.set(k) = xi_phi.val_grid_point(1, k, nt-1, 0) * ratio ;
211  }
212 
213 
214  // Solution of the Killing vector to the pole angle
215  // ------------------------------------------------
216 
217  double dt = 0.5 * M_PI / double(nt-1) ; // Angular step
218 
219  // Killing vector to the polar angle
220  Tbl xi_th(nt, np) ;
221  xi_th.set_etat_qcq() ;
222  Tbl xi_ph(nt, np) ;
223  xi_ph.set_etat_qcq() ;
224  Tbl xi_ll(nt, np) ;
225  xi_ll.set_etat_qcq() ;
226 
227  // Values on the equator
228  for (int i=0; i<np; i++) {
229 
230  xi_th.set(nt-1, i) = xi_t(i) ;
231  xi_ph.set(nt-1, i) = xi_p(i) ;
232  xi_ll.set(nt-1, i) = xi_l(i) ;
233 
234  }
235 
236  for (int i=0; i<np; i++) {
237 
238  // Value at theta=pi/2, phi=phi_i
239  xi_ini.set(0) = xi_t(i) ;
240  xi_ini.set(1) = xi_p(i) ;
241  xi_ini.set(2) = xi_l(i) ;
242 
243  double pp = double(i) * dp ;
244  double tt_0 = theta_i ; // polar angle theta
245 
246  for (int j=1; j<nt; j++) {
247 
248  xi = runge_kutta_theta(xi_ini, tt_0, pp, nrk_theta) ;
249 
250  xi_th.set(nt-1-j, i) = xi(0) ;
251  xi_ph.set(nt-1-j, i) = xi(1) ;
252  xi_ll.set(nt-1-j, i) = xi(2) ;
253 
254  // New data for the nxt step
255  // -------------------------
256  tt_0 -= dt ; // New angle
257  xi_ini = xi ;
258 
259  } // End of the loop to the polar direction
260 
261  } // End of the loop to the azimuhtal direction
262 
263  cout << "Hole_bhns::killing_vect :" << endl ;
264  cout.precision(16) ;
265  cout << " xi_ph(nt-1,0) : " << xi_ph(nt-1,0) << endl ;
266  cout << " xi_ph(0,0) : " << xi_ph(0,0) << endl ;
267  cout << " xi_ph(0,np/4) : " << xi_ph(0,np/4) << endl ;
268  cout << " xi_ph(0,np/2) : " << xi_ph(0,np/2) << endl ;
269  cout << " xi_ph(0,3np/4) : " << xi_ph(0,3*np/4) << endl ;
270 
271 
272  // Construction of the Killing vector
273  // ----------------------------------
274 
275  // Definition of the vector is at the top of this code
276  killing.set_etat_qcq() ;
277  killing.set(1) = 0.5 ; // initialization of L instead of
278  // the r component
279  killing.set(2) = 0.5 ; // initialization of the theta component
280  killing.set(3) = 0.5 ; // initialization of the phi component
281 
282  for (int l=0; l<2; l++) {
283  for (int i=0; i<nr; i++) {
284  for (int j=0; j<nt; j++) {
285  for (int k=0; k<np; k++) {
286  (killing.set(1)).set_grid_point(l, k, j, i) = xi_ll(j, k) ;
287  (killing.set(2)).set_grid_point(l, k, j, i) = xi_th(j, k) ;
288  (killing.set(3)).set_grid_point(l, k, j, i) = xi_ph(j, k) ;
289  }
290  }
291  }
292  }
293  killing.std_spectral_base() ;
294 
295  // Check the normalization
296  // -----------------------
297 
298  double check_norm1 = 0. ;
299  double check_norm2 = 0. ;
300  source_phi = pow(confo_tot, 2.) * rr * st / killing(3) ;
301  source_phi.std_spectral_base() ;
302 
303  for (int i=0; i<mm; i++) {
304 
305  t1 = hh * double(4*i) ;
306  t2 = hh * double(4*i+1) ;
307  t3 = hh * double(4*i+2) ;
308  t4 = hh * double(4*i+3) ;
309  t5 = hh * double(4*i+4) ;
310 
311  check_norm1 += (hh/45.) *
312  ( 14.*source_phi.val_point(rah,M_PI/4.,t1)
313  + 64.*source_phi.val_point(rah,M_PI/4.,t2)
314  + 24.*source_phi.val_point(rah,M_PI/4.,t3)
315  + 64.*source_phi.val_point(rah,M_PI/4.,t4)
316  + 14.*source_phi.val_point(rah,M_PI/4.,t5) ) ;
317 
318  check_norm2 += (hh/45.) *
319  ( 14.*source_phi.val_point(rah,M_PI/8.,t1)
320  + 64.*source_phi.val_point(rah,M_PI/8.,t2)
321  + 24.*source_phi.val_point(rah,M_PI/8.,t3)
322  + 64.*source_phi.val_point(rah,M_PI/8.,t4)
323  + 14.*source_phi.val_point(rah,M_PI/8.,t5) ) ;
324 
325  }
326 
327  cout.precision(16) ;
328  cout << "Hole_bhns:: t_f for M_PI/4 = " << check_norm1 / M_PI
329  << " M_PI" << endl ;
330  cout << "Hole_bhns:: t_f for M_PI/8 = " << check_norm2 / M_PI
331  << " M_PI" << endl ;
332 
333  // Checking the accuracy of the Killing vector
334  // -------------------------------------------
335 
336  // xi^i D_i L = 0
337  Scalar dldt(mp) ;
338  dldt = killing(1).dsdt() ;
339 
340  Scalar dldp(mp) ;
341  dldp = killing(1).stdsdp() ;
342 
343  Scalar xidl(mp) ;
344  xidl = killing(2) * dldt + killing(3) * dldp ;
345  xidl.std_spectral_base() ;
346 
347  double xidl_error1 = 0. ;
348  double xidl_error2 = 0. ;
349  double xidl_norm = 0. ;
350 
351  for (int j=0; j<nt; j++) {
352  for (int k=0; k<np/2; k++) {
353  xidl_error1 += xidl.val_grid_point(1, k, j, 0) ;
354  }
355  }
356 
357  for (int j=0; j<nt; j++) {
358  for (int k=np/2; k<np; k++) {
359  xidl_error2 += xidl.val_grid_point(1, k, j, 0) ;
360  }
361  }
362 
363  for (int j=0; j<nt; j++) {
364  for (int k=0; k<np; k++) {
365  xidl_norm += fabs(xidl.val_grid_point(1, k, j, 0)) ;
366  }
367  }
368 
369  cout.precision(6) ;
370  cout << "Relative error in the 1st condition : "
371  << xidl_error1 / xidl_norm / double(nt) / double(np) * 2.
372  << " "
373  << xidl_error2 / xidl_norm / double(nt) / double(np) * 2.
374  << " "
375  << (xidl_error1+xidl_error2)/xidl_norm/double(nt)/double(np)
376  << endl ;
377 
378  // D_i xi^i = 0
379  Scalar xitst(mp) ;
380  xitst = pow(confo_tot, 2.) * rr * st * killing(2) ;
381  xitst.std_spectral_base() ;
382 
383  Scalar dxitstdt(mp) ;
384  dxitstdt = xitst.dsdt() ;
385 
386  Scalar xip(mp) ;
387  xip = pow(confo_tot, 2.) * rr * killing(3) ;
388  xip.std_spectral_base() ;
389 
390  Scalar dxipdp(mp) ;
391  dxipdp = xip.stdsdp() ;
392 
393  Scalar dxi(mp) ;
394  dxi = dxitstdt + st * dxipdp ;
395  dxi.std_spectral_base() ;
396 
397  double dxi_error1 = 0. ;
398  double dxi_error2 = 0. ;
399  double dxi_norm = 0. ;
400 
401  for (int j=0; j<nt; j++) {
402  for (int k=0; k<np/2; k++) {
403  dxi_error1 += dxi.val_grid_point(1, k, j, 0) ;
404  }
405  }
406 
407  for (int j=0; j<nt; j++) {
408  for (int k=np/2; k<np; k++) {
409  dxi_error2 += dxi.val_grid_point(1, k, j, 0) ;
410  }
411  }
412 
413  for (int j=0; j<nt; j++) {
414  for (int k=0; k<np; k++) {
415  dxi_norm += fabs(dxi.val_grid_point(1, k, j, 0)) ;
416  }
417  }
418 
419  cout << "Relative error in the 2nd condition : "
420  << dxi_error1 / dxi_norm / double(nt) / double(np) * 2.
421  << " "
422  << dxi_error2 / dxi_norm / double(nt) / double(np) * 2.
423  << " "
424  << (dxi_error1+dxi_error2)/dxi_norm/double(nt)/double(np)
425  << endl ;
426 
427  // e^{ij} D_i \xi_j = 2L
428  Scalar xipst(mp) ;
429  xipst = pow(confo_tot, 2.) * rr * st * killing(3) ;
430  xipst.std_spectral_base() ;
431 
432  Scalar dxipstdt(mp) ;
433  dxipstdt = xipst.dsdt() ;
434 
435  Scalar xit(mp) ;
436  xit = pow(confo_tot, 2.) * rr * killing(2) ;
437  xit.std_spectral_base() ;
438 
439  Scalar dxitdp(mp) ;
440  dxitdp = xit.stdsdp() ;
441 
442  Scalar dxi2l(mp) ;
443  dxi2l = dxipstdt - st * dxitdp
444  - 2. * killing(1) * pow(confo_tot, 4.) * rr * rr * st ;
445  dxi2l.std_spectral_base() ;
446 
447  double dxi2l_error1 = 0. ;
448  double dxi2l_error2 = 0. ;
449  double dxi2l_norm = 0. ;
450 
451  for (int j=0; j<nt; j++) {
452  for (int k=0; k<np/2; k++) {
453  dxi2l_error1 += dxi2l.val_grid_point(1, k, j, 0) ;
454  }
455  }
456 
457  for (int j=0; j<nt; j++) {
458  for (int k=np/2; k<np; k++) {
459  dxi2l_error2 += dxi2l.val_grid_point(1, k, j, 0) ;
460  }
461  }
462 
463  for (int j=0; j<nt; j++) {
464  for (int k=0; k<np; k++) {
465  dxi2l_norm += fabs(dxi2l.val_grid_point(1, k, j, 0)) ;
466  }
467  }
468 
469  cout << "Relative error in the 3rd condition : "
470  << dxi2l_error1 / dxi2l_norm / double(nt) / double(np) * 2.
471  << " "
472  << dxi2l_error2 / dxi2l_norm / double(nt) / double(np) * 2.
473  << " "
474  << (dxi2l_error1+dxi2l_error2)/dxi2l_norm/double(nt)/double(np)
475  << endl ;
476 
477  cout.precision(4) ;
478 
479  } // End of the loop for isotropic coordinates
480 
481  return killing ;
482 
483 }
484 }
Map & mp
Mapping associated with the black hole.
Definition: blackhole.h:80
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:795
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Scalar & dsdt() const
Returns of *this .
Definition: scalar_deriv.C:208
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:777
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
bool kerrschild
true for a Kerr-Schild background, false for a conformally flat background
Definition: blackhole.h:85
Tbl runge_kutta_theta(const Tbl &xi_i, const double &theta_i, const double &phi, const int &nrk) const
Compute a forth-order Runge-Kutta integration to the theta direction for the solution of the Killing ...
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
Tensor field of valence 1.
Definition: vector.h:188
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:643
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
Coord sint
Definition: map.h:733
int get_ndim() const
Gives the number of dimensions (ie dim.ndim)
Definition: tbl.h:420
Vector killing_vect(const Tbl &xi_i, const double &phi_i, const double &theta_i, const int &nrk_phi, const int &nrk_theta) const
Compute the Killing vector of a black hole normalized so that its affine length is 2 M_PI...
virtual double rad_ah() const
Radius of the apparent horizon.
const Scalar & stdsdp() const
Returns of *this .
Definition: scalar_deriv.C:238
int get_dim(int i) const
Gives the i-th dimension (ie dim.dim[i])
Definition: tbl.h:423
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
Definition: scalar.C:896
Tbl runge_kutta_phi(const Tbl &xi_i, const double &phi_i, const int &nrk) const
Compute a forth-order Runge-Kutta integration to the phi direction for the solution of the Killing ve...
Scalar confo_tot
Total conformal factor.
Definition: hole_bhns.h:169
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
Multi-domain grid.
Definition: grilles.h:279
double & set_grid_point(int l, int k, int j, int i)
Setting the value of the field at a given grid point.
Definition: scalar.h:690
Basic array class.
Definition: tbl.h:164
Coord r
r coordinate centered on the grid
Definition: map.h:730