et_bin_bhns_extr_kinema.C

00001 /*
00002  *  Method Et_bin_bhns_extr::kinematics_extr_ks
00003  *  and Et_bin_bhns_extr::kinematics_extr_cf
00004  *
00005  *    (see file et_bin_bhns_extr.h for documentation).
00006  *
00007  */
00008 
00009 /*
00010  *   Copyright (c) 2004-2005 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 et_bin_bhns_extr_kinema_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_kinema.C,v 1.2 2005/02/28 23:15:09 k_taniguchi Exp $" ;
00030 
00031 /*
00032  * $Id: et_bin_bhns_extr_kinema.C,v 1.2 2005/02/28 23:15:09 k_taniguchi Exp $
00033  * $Log: et_bin_bhns_extr_kinema.C,v $
00034  * Revision 1.2  2005/02/28 23:15:09  k_taniguchi
00035  * Modification to include the case of the conformally flat background metric
00036  *
00037  * Revision 1.1  2004/11/30 20:49:58  k_taniguchi
00038  * *** empty log message ***
00039  *
00040  *
00041  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_bhns_extr_kinema.C,v 1.2 2005/02/28 23:15:09 k_taniguchi Exp $
00042  *
00043  */
00044 
00045 // Lorene headers
00046 #include "et_bin_bhns_extr.h"
00047 #include "etoile.h"
00048 #include "coord.h"
00049 #include "unites.h"
00050 
00051 void Et_bin_bhns_extr::kinematics_extr(double omega, const double& mass,
00052                        const double& sepa) {
00053 
00054   using namespace Unites ;
00055 
00056     if (kerrschild) {
00057 
00058         int nz = mp.get_mg()->get_nzone() ;
00059     int nzm1 = nz - 1 ;
00060 
00061     // --------------------
00062     // Computation of B^i/N
00063     // --------------------
00064 
00065     //  1/ Computation of  - omega m^i
00066 
00067     const Coord& xa = mp.xa ;
00068     const Coord& ya = mp.ya ;
00069 
00070     bsn.set_etat_qcq() ;
00071 
00072     bsn.set(0) = omega * ya ;
00073     bsn.set(1) = - omega * xa ;
00074     bsn.set(2) = 0 ;
00075 
00076     bsn.annule(nzm1, nzm1) ;    // set to zero in the ZEC
00077 
00078     //  2/ Addition of shift and division by lapse
00079 
00080     bsn = ( bsn + shift ) / nnn ;
00081 
00082     bsn.annule(nzm1, nzm1) ;    // set to zero in the ZEC
00083     bsn.set_std_base() ;   // set the bases for spectral expansions
00084 
00085     //-------------------------
00086     // Centrifugal potentatial
00087     //-------------------------
00088 
00089     const Coord& xx = mp.x ;
00090     const Coord& yy = mp.y ;
00091     const Coord& zz = mp.z ;
00092 
00093     Tenseur r_bh(mp) ;
00094     r_bh.set_etat_qcq() ;
00095     r_bh.set() = pow( (xx+sepa)*(xx+sepa) + yy*yy + zz*zz, 0.5) ;
00096     r_bh.set_std_base() ;
00097 
00098     Tenseur xx_cov(mp, 1, COV, ref_triad) ;
00099     xx_cov.set_etat_qcq() ;
00100     xx_cov.set(0) = xx + sepa ;
00101     xx_cov.set(1) = yy ;
00102     xx_cov.set(2) = zz ;
00103     xx_cov.set_std_base() ;
00104 
00105     Tenseur xsr_cov(mp, 1, COV, ref_triad) ;
00106     xsr_cov = xx_cov / r_bh ;
00107     xsr_cov.set_std_base() ;
00108 
00109     Tenseur msr(mp) ;
00110     msr = ggrav * mass / r_bh ;
00111     msr.set_std_base() ;
00112 
00113     if (relativistic) {
00114 
00115         // Lorentz factor between the co-orbiting observer and
00116         // the Eulerian one
00117 
00118         Tenseur tmp1(mp) ;
00119         tmp1.set_etat_qcq() ;
00120         tmp1.set() = 0 ;
00121         tmp1.set_std_base() ;
00122 
00123         for (int i=0; i<3; i++)
00124             tmp1.set() += xsr_cov(i) % bsn(i) ;
00125 
00126         tmp1.set_std_base() ;
00127 
00128         Tenseur tmp2 = 2.*msr % tmp1 % tmp1 ;
00129         tmp2.set_std_base() ;
00130 
00131         for (int i=0; i<3; i++)
00132             tmp2.set() += bsn(i) % bsn(i) ;
00133 
00134         tmp2 = a_car % tmp2 ;
00135 
00136         Tenseur gam0 = 1 / sqrt( 1 - tmp2 ) ;
00137 
00138         pot_centri = - log( gam0 ) ;
00139 
00140     }
00141     else {
00142         cout << "BH-NS binary system should be relativistic !!!" << endl ;
00143         abort() ;
00144     }
00145 
00146     pot_centri.annule(nzm1, nzm1) ; // set to zero in the external domain
00147     pot_centri.set_std_base() ;   // set the bases for spectral expansions
00148 
00149     // The derived quantities are obsolete
00150     // -----------------------------------
00151 
00152     Etoile_bin::del_deriv() ;
00153 
00154     }
00155     else {
00156 
00157         int nz = mp.get_mg()->get_nzone() ;
00158     int nzm1 = nz - 1 ;
00159 
00160     // --------------------
00161     // Computation of B^i/N
00162     // --------------------
00163 
00164     //  1/ Computation of  - omega m^i
00165 
00166     const Coord& xa = mp.xa ;
00167     const Coord& ya = mp.ya ;
00168 
00169     bsn.set_etat_qcq() ;
00170 
00171     bsn.set(0) = omega * ya ;
00172     bsn.set(1) = - omega * xa ;
00173     bsn.set(2) = 0 ;
00174 
00175     bsn.annule(nzm1, nzm1) ;    // set to zero in the ZEC
00176 
00177     //  2/ Addition of shift and division by lapse
00178 
00179     bsn = ( bsn + shift ) / nnn ;
00180 
00181     bsn.annule(nzm1, nzm1) ;    // set to zero in the ZEC
00182     bsn.set_std_base() ;   // set the bases for spectral expansions
00183 
00184     //-------------------------
00185     // Centrifugal potentatial
00186     //-------------------------
00187 
00188     if (relativistic) {
00189 
00190         // Lorentz factor between the co-orbiting observer and
00191         // the Eulerian one
00192 
00193         Tenseur tmp(mp) ;
00194         tmp.set_etat_qcq() ;
00195         tmp.set() = 0. ;
00196         tmp.set_std_base() ;
00197 
00198         for (int i=0; i<3; i++)
00199             tmp.set() += bsn(i) % bsn(i) ;
00200 
00201         tmp = a_car % tmp ;
00202 
00203         Tenseur gam0 = 1 / sqrt( 1 - tmp ) ;
00204 
00205         pot_centri = - log( gam0 ) ;
00206 
00207     }
00208     else {
00209         cout << "BH-NS binary system should be relativistic !!!" << endl ;
00210         abort() ;
00211     }
00212 
00213     pot_centri.annule(nzm1, nzm1) ; // set to zero in the external domain
00214     pot_centri.set_std_base() ;   // set the bases for spectral expansions
00215 
00216     // The derived quantities are obsolete
00217     // -----------------------------------
00218 
00219     Etoile_bin::del_deriv() ;
00220 
00221     }
00222 
00223 }

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