LORENE
blackhole_killing.C
1/*
2 * Methods of class Black_hole to compute Killing vectors
3 *
4 * (see file blackhole.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 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
28char blackhole_killing_C[] = "$Header: /cvsroot/Lorene/C++/Source/Black_hole/blackhole_killing.C,v 1.4 2014/10/13 08:52:46 j_novak Exp $" ;
29
30/*
31 * $Id: blackhole_killing.C,v 1.4 2014/10/13 08:52:46 j_novak Exp $
32 * $Log: blackhole_killing.C,v $
33 * Revision 1.4 2014/10/13 08:52:46 j_novak
34 * Lorene classes and functions now belong to the namespace Lorene.
35 *
36 * Revision 1.3 2014/10/06 15:13:02 j_novak
37 * Modified #include directives to use c++ syntax.
38 *
39 * Revision 1.2 2008/07/02 20:45:07 k_taniguchi
40 * A bug removed.
41 *
42 * Revision 1.1 2008/05/15 19:33:12 k_taniguchi
43 * *** empty log message ***
44 *
45 *
46 * $Header: /cvsroot/Lorene/C++/Source/Black_hole/blackhole_killing.C,v 1.4 2014/10/13 08:52:46 j_novak Exp $
47 *
48 */
49
50// C++ headers
51//#include <>
52
53// C headers
54#include <cmath>
55
56// Lorene headers
57#include "blackhole.h"
58#include "unites.h"
59#include "utilitaires.h"
60
61 //---------------------------------------------//
62 // Killing vectors on the AH //
63 //---------------------------------------------//
64
65namespace Lorene {
66Vector Black_hole::killing_vect_bh(const Tbl& xi_i, const double& phi_i,
67 const double& theta_i, const int& nrk_phi,
68 const int& nrk_theta) const {
69
70 using namespace Unites ;
71
72 assert(xi_i.get_ndim() == 1) ;
73 assert(xi_i.get_dim(0) == 3) ;
74
75 const Mg3d* mg = mp.get_mg() ;
76 int nr = mg->get_nr(1) ;
77 int nt = mg->get_nt(1) ;
78 int np = mg->get_np(1) ;
79
80 // Vector which is returned to the main code
81 // Spherical basis, covariant
82 Vector killing(mp, COV, mp.get_bvect_spher()) ;
83
84 if (kerrschild) {
85
86 cout << "Not yet prepared!!!" << endl ;
87 abort() ;
88
89 }
90 else { // Isotropic coordinates
91
92 // Solution of the Killing vector on the equator
93 // ---------------------------------------------
94
95 double dp = 2. * M_PI / double(np) ; // Angular step
96
97 // Killing vector on the equator
98 // np+1 is for the check of xi(phi=0)= xi(phi=2pi)
99 Tbl xi_t(np+1) ;
100 xi_t.set_etat_qcq() ;
101 Tbl xi_p(np+1) ;
102 xi_p.set_etat_qcq() ;
103 Tbl xi_l(np+1) ;
104 xi_l.set_etat_qcq() ;
105
106 xi_t.set(0) = xi_i(0) ;
107 xi_p.set(0) = xi_i(1) ;
108 xi_l.set(0) = xi_i(2) ;
109
110 Tbl xi(3) ;
111 xi.set_etat_qcq() ;
112
113 Tbl xi_ini(3) ;
114 xi_ini.set_etat_qcq() ;
115 xi_ini.set(0) = xi_i(0) ;
116 xi_ini.set(1) = xi_i(1) ;
117 xi_ini.set(2) = xi_i(2) ;
118
119 double pp_0 = phi_i ; // azimuthal angle phi
120
121 for (int i=1; i<np+1; i++) {
122
123 xi = runge_kutta_phi_bh(xi_ini, pp_0, nrk_phi) ;
124
125 xi_t.set(i) = xi(0) ;
126 xi_p.set(i) = xi(1) ;
127 xi_l.set(i) = xi(2) ;
128
129 // New data for the next step
130 // -------------------------
131 pp_0 += dp ; // New angle
132 xi_ini = xi ;
133
134 }
135
136 /*
137 for (int i=0; i<np+1; i++) {
138
139 cout << "xi_t xi_p xi_l" << endl ;
140 cout << xi_t(i) << " " << xi_p(i) << " " << xi_l(i) << endl ;
141
142 }
143 arrete() ;
144 */
145
146 // Normalization of the Killing vector
147 // -----------------------------------
148
149 // Initialization of the Killing vector to the phi direction
150 Scalar xi_phi(mp) ;
151 xi_phi = 0.5 ; // If we use "1." for the initialization value,
152 // the state flag becomes "ETATUN" which does not
153 // work for "set_grid_point".
154
155 for (int k=0; k<np; k++) {
156 xi_phi.set_grid_point(0, k, nt-1, nr-1) = xi_p(k) ;
157 xi_phi.set_grid_point(1, k, nt-1, 0) = xi_p(k) ;
158 }
159 xi_phi.std_spectral_base() ;
160 /*
161 for (int l=0; l<2; l++) {
162 for (int k=0; k<np; k++) {
163 for (int j=0; j<nt; j++) {
164 for (int i=0; i<nr; i++) {
165 cout << "(l,k,j,i)=" << l << "," << k << "," << j
166 << "," << i << ": "
167 << xi_phi.val_grid_point(l,k,j,i) << endl ;
168 }
169 arrete() ;
170 }
171 arrete() ;
172 }
173 arrete() ;
174 }
175 */
176
177 // Normalization of the Killing vector
178 Scalar rr(mp) ;
179 rr = mp.r ;
180 rr.std_spectral_base() ;
181
182 Scalar st(mp) ;
183 st = mp.sint ;
184 st.std_spectral_base() ;
185
186 Scalar source_phi(mp) ;
187 source_phi = pow(confo, 2.) * rr * st / xi_phi ;
188 source_phi.std_spectral_base() ;
189
190 double rah = rad_ah() ;
191
192 int nn = 1000 ;
193 double hh = 2. * M_PI / double(nn) ;
194 double integ = 0. ;
195
196 int mm ;
197 double t1, t2, t3, t4, t5 ;
198
199 // Boole's Rule (Newton-Cotes Integral) for integration
200 // ----------------------------------------------------
201
202 assert(nn%4 == 0) ;
203 mm = nn/4 ;
204
205 for (int i=0; i<mm; i++) {
206
207 t1 = hh * double(4*i) ;
208 t2 = hh * double(4*i+1) ;
209 t3 = hh * double(4*i+2) ;
210 t4 = hh * double(4*i+3) ;
211 t5 = hh * double(4*i+4) ;
212
213 integ += (hh/45.) * (14.*source_phi.val_point(rah,M_PI/2.,t1)
214 + 64.*source_phi.val_point(rah,M_PI/2.,t2)
215 + 24.*source_phi.val_point(rah,M_PI/2.,t3)
216 + 64.*source_phi.val_point(rah,M_PI/2.,t4)
217 + 14.*source_phi.val_point(rah,M_PI/2.,t5)
218 ) ;
219
220 }
221
222 cout << "Black_hole:: t_f = " << integ << endl ;
223 double ratio = 0.5 * integ / M_PI ;
224
225 cout << "Black_hole:: t_f / 2M_PI = " << ratio << endl ;
226
227 for (int k=0; k<np; k++) {
228 xi_p.set(k) = xi_phi.val_grid_point(1, k, nt-1, 0) * ratio ;
229 }
230
231 /*
232 for (int k=0; k<np; k++) {
233 cout << "Normalized xi_p" << "(" << k << ") :" << xi_p(k) << endl ;
234 }
235 */
236
237 // Solution of the Killing vector to the pole angle
238 // ------------------------------------------------
239
240 double dt = 0.5 * M_PI / double(nt-1) ; // Angular step
241
242 // Killing vector to the polar angle
243 Tbl xi_th(nt, np) ;
244 xi_th.set_etat_qcq() ;
245 Tbl xi_ph(nt, np) ;
246 xi_ph.set_etat_qcq() ;
247 Tbl xi_ll(nt, np) ;
248 xi_ll.set_etat_qcq() ;
249
250 // Values on the equator
251 for (int i=0; i<np; i++) {
252
253 xi_th.set(nt-1, i) = xi_t(i) ;
254 xi_ph.set(nt-1, i) = xi_p(i) ;
255 xi_ll.set(nt-1, i) = xi_l(i) ;
256
257 }
258
259 for (int i=0; i<np; i++) {
260
261 // Value at theta=pi/2, phi=phi_i
262 xi_ini.set(0) = xi_t(i) ;
263 xi_ini.set(1) = xi_p(i) ;
264 xi_ini.set(2) = xi_l(i) ;
265
266 double pp = double(i) * dp ;
267 double tt_0 = theta_i ; // polar angle theta
268
269 for (int j=1; j<nt; j++) {
270
271 xi = runge_kutta_theta_bh(xi_ini, tt_0, pp, nrk_theta) ;
272
273 xi_th.set(nt-1-j, i) = xi(0) ;
274 xi_ph.set(nt-1-j, i) = xi(1) ;
275 xi_ll.set(nt-1-j, i) = xi(2) ;
276
277 // New data for the nxt step
278 // -------------------------
279 tt_0 -= dt ; // New angle
280 xi_ini = xi ;
281
282 } // End of the loop to the polar direction
283
284 } // End of the loop to the azimuhtal direction
285
286
287 // Construction of the Killing vector
288 // ----------------------------------
289
290 // Definition of the vector is at the top of this code
291 killing.set_etat_qcq() ;
292 killing.set(1) = 0. ; // r component
293 killing.set(2) = 0.5 ; // initialization of the theta component
294 killing.set(3) = 0.5 ; // initialization of the phi component
295
296 for (int l=0; l<2; l++) {
297 for (int i=0; i<nr; i++) {
298 for (int j=0; j<nt; j++) {
299 for (int k=0; k<np; k++) {
300 (killing.set(2)).set_grid_point(l, k, j, i) = xi_th(j, k) ;
301 (killing.set(3)).set_grid_point(l, k, j, i) = xi_ph(j, k) ;
302 }
303 }
304 }
305 }
306 killing.std_spectral_base() ;
307
308 // Check the normalization
309 // -----------------------
310
311 double check_norm = 0. ;
312 source_phi = pow(confo, 2.) * rr * st / killing(3) ;
313 source_phi.std_spectral_base() ;
314
315 for (int i=0; i<mm; i++) {
316
317 t1 = hh * double(4*i) ;
318 t2 = hh * double(4*i+1) ;
319 t3 = hh * double(4*i+2) ;
320 t4 = hh * double(4*i+3) ;
321 t5 = hh * double(4*i+4) ;
322
323 check_norm += (hh/45.) *
324 ( 14.*source_phi.val_point(rah,M_PI/4.,t1)
325 + 64.*source_phi.val_point(rah,M_PI/4.,t2)
326 + 24.*source_phi.val_point(rah,M_PI/4.,t3)
327 + 64.*source_phi.val_point(rah,M_PI/4.,t4)
328 + 14.*source_phi.val_point(rah,M_PI/4.,t5) ) ;
329
330 }
331
332 cout << "Black_hole:: t_f for M_PI/4 = " << check_norm / M_PI
333 << " M_PI" << endl ;
334
335 } // End of the loop for isotropic coordinates
336
337 return killing ;
338
339}
340}
Tbl runge_kutta_phi_bh(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...
Map & mp
Mapping associated with the black hole.
Definition blackhole.h:80
virtual double rad_ah() const
Radius of the apparent horizon.
Tbl runge_kutta_theta_bh(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 ...
Vector killing_vect_bh(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.
bool kerrschild
true for a Kerr-Schild background, false for a conformally flat background
Definition blackhole.h:85
Scalar confo
Conformal factor generated by the black hole.
Definition blackhole.h:118
Coord sint
Definition map.h:721
Coord r
r coordinate centered on the grid
Definition map.h:718
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition map.h:783
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
Multi-domain grid.
Definition grilles.h:273
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:462
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:457
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:452
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:387
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
Definition scalar.C:784
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:637
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:684
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:890
Basic array class.
Definition tbl.h:161
int get_ndim() const
Gives the number of dimensions (ie dim.ndim)
Definition tbl.h:400
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tbl.C:361
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition tbl.h:281
int get_dim(int i) const
Gives the i-th dimension (ie dim.dim[i])
Definition tbl.h:403
Tensor field of valence 1.
Definition vector.h:188
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition vector.C:316
Scalar & set(int)
Read/write access to a component.
Definition vector.C:296
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:348
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition tensor.C:481
Lorene prototypes.
Definition app_hor.h:64
Standard units of space, time and mass.