blackhole_extr_curv.C

00001 /*
00002  *  Method of class Black_hole to compute the extrinsic curvature tensor
00003  *
00004  *    (see file blackhole.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2005-2007 Keisuke Taniguchi
00010  *
00011  *   This file is part of LORENE.
00012  *
00013  *   LORENE is free software; you can redistribute it and/or modify
00014  *   it under the terms of the GNU General Public License version 2
00015  *   as published by the Free Software Foundation.
00016  *
00017  *   LORENE is distributed in the hope that it will be useful,
00018  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00019  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00020  *   GNU General Public License for more details.
00021  *
00022  *   You should have received a copy of the GNU General Public License
00023  *   along with LORENE; if not, write to the Free Software
00024  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00025  *
00026  */
00027 
00028 char blackhole_extr_curv_C[] = "$Header: /cvsroot/Lorene/C++/Source/Black_hole/blackhole_extr_curv.C,v 1.2 2008/05/15 19:27:14 k_taniguchi Exp $" ;
00029 
00030 /*
00031  * $Id: blackhole_extr_curv.C,v 1.2 2008/05/15 19:27:14 k_taniguchi Exp $
00032  * $Log: blackhole_extr_curv.C,v $
00033  * Revision 1.2  2008/05/15 19:27:14  k_taniguchi
00034  * Change of some parameters.
00035  *
00036  * Revision 1.1  2007/06/22 01:19:32  k_taniguchi
00037  * *** empty log message ***
00038  *
00039  *
00040  * $Header: /cvsroot/Lorene/C++/Source/Black_hole/blackhole_extr_curv.C,v 1.2 2008/05/15 19:27:14 k_taniguchi Exp $
00041  *
00042  */
00043 
00044 // C++ headers
00045 //#include <>
00046 
00047 // C headers
00048 #include <math.h>
00049 
00050 // Lorene headers
00051 #include "blackhole.h"
00052 #include "utilitaires.h"
00053 #include "unites.h"
00054 
00055 void Black_hole::extr_curv_bh() {
00056 
00057     // Fundamental constants and units
00058     // -------------------------------
00059     using namespace Unites ;
00060 
00061     if (kerrschild) {
00062 
00063         double mass = ggrav * mass_bh ;
00064 
00065     Scalar rr(mp) ;
00066     rr = mp.r ;
00067     rr.std_spectral_base() ;
00068     Scalar st(mp) ;
00069     st = mp.sint ;
00070     st.std_spectral_base() ;
00071     Scalar ct(mp) ;
00072     ct = mp.cost ;
00073     ct.std_spectral_base() ;
00074     Scalar sp(mp) ;
00075     sp = mp.sinp ;
00076     sp.std_spectral_base() ;
00077     Scalar cp(mp) ;
00078     cp = mp.cosp ;
00079     cp.std_spectral_base() ;
00080 
00081     Vector ll(mp, CON, mp.get_bvect_cart()) ;
00082     ll.set_etat_qcq() ;
00083     ll.set(1) = st * cp ;
00084     ll.set(2) = st * sp ;
00085     ll.set(3) = ct ;
00086     ll.std_spectral_base() ;
00087 
00088 
00089     // Computation of \tilde{A}^{ij}
00090     // -----------------------------
00091 
00092     Scalar divshift(mp) ;
00093     divshift = shift_rs(1).dsdx() + shift_rs(2).dsdy()
00094       + shift_rs(3).dsdz() ;
00095     divshift.std_spectral_base() ;
00096 
00097     Sym_tensor flat_taij(mp, CON, mp.get_bvect_cart()) ;
00098     flat_taij.set_etat_qcq() ;
00099 
00100     for (int i=1; i<=3; i++) {
00101         for (int j=1; j<=3; j++) {
00102             flat_taij.set(i,j) = shift_rs(j).deriv(i)
00103           + shift_rs(i).deriv(j)
00104           - 2. * divshift * flat.con()(i,j) / 3. ;
00105         }
00106     }
00107 
00108     flat_taij.std_spectral_base() ;
00109 
00110     Scalar lapse_bh(mp) ;
00111     lapse_bh = 1. / sqrt(1. + 2. * mass / rr) ;
00112     lapse_bh.std_spectral_base() ;
00113     lapse_bh.annule_domain(0) ;
00114     lapse_bh.raccord(1) ;
00115 
00116     Sym_tensor curv_taij(mp, CON, mp.get_bvect_cart()) ;
00117     curv_taij.set_etat_qcq() ;
00118 
00119     for (int i=1; i<=3; i++) {
00120         for (int j=1; j<=3; j++) {
00121             curv_taij.set(i,j) = -2. * lapse_bh * lapse_bh * mass
00122           * (ll(i)*(shift_rs(j).dsdr()) + ll(j)*(shift_rs(i).dsdr())
00123              - 2. * ll(i) * ll(j) * divshift / 3.) / rr ;
00124         }
00125     }
00126 
00127     curv_taij.std_spectral_base() ;
00128 
00129     Sym_tensor resi_taij(mp, CON, mp.get_bvect_cart()) ;
00130     resi_taij.set_etat_qcq() ;
00131 
00132     for (int i=1; i<=3; i++) {
00133         for (int j=1; j<=3; j++) {
00134             resi_taij.set(i,j) = 2. * lapse_bh * lapse_bh * mass
00135           * ( ll(i) * shift_rs(j) + ll(j) * shift_rs(i)
00136               + ( flat.con()(i,j)
00137               - lapse_bh*lapse_bh*(9.+14.*mass/rr)*ll(i)*ll(j) )
00138               * ( ll(1)*shift_rs(1)+ll(2)*shift_rs(2)
00139               +ll(3)*shift_rs(3) )/ 3. )
00140           / rr / rr ;
00141         }
00142     }
00143 
00144     resi_taij.std_spectral_base() ;
00145     resi_taij.inc_dzpuis(2) ;
00146 
00147     taij_rs = 0.5 * pow(confo, 7.)
00148       * (flat_taij + curv_taij + resi_taij) / lapconf ;
00149 
00150     taij_rs.std_spectral_base() ;
00151     taij_rs.annule_domain(0) ;
00152 
00153     Sym_tensor taij_bh(mp, CON, mp.get_bvect_cart()) ;
00154     taij_bh.set_etat_qcq() ;
00155 
00156     for (int i=1; i<=3; i++) {
00157         for (int j=1; j<=3; j++) {
00158             taij_bh.set(i,j) = 2.*pow(lapse_bh,6.)*mass*(2.+3.*mass/rr)
00159           *( (1.+2.*mass/rr) * flat.con()(i,j)
00160              - (3.+2.*mass/rr) * ll(i) * ll(j) )
00161           *pow(confo, 7.)/lapconf/3./rr/rr ;
00162         }
00163     }
00164 
00165     taij_bh.std_spectral_base() ;
00166     taij_bh.inc_dzpuis(2) ;
00167     taij_bh.annule_domain(0) ;
00168 
00169     taij = taij_rs + taij_bh ;
00170     taij.std_spectral_base() ;
00171     taij.annule_domain(0) ;
00172 
00173     /*
00174     Sym_tensor taij_ks_con(mp, CON, mp.get_bvect_cart()) ;
00175     taij_ks_con.set_etat_qcq() ;
00176 
00177     for (int i=1; i<=3; i++) {
00178         for (int j=1; j<=3; j++) {
00179             taij_ks_con.set(i,j) = 2.*pow(lap_bh2,2.5)*mass
00180           * (2.+3.*mass/rr)
00181           * ( (1.+2.*mass/rr)*flat.con()(i,j)
00182               - (3.+2.*mass/rr)*ll(i)*ll(j) ) / 3. / rr / rr ;
00183         }
00184     }
00185     taij_ks_con.std_spectral_base() ;
00186     taij_ks_con.annule_domain(0) ;
00187     taij_ks_con.inc_dzpuis(2) ;
00188 
00189     cout << "taij(1,1) - taij_ks_con(1,1) : " << endl ;
00190     cout << taij(1,1) - taij_ks_con(1,1) << endl ;
00191     arrete() ;
00192 
00193     cout << "taij_ks_con(1,1) : " << endl ;
00194     cout << taij_ks_con(1,1) << endl ;
00195     arrete() ;
00196 
00197     cout << "taij(1,1) : " << endl ;
00198     cout << taij(1,1) << endl ;
00199     arrete() ;
00200     */
00201 
00202     // Computation of \tilde{A}^{ij} \tilde{A}_{ij}
00203     // --------------------------------------------
00204 
00205     Sym_tensor flat_dshift(mp, COV, mp.get_bvect_cart()) ;
00206     flat_dshift.set_etat_qcq() ;
00207 
00208     for (int i=1; i<=3; i++) {
00209         for (int j=1; j<=3; j++) {
00210             flat_dshift.set(i,j) = flat.cov()(j,1)*(shift_rs(1).deriv(i))
00211           + flat.cov()(j,2)*(shift_rs(2).deriv(i))
00212           + flat.cov()(j,3)*(shift_rs(3).deriv(i))
00213           + flat.cov()(i,1)*(shift_rs(1).deriv(j))
00214           + flat.cov()(i,2)*(shift_rs(2).deriv(j))
00215           + flat.cov()(i,3)*(shift_rs(3).deriv(j))
00216           - 2. * divshift * flat.cov()(i,j) / 3. ;
00217         }
00218     }
00219 
00220     flat_dshift.std_spectral_base() ;
00221 
00222     Sym_tensor curv_dshift(mp, COV, mp.get_bvect_cart()) ;
00223     curv_dshift.set_etat_qcq() ;
00224 
00225     for (int i=1; i<=3; i++) {
00226         for (int j=1; j<=3; j++) {
00227             curv_dshift.set(i,j) = 2. * mass
00228           *( ll(j) *( ll(1)*(shift_rs(1).deriv(i))
00229                   + ll(2)*(shift_rs(2).deriv(i))
00230                   + ll(3)*(shift_rs(3).deriv(i)) )
00231              + ll(i) *( ll(1)*(shift_rs(1).deriv(j))
00232                 + ll(2)*(shift_rs(2).deriv(j))
00233                 + ll(3)*(shift_rs(3).deriv(j)) )
00234              - 2. * divshift * ll(i) * ll(j) / 3. ) / rr ;
00235         }
00236     }
00237 
00238     curv_dshift.std_spectral_base() ;
00239 
00240     Sym_tensor tmp1(mp, COV, mp.get_bvect_cart()) ;
00241     tmp1.set_etat_qcq() ;
00242 
00243     for (int i=1; i<=3; i++) {
00244         for (int j=1; j<=3; j++) {
00245             tmp1.set(i,j) = 2. * mass
00246           *(ll(j)*( (flat.cov()(i,1)+2.*mass*ll(i)*ll(1)/rr)
00247                 *shift_rs(1)
00248                 + (flat.cov()(i,2)+2.*mass*ll(i)*ll(2)/rr)
00249                 *shift_rs(2)
00250                 + (flat.cov()(i,3)+2.*mass*ll(i)*ll(3)/rr)
00251                 *shift_rs(3)
00252                 )
00253             + ll(i)*( (flat.cov()(j,1)+2.*mass*ll(j)*ll(1)/rr)
00254                   *shift_rs(1)
00255                   + (flat.cov()(j,2)+2.*mass*ll(j)*ll(2)/rr)
00256                   *shift_rs(2)
00257                   + (flat.cov()(j,3)+2.*mass*ll(j)*ll(3)/rr)
00258                   *shift_rs(3) )
00259             ) / rr / rr ;
00260         }
00261     }
00262     tmp1.std_spectral_base() ;
00263     tmp1.inc_dzpuis(2) ;
00264 
00265     Sym_tensor tmp2(mp, COV, mp.get_bvect_cart()) ;
00266     tmp2.set_etat_qcq() ;
00267 
00268     for (int i=1; i<=3; i++) {
00269         for (int j=1; j<=3; j++) {
00270             tmp2.set(i,j) = 2. * mass * lapse_bh * lapse_bh
00271           * (ll(1)*shift_rs(1)+ll(2)*shift_rs(2)+ll(3)*shift_rs(3))
00272           * (flat.cov()(i,j)
00273              - (9.+28.*mass/rr+24.*mass*mass/rr/rr)*ll(i)*ll(j))
00274           / 3. / rr / rr ;
00275         }
00276     }
00277     tmp2.std_spectral_base() ;
00278     tmp2.inc_dzpuis(2) ;
00279 
00280     Sym_tensor taij_rs_down(mp, COV, mp.get_bvect_cart()) ;
00281     taij_rs_down.set_etat_qcq() ;
00282 
00283     taij_rs_down = 0.5 * pow(confo, 7.)
00284       * (flat_dshift + curv_dshift + tmp1 + tmp2) / lapconf ;
00285 
00286     taij_rs_down.std_spectral_base() ;
00287     taij_rs_down.annule_domain(0) ;
00288 
00289     Sym_tensor taij_bh_down(mp, COV, mp.get_bvect_cart()) ;
00290     taij_bh_down.set_etat_qcq() ;
00291 
00292     for (int i=1; i<=3; i++) {
00293         for (int j=1; j<=3; j++) {
00294           taij_bh_down.set(i,j) = 2.*pow(lapse_bh,4.)*mass*(2.+3.*mass/rr)
00295         *pow(confo,7.)*(flat.cov()(i,j)-(3.+4.*mass/rr)*ll(i)*ll(j))
00296         /lapconf/3./rr/rr ;
00297         }
00298     }
00299 
00300     taij_bh_down.std_spectral_base() ;
00301     taij_bh_down.inc_dzpuis(2) ;
00302     taij_bh_down.annule_domain(0) ;
00303 
00304     /*
00305     Sym_tensor taij_ks_cov(mp, COV, mp.get_bvect_cart()) ;
00306     taij_ks_cov.set_etat_qcq() ;
00307 
00308     for (int i=1; i<=3; i++) {
00309         for (int j=1; j<=3; j++) {
00310           taij_ks_cov.set(i,j) = 2.*pow(lap_bh2,1.5)*mass * (2.+3.*mass/rr)
00311         * ( flat.cov()(i,j)
00312             - (3.+4.*mass/rr)*ll(i)*ll(j) ) / 3. / rr / rr ;
00313         }
00314     }
00315     taij_ks_cov.std_spectral_base() ;
00316     taij_ks_cov.annule_domain(0) ;
00317     taij_ks_cov.inc_dzpuis(2) ;
00318 
00319     cout << "taij_down(1,1) - taij_ks_cov(1,1) : " << endl ;
00320     cout << taij_down(1,1) - taij_ks_cov(1,1) << endl ;
00321     arrete() ;
00322 
00323     cout << "taij_ks_cov(1,1) : " << endl ;
00324     cout << taij_ks_cov(1,1) << endl ;
00325     arrete() ;
00326 
00327     cout << "taij_down(1,1) : " << endl ;
00328     cout << taij_down(1,1) << endl ;
00329     arrete() ;
00330     */
00331 
00332     Scalar taij_quad_rsrs(mp) ;
00333     taij_quad_rsrs = 0. ;
00334 
00335     for (int i=1; i<=3; i++) {
00336         for (int j=1; j<=3; j++) {
00337             taij_quad_rsrs += taij_rs_down(i,j) * taij_rs(i,j) ;
00338         }
00339     }
00340     taij_quad_rsrs.std_spectral_base() ;
00341 
00342     Scalar taij_quad_rsbh1(mp) ;
00343     taij_quad_rsbh1 = 0. ;
00344 
00345     for (int i=1; i<=3; i++) {
00346         for (int j=1; j<=3; j++) {
00347             taij_quad_rsbh1 += taij_rs_down(i,j) * taij_bh(i,j) ;
00348         }
00349     }
00350     taij_quad_rsbh1.std_spectral_base() ;
00351 
00352     Scalar taij_quad_rsbh2(mp) ;
00353     taij_quad_rsbh2 = 0. ;
00354 
00355     for (int i=1; i<=3; i++) {
00356         for (int j=1; j<=3; j++) {
00357             taij_quad_rsbh2 += taij_bh_down(i,j) * taij_rs(i,j) ;
00358         }
00359     }
00360     taij_quad_rsbh2.std_spectral_base() ;
00361 
00362     taij_quad_rs = taij_quad_rsrs + taij_quad_rsbh1 + taij_quad_rsbh2 ;
00363     taij_quad_rs.std_spectral_base() ;
00364 
00365     Scalar taij_quad_bh(mp) ;
00366     taij_quad_bh = 8.*pow(lapse_bh,10.)*mass*mass*(2.+3.*mass/rr)
00367       *(2.+3.*mass/rr)*pow(confo,12.)/3./pow(rr,4.)/lapconf/lapconf ;
00368     taij_quad_bh.std_spectral_base() ;
00369     taij_quad_bh.inc_dzpuis(4) ;
00370 
00371     taij_quad = taij_quad_rs + taij_quad_bh ;
00372 
00373     taij_quad.std_spectral_base() ;
00374 
00375     /*
00376     Scalar taij_quad_ks(mp) ;
00377     taij_quad_ks = 8. * pow(lap_bh2,3.) * mass * mass * (2.+3.*mass/rr)
00378       * (2.+3.*mass/rr) / 3. / pow(rr, 4.) ;
00379     taij_quad_ks.std_spectral_base() ;
00380     taij_quad_ks.annule_domain(0) ;
00381     taij_quad_ks.inc_dzpuis(4) ;
00382 
00383     cout << "taij_quad - taij_quad_ks : " << endl ;
00384     cout << taij_quad - taij_quad_ks << endl ;
00385     arrete() ;
00386     */
00387     }
00388     else {  // Isotropic coordinates with the maximal slicing
00389 
00390         // Computation of \tilde{A}^{ij}
00391     // -----------------------------
00392 
00393     Scalar divshift(mp) ;
00394     divshift = shift(1).dsdx() + shift(2).dsdy() + shift(3).dsdz() ;
00395     divshift.std_spectral_base() ;
00396 
00397     Sym_tensor flat_taij(mp, CON, mp.get_bvect_cart()) ;
00398     flat_taij.set_etat_qcq() ;
00399 
00400     for (int i=1; i<=3; i++) {
00401         for (int j=1; j<=3; j++) {
00402             flat_taij.set(i,j) = shift(j).deriv(i) + shift(i).deriv(j)
00403           - 2. * divshift * flat.con()(i,j) / 3. ;
00404         }
00405     }
00406 
00407     flat_taij.std_spectral_base() ;
00408 
00409     taij = 0.5 * pow(confo, 7.) * flat_taij / lapconf ;
00410 
00411     taij.std_spectral_base() ;
00412     taij.annule_domain(0) ;
00413 
00414 
00415     // Computation of \tilde{A}^{ij} \tilde{A}_{ij}
00416     // --------------------------------------------
00417 
00418     Sym_tensor flat_dshift(mp, COV, mp.get_bvect_cart()) ;
00419     flat_dshift.set_etat_qcq() ;
00420 
00421     for (int i=1; i<=3; i++) {
00422         for (int j=1; j<=3; j++) {
00423             flat_dshift.set(i,j) = flat.cov()(j,1)*(shift(1).deriv(i))
00424           + flat.cov()(j,2)*(shift(2).deriv(i))
00425           + flat.cov()(j,3)*(shift(3).deriv(i))
00426           + flat.cov()(i,1)*(shift(1).deriv(j))
00427           + flat.cov()(i,2)*(shift(2).deriv(j))
00428           + flat.cov()(i,3)*(shift(3).deriv(j))
00429           - 2. * divshift * flat.cov()(i,j) / 3. ;
00430         }
00431     }
00432 
00433     Sym_tensor taij_down(mp, COV, mp.get_bvect_cart()) ;
00434     taij_down.set_etat_qcq() ;
00435 
00436     for (int i=1; i<=3; i++) {
00437         for (int j=1; j<=3; j++) {
00438           taij_down.set(i,j) = 0.5 * pow(confo, 7.) * flat_dshift(i,j)
00439           / lapconf ;
00440         }
00441     }
00442 
00443     taij_down.std_spectral_base() ;
00444     taij_down.annule_domain(0) ;
00445 
00446     taij_quad = 0. ;
00447 
00448     for (int i=1; i<=3; i++) {
00449         for (int j=1; j<=3; j++) {
00450             taij_quad += taij_down(i,j) * taij(i,j) ;
00451         }
00452     }
00453     taij_quad.std_spectral_base() ;
00454 
00455     }
00456 
00457 }

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