hole_bhns_rk_theta.C

00001 /*
00002  *  Methods of class Hole_bhns to compute a forth-order Runge-Kutta
00003  *  integration to the theta direction for the solution of the Killing vectors
00004  *
00005  *    (see file hole_bhns.h for documentation).
00006  *
00007  */
00008 
00009 /*
00010  *   Copyright (c) 2006-2007 Keisuke Taniguchi
00011  *
00012  *   This file is part of LORENE.
00013  *
00014  *   LORENE is free software; you can redistribute it and/or modify
00015  *   it under the terms of the GNU General Public License version 2
00016  *   as published by the Free Software Foundation.
00017  *
00018  *   LORENE is distributed in the hope that it will be useful,
00019  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00020  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00021  *   GNU General Public License for more details.
00022  *
00023  *   You should have received a copy of the GNU General Public License
00024  *   along with LORENE; if not, write to the Free Software
00025  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00026  *
00027  */
00028 
00029 char hole_bhns_rk_theta_C[] = "$Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_rk_theta.C,v 1.2 2008/07/02 20:49:05 k_taniguchi Exp $" ;
00030 
00031 /*
00032  * $Id: hole_bhns_rk_theta.C,v 1.2 2008/07/02 20:49:05 k_taniguchi Exp $
00033  * $Log: hole_bhns_rk_theta.C,v $
00034  * Revision 1.2  2008/07/02 20:49:05  k_taniguchi
00035  * Typos removed.
00036  *
00037  * Revision 1.1  2008/05/15 19:11:01  k_taniguchi
00038  * *** empty log message ***
00039  *
00040  *
00041  * $Header: /cvsroot/Lorene/C++/Source/Hole_bhns/hole_bhns_rk_theta.C,v 1.2 2008/07/02 20:49:05 k_taniguchi Exp $
00042  *
00043  */
00044 
00045 // C++ headers
00046 //#include <>
00047 
00048 // C headers
00049 #include <math.h>
00050 
00051 // Lorene headers
00052 #include "hole_bhns.h"
00053 #include "unites.h"
00054 #include "utilitaires.h"
00055 
00056           //---------------------------------------------------//
00057           //     Forth-order Runge-Kutta to the polar angle    //
00058           //---------------------------------------------------//
00059 
00060 Tbl Hole_bhns::runge_kutta_theta(const Tbl& xi_i, const double& theta_i,
00061                  const double& phi,
00062                  const int& nrk_theta) const {
00063 
00064     using namespace Unites ;
00065 
00066     const Mg3d* mg = mp.get_mg() ;
00067     int nt = mg->get_nt(1) ;
00068 
00069     Tbl xi_f(3) ;  // xi_f(0)=xi_hat{theta}, xi_f(1)=xi_hat{phi}, xi_f(2)=L
00070     xi_f.set_etat_qcq() ;
00071 
00072     if (kerrschild) {
00073 
00074       cout << "Not yet prepared!!!" << endl ;
00075       abort() ;
00076 
00077     }
00078     else {  // Isotropic coordinates
00079 
00080       // Initial data at phi on the equator
00081       double xi_t0 = xi_i(0) ;  // xi_hat{theta}
00082       double xi_p0 = xi_i(1) ;  // xi_hat{phi}
00083       double xi_l0 = xi_i(2) ;  // L
00084       double theta0 = theta_i ;
00085 
00086       double dt = - 0.5 * M_PI / double(nt-1) / double(nrk_theta) ;
00087       // Compute from M_PI/2 to 0
00088 
00089       double rah = rad_ah() ;
00090 
00091       Scalar dlnconfo(mp) ;
00092       dlnconfo = confo_tot.stdsdp() / confo_tot ;
00093       dlnconfo.std_spectral_base() ;
00094 
00095       Scalar laplnconfo(mp) ;
00096       laplnconfo = confo_tot.lapang() / confo_tot ;
00097       laplnconfo.std_spectral_base() ;
00098 
00099       Scalar confo2(mp) ;
00100       confo2 = confo_tot * confo_tot ;
00101       confo2.std_spectral_base() ;
00102 
00103       double xi_t1, xi_t2, xi_t3, xi_t4, xi_tf ;
00104       double xi_p1, xi_p2, xi_p3, xi_p4, xi_pf ;
00105       double xi_l1, xi_l2, xi_l3, xi_l4, xi_lf ;
00106       double f1, f2, f3, f4 ;
00107       double g1, g2, g3, g4 ;
00108       double h1, h2, h3, h4 ;
00109 
00110       // Forth-order Runge-Kutta
00111       // (nrk_theta times steps between two collocation points)
00112       // ------------------------------------------------------
00113 
00114       for (int i=0; i<nrk_theta; i++) {
00115 
00116     // First
00117     f1 = -2. * xi_p0 * dlnconfo.val_point(rah, theta0, phi) ;
00118     g1 = xi_l0 * rah * confo2.val_point(rah, theta0, phi)
00119       + 2. * xi_t0 * dlnconfo.val_point(rah, theta0, phi) ;
00120     h1 = - (1. - 2.*laplnconfo.val_point(rah, theta0, phi)) * xi_p0
00121       / rah / confo2.val_point(rah, theta0, phi) ;
00122 
00123     xi_t1 = dt * f1 ;
00124     xi_p1 = dt * g1 ;
00125     xi_l1 = dt * h1 ;
00126 
00127     // Second
00128     f2 = -2. * (xi_p0+0.5*xi_p1)
00129       * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
00130     g2 = (xi_l0+0.5*xi_l1) * rah
00131       * confo2.val_point(rah, theta0+0.5*dt, phi)
00132       + 2. * (xi_t0+0.5*xi_t1)
00133       * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
00134     h2 = - (1. - 2.*laplnconfo.val_point(rah, theta0+0.5*dt, phi))
00135       * (xi_p0+0.5*xi_p1) / rah
00136       / confo2.val_point(rah, theta0+0.5*dt, phi) ;
00137 
00138     xi_t2 = dt * f2 ;
00139     xi_p2 = dt * g2 ;
00140     xi_l2 = dt * h2 ;
00141 
00142     // Third
00143     f3 = -2. * (xi_p0+0.5*xi_p2)
00144       * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
00145     g3 = (xi_l0+0.5*xi_l2) * rah
00146       * confo2.val_point(rah, theta0+0.5*dt, phi)
00147       + 2. * (xi_t0+0.5*xi_t2)
00148       * dlnconfo.val_point(rah, theta0+0.5*dt, phi) ;
00149     h3 = - (1. - 2.*laplnconfo.val_point(rah, theta0+0.5*dt, phi))
00150       * (xi_p0+0.5*xi_p2) / rah
00151       / confo2.val_point(rah, theta0+0.5*dt, phi) ;
00152 
00153     xi_t3 = dt * f3 ;
00154     xi_p3 = dt * g3 ;
00155     xi_l3 = dt * h3 ;
00156 
00157     // Forth
00158     f4 = -2. * (xi_p0+xi_p3) * dlnconfo.val_point(rah, theta0+dt, phi) ;
00159     g4 = (xi_l0+xi_l3) * rah * confo2.val_point(rah, theta0+dt, phi)
00160       + 2. * (xi_t0+xi_t3) * dlnconfo.val_point(rah, theta0+dt, phi) ;
00161     h4 = - (1. - 2.*laplnconfo.val_point(rah, theta0+dt, phi))
00162       * (xi_p0+xi_p3) / rah / confo2.val_point(rah, theta0+dt, phi) ;
00163 
00164     xi_t4 = dt * f4 ;
00165     xi_p4 = dt * g4 ;
00166     xi_l4 = dt * h4 ;
00167 
00168     // Final results
00169     // -------------
00170     xi_tf = xi_t0 + (xi_t1 + 2.*xi_t2 + 2.*xi_t3 + xi_t4) / 6. ;
00171     xi_pf = xi_p0 + (xi_p1 + 2.*xi_p2 + 2.*xi_p3 + xi_p4) / 6. ;
00172     xi_lf = xi_l0 + (xi_l1 + 2.*xi_l2 + 2.*xi_l3 + xi_l4) / 6. ;
00173 
00174     // Final results are put into the initial data
00175     // in order for the next step
00176     // -------------------------------------------
00177     xi_t0 = xi_tf ;
00178     xi_p0 = xi_pf ;
00179     xi_l0 = xi_lf ;
00180 
00181       } // End of the loop
00182 
00183       xi_f.set(0) = xi_tf ;
00184       xi_f.set(1) = xi_pf ;
00185       xi_f.set(2) = xi_lf ;
00186 
00187     }
00188 
00189     return xi_f ;
00190 
00191 }

Generated on Tue Feb 7 01:35:17 2012 for LORENE by  doxygen 1.4.6