LORENE
hole_bhns_rk_theta.C
1 /*
2  * Methods of class Hole_bhns to compute a forth-order Runge-Kutta
3  * integration to the theta direction for the solution of the Killing vectors
4  *
5  * (see file hole_bhns.h for documentation).
6  *
7  */
8 
9 /*
10  * Copyright (c) 2006-2007 Keisuke Taniguchi
11  *
12  * This file is part of LORENE.
13  *
14  * LORENE is free software; you can redistribute it and/or modify
15  * it under the terms of the GNU General Public License version 2
16  * as published by the Free Software Foundation.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 char hole_bhns_rk_theta_C[] = "$Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_rk_theta.C,v 1.4 2014/10/13 08:53:00 j_novak Exp $" ;
30 
31 /*
32  * $Id: hole_bhns_rk_theta.C,v 1.4 2014/10/13 08:53:00 j_novak Exp $
33  * $Log: hole_bhns_rk_theta.C,v $
34  * Revision 1.4 2014/10/13 08:53:00 j_novak
35  * Lorene classes and functions now belong to the namespace Lorene.
36  *
37  * Revision 1.3 2014/10/06 15:13:10 j_novak
38  * Modified #include directives to use c++ syntax.
39  *
40  * Revision 1.2 2008/07/02 20:49:05 k_taniguchi
41  * Typos removed.
42  *
43  * Revision 1.1 2008/05/15 19:11:01 k_taniguchi
44  * *** empty log message ***
45  *
46  *
47  * $Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_rk_theta.C,v 1.4 2014/10/13 08:53:00 j_novak Exp $
48  *
49  */
50 
51 // C++ headers
52 //#include <>
53 
54 // C headers
55 #include <cmath>
56 
57 // Lorene headers
58 #include "hole_bhns.h"
59 #include "unites.h"
60 #include "utilitaires.h"
61 
62  //---------------------------------------------------//
63  // Forth-order Runge-Kutta to the polar angle //
64  //---------------------------------------------------//
65 
66 namespace Lorene {
67 Tbl Hole_bhns::runge_kutta_theta(const Tbl& xi_i, const double& theta_i,
68  const double& phi,
69  const int& nrk_theta) const {
70 
71  using namespace Unites ;
72 
73  const Mg3d* mg = mp.get_mg() ;
74  int nt = mg->get_nt(1) ;
75 
76  Tbl xi_f(3) ; // xi_f(0)=xi_hat{theta}, xi_f(1)=xi_hat{phi}, xi_f(2)=L
77  xi_f.set_etat_qcq() ;
78 
79  if (kerrschild) {
80 
81  cout << "Not yet prepared!!!" << endl ;
82  abort() ;
83 
84  }
85  else { // Isotropic coordinates
86 
87  // Initial data at phi on the equator
88  double xi_t0 = xi_i(0) ; // xi_hat{theta}
89  double xi_p0 = xi_i(1) ; // xi_hat{phi}
90  double xi_l0 = xi_i(2) ; // L
91  double theta0 = theta_i ;
92 
93  double dt = - 0.5 * M_PI / double(nt-1) / double(nrk_theta) ;
94  // Compute from M_PI/2 to 0
95 
96  double rah = rad_ah() ;
97 
98  Scalar dlnconfo(mp) ;
99  dlnconfo = confo_tot.stdsdp() / confo_tot ;
100  dlnconfo.std_spectral_base() ;
101 
102  Scalar laplnconfo(mp) ;
103  laplnconfo = confo_tot.lapang() / confo_tot ;
104  laplnconfo.std_spectral_base() ;
105 
106  Scalar confo2(mp) ;
107  confo2 = confo_tot * confo_tot ;
108  confo2.std_spectral_base() ;
109 
110  double xi_t1, xi_t2, xi_t3, xi_t4, xi_tf ;
111  double xi_p1, xi_p2, xi_p3, xi_p4, xi_pf ;
112  double xi_l1, xi_l2, xi_l3, xi_l4, xi_lf ;
113  double f1, f2, f3, f4 ;
114  double g1, g2, g3, g4 ;
115  double h1, h2, h3, h4 ;
116 
117  // Forth-order Runge-Kutta
118  // (nrk_theta times steps between two collocation points)
119  // ------------------------------------------------------
120 
121  for (int i=0; i<nrk_theta; i++) {
122 
123  // First
124  f1 = -2. * xi_p0 * dlnconfo.val_point(rah, theta0, phi) ;
125  g1 = xi_l0 * rah * confo2.val_point(rah, theta0, phi)
126  + 2. * xi_t0 * dlnconfo.val_point(rah, theta0, phi) ;
127  h1 = - (1. - 2.*laplnconfo.val_point(rah, theta0, phi)) * xi_p0
128  / rah / confo2.val_point(rah, theta0, phi) ;
129 
130  xi_t1 = dt * f1 ;
131  xi_p1 = dt * g1 ;
132  xi_l1 = dt * h1 ;
133 
134  // Second
135  f2 = -2. * (xi_p0+0.5*xi_p1)
136  * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
137  g2 = (xi_l0+0.5*xi_l1) * rah
138  * confo2.val_point(rah, theta0+0.5*dt, phi)
139  + 2. * (xi_t0+0.5*xi_t1)
140  * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
141  h2 = - (1. - 2.*laplnconfo.val_point(rah, theta0+0.5*dt, phi))
142  * (xi_p0+0.5*xi_p1) / rah
143  / confo2.val_point(rah, theta0+0.5*dt, phi) ;
144 
145  xi_t2 = dt * f2 ;
146  xi_p2 = dt * g2 ;
147  xi_l2 = dt * h2 ;
148 
149  // Third
150  f3 = -2. * (xi_p0+0.5*xi_p2)
151  * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
152  g3 = (xi_l0+0.5*xi_l2) * rah
153  * confo2.val_point(rah, theta0+0.5*dt, phi)
154  + 2. * (xi_t0+0.5*xi_t2)
155  * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
156  h3 = - (1. - 2.*laplnconfo.val_point(rah, theta0+0.5*dt, phi))
157  * (xi_p0+0.5*xi_p2) / rah
158  / confo2.val_point(rah, theta0+0.5*dt, phi) ;
159 
160  xi_t3 = dt * f3 ;
161  xi_p3 = dt * g3 ;
162  xi_l3 = dt * h3 ;
163 
164  // Forth
165  f4 = -2. * (xi_p0+xi_p3) * dlnconfo.val_point(rah, theta0+dt, phi) ;
166  g4 = (xi_l0+xi_l3) * rah * confo2.val_point(rah, theta0+dt, phi)
167  + 2. * (xi_t0+xi_t3) * dlnconfo.val_point(rah, theta0+dt, phi) ;
168  h4 = - (1. - 2.*laplnconfo.val_point(rah, theta0+dt, phi))
169  * (xi_p0+xi_p3) / rah / confo2.val_point(rah, theta0+dt, phi) ;
170 
171  xi_t4 = dt * f4 ;
172  xi_p4 = dt * g4 ;
173  xi_l4 = dt * h4 ;
174 
175  // Final results
176  // -------------
177  xi_tf = xi_t0 + (xi_t1 + 2.*xi_t2 + 2.*xi_t3 + xi_t4) / 6. ;
178  xi_pf = xi_p0 + (xi_p1 + 2.*xi_p2 + 2.*xi_p3 + xi_p4) / 6. ;
179  xi_lf = xi_l0 + (xi_l1 + 2.*xi_l2 + 2.*xi_l3 + xi_l4) / 6. ;
180 
181  // Final results are put into the initial data
182  // in order for the next step
183  // -------------------------------------------
184  xi_t0 = xi_tf ;
185  xi_p0 = xi_pf ;
186  xi_l0 = xi_lf ;
187 
188  } // End of the loop
189 
190  xi_f.set(0) = xi_tf ;
191  xi_f.set(1) = xi_pf ;
192  xi_f.set(2) = xi_lf ;
193 
194  }
195 
196  return xi_f ;
197 
198 }
199 }
Lorene::Tbl::set
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:281
Lorene::Scalar::val_point
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
Lorene::Mg3d
Multi-domain grid.
Definition: grilles.h:273
Lorene::Mg3d::get_nt
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:457
Lorene
Lorene prototypes.
Definition: app_hor.h:64
Lorene::Hole_bhns::runge_kutta_theta
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 ...
Definition: hole_bhns_rk_theta.C:67
Lorene::Map::get_mg
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:765
Lorene::Black_hole::mp
Map & mp
Mapping associated with the black hole.
Definition: blackhole.h:80
Lorene::Scalar
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:387
Unites
Standard units of space, time and mass.
Lorene::Tbl::set_etat_qcq
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:361
Lorene::Scalar::lapang
const Scalar & lapang() const
Returns the angular Laplacian of *this , where .
Definition: scalar_deriv.C:461
Lorene::Tbl
Basic array class.
Definition: tbl.h:161
Lorene::Scalar::std_spectral_base
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
Lorene::Hole_bhns::confo_tot
Scalar confo_tot
Total conformal factor.
Definition: hole_bhns.h:169
Lorene::Black_hole::rad_ah
virtual double rad_ah() const
Radius of the apparent horizon.
Definition: blackhole_global.C:505
Lorene::Scalar::stdsdp
const Scalar & stdsdp() const
Returns of *this .
Definition: scalar_deriv.C:238
Lorene::Black_hole::kerrschild
bool kerrschild
true for a Kerr-Schild background, false for a conformally flat background
Definition: blackhole.h:85