bin_ns_bh_orbit.C

00001 /*
00002  *  Method of class Bin_ns_bh to compute the orbital angular velocity
00003  *  {\tt omega}
00004  *
00005  *    (see file bin_ns_bh.h for documentation).
00006  *
00007  */
00008 
00009 /*
00010  *   Copyright (c) 2003 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 bin_ns_bh_orbit_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bin_ns_bh/bin_ns_bh_orbit.C,v 1.4 2004/06/09 07:26:16 k_taniguchi Exp $" ;
00030 
00031 /*
00032  * $Id: bin_ns_bh_orbit.C,v 1.4 2004/06/09 07:26:16 k_taniguchi Exp $
00033  * $Log: bin_ns_bh_orbit.C,v $
00034  * Revision 1.4  2004/06/09 07:26:16  k_taniguchi
00035  * Minor changes.
00036  *
00037  * Revision 1.3  2004/06/09 06:20:11  k_taniguchi
00038  * Set the standard basis for some Cmp.
00039  *
00040  * Revision 1.2  2004/03/25 10:28:58  j_novak
00041  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00042  *
00043  * Revision 1.1  2003/10/24 16:57:43  k_taniguchi
00044  * Method for the calculation of the orbital angular velocity
00045  *
00046  *
00047  * $Header: /cvsroot/Lorene/C++/Source/Bin_ns_bh/bin_ns_bh_orbit.C,v 1.4 2004/06/09 07:26:16 k_taniguchi Exp $
00048  *
00049  */
00050 
00051 // C headers
00052 #include <math.h>
00053 
00054 // Lorene headers
00055 #include "bin_ns_bh.h"
00056 #include "eos.h"
00057 #include "param.h"
00058 #include "utilitaires.h"
00059 #include "unites.h"
00060 
00061 double  fonc_bin_ns_bh_orbit(double , const Param& ) ;
00062 
00063 //*************************************************************************
00064 
00065 void Bin_ns_bh::orbit_omega(double fact_omeg_min, double fact_omeg_max) {
00066 
00067   using namespace Unites ; 
00068 
00069     //------------------------------------------------------------------
00070     // Evaluation of various quantities at the center of a neutron star
00071     //------------------------------------------------------------------
00072 
00073     double dnulg, asn2, dasn2, ny, dny, npn, dnpn ;
00074 
00075     const Map& mp = star.get_mp() ;
00076 
00077     const Cmp& loggam = star.get_loggam()() ;
00078     const Cmp& nnn = star.get_nnn()() ;
00079     const Cmp& confpsi = star.get_confpsi()() ;
00080     const Tenseur& shift = star.get_shift() ;
00081 
00082     Cmp confpsi_q = pow(confpsi, 4.) ;
00083     confpsi_q.std_base_scal() ;
00084 
00085     //-----------------------------------------------------------------
00086     // Calculation of d/dX( nu + ln(Gamma) ) at the center of the star
00087     //  ---> dnulg
00088     //-----------------------------------------------------------------
00089 
00090     // Factor for the coordinate transformation x --> X :
00091     double factx ;
00092     if (fabs(mp.get_rot_phi()) < 1.e-14) {
00093     factx = 1. ;
00094     }
00095     else {
00096     if (fabs(mp.get_rot_phi() - M_PI) < 1.e-14) {
00097         factx = - 1. ;
00098     }
00099     else {
00100         cout << "Bin_ns_bh::orbit_omega : unknown value of rot_phi !"
00101          << endl ;
00102         abort() ;
00103     }
00104     }
00105 
00106     Cmp tmp = log( nnn ) + loggam ;
00107     tmp.std_base_scal() ;
00108 
00109     // ... gradient
00110     dnulg = factx * tmp.dsdx()(0, 0, 0, 0) ;
00111 
00112     //------------------------------------------------------------
00113     // Calculation of A^2/N^2 at the center of the star ---> asn2
00114     //------------------------------------------------------------
00115     double nc = nnn(0, 0, 0, 0) ;
00116     double a2c = confpsi_q(0, 0, 0, 0) ;
00117     asn2 = a2c / (nc * nc) ;
00118 
00119     if ( star.is_relativistic() ) {
00120 
00121     //------------------------------------------------------------------
00122     // Calculation of d/dX(A^2/N^2) at the center of the star ---> dasn
00123     //------------------------------------------------------------------
00124     double da2c = factx * confpsi_q.dsdx()(0, 0, 0, 0) ;
00125     double dnc =  factx * nnn.dsdx()(0, 0, 0, 0) ;
00126 
00127     dasn2 = ( da2c - 2 * a2c / nc * dnc ) / (nc*nc) ;
00128 
00129     //------------------------------------------------------
00130     // Calculation of N^Y at the center of the star ---> ny
00131     //------------------------------------------------------
00132     ny = shift(1)(0, 0, 0, 0) ;
00133 
00134     //-----------------------------------------------------------
00135     // Calculation of dN^Y/dX at the center of the star ---> dny
00136     //-----------------------------------------------------------
00137     dny = factx * shift(1).dsdx()(0, 0, 0, 0) ;
00138 
00139     //--------------------------------------------
00140     // Calculation of (N^X)^2 + (N^Y)^2 + (N^Z)^2
00141     //  at the center of the star ---> npn
00142     //--------------------------------------------
00143     tmp = flat_scalar_prod(shift, shift)() ;
00144     npn = tmp(0, 0, 0, 0) ;
00145 
00146     //----------------------------------------------------
00147     // Calculation of d/dX( (N^X)^2 + (N^Y)^2 + (N^Z)^2 )
00148     //  at the center of the star ---> dnpn
00149     //----------------------------------------------------
00150     dnpn = factx * tmp.dsdx()(0, 0, 0, 0) ;
00151 
00152     }  // Finish of the relativistic case
00153     else {
00154     cout << "Bin_ns_bh::orbit_omega : "
00155          << "It should be the relativistic calculation !" << endl ;
00156     abort() ;
00157     }
00158 
00159     cout << "Bin_ns_bh::orbit_omega: central d(nu+log(Gam))/dX : " 
00160      << dnulg << endl ; 
00161     cout << "Bin_ns_bh::orbit_omega: central A^2/N^2 : " << asn2 << endl ; 
00162     cout << "Bin_ns_bh::orbit_omega: central d(A^2/N^2)/dX : "
00163      << dasn2 << endl ; 
00164     cout << "Bin_ns_bh::orbit_omega: central N^Y : " << ny << endl ; 
00165     cout << "Bin_ns_bh::orbit_omega: central dN^Y/dX : " << dny << endl ; 
00166     cout << "Bin_ns_bh::orbit_omega: central N.N : " << npn << endl ; 
00167     cout << "Bin_ns_bh::orbit_omega: central d(N.N)/dX : "
00168      << dnpn << endl ; 
00169 
00170     //------------------------------------------------------
00171     // Start of calculation of the orbital angular velocity
00172     //------------------------------------------------------
00173     int relat = ( star.is_relativistic() ) ? 1 : 0 ;
00174 
00175     Param parf ;
00176     parf.add_int(relat) ;
00177     parf.add_double( (star.get_mp()).get_ori_x(), 0) ;
00178     parf.add_double( dnulg, 1) ;
00179     parf.add_double( asn2, 2) ;
00180     parf.add_double( dasn2, 3) ;
00181     parf.add_double( ny, 4) ;
00182     parf.add_double( dny, 5) ;
00183     parf.add_double( npn, 6) ;
00184     parf.add_double( dnpn, 7) ;
00185     parf.add_double( x_axe, 8) ;
00186 
00187     double omega1 = fact_omeg_min * omega ;
00188     double omega2 = fact_omeg_max * omega ;
00189 
00190     cout << "Bin_ns_bh::orbit_omega: omega1, omega2 [rad/s] : "
00191      << omega1 * f_unit << "  " << omega2 * f_unit << endl ;
00192 
00193     // Search for the various zeros in the interval [omega1,omega2]
00194     // ------------------------------------------------------------
00195     int nsub = 50 ;  // total number of subdivisions of the interval
00196     Tbl* azer = 0x0 ;
00197     Tbl* bzer = 0x0 ;
00198     zero_list(fonc_bin_ns_bh_orbit, parf, omega1, omega2, nsub,
00199           azer, bzer) ;
00200 
00201     // Search for the zero closest to the previous value of omega
00202     // ----------------------------------------------------------
00203     double omeg_min, omeg_max ;
00204     int nzer = azer->get_taille() ; // number of zeros found by zero_list
00205     cout << "Bin_ns_bh:orbit_omega : " << nzer <<
00206     "zero(s) found in the interval [omega1,  omega2]." << endl ;
00207     cout << "omega, omega1, omega2 : " << omega << "  " << omega1
00208      << "  " << omega2 << endl ; 
00209     cout << "azer : " << *azer << endl ;
00210     cout << "bzer : " << *bzer << endl ;
00211 
00212     if (nzer == 0) {
00213     cout << "Bin_ns_bh::orbit_omega: WARNING : "
00214          << "no zero detected in the interval" << endl
00215          << "   [" << omega1 * f_unit << ", " 
00216          << omega2 * f_unit << "]  rad/s  !" << endl ;
00217     omeg_min = omega1 ;
00218     omeg_max = omega2 ;
00219     }
00220     else {
00221     double dist_min = fabs(omega2 - omega1) ;
00222     int i_dist_min = -1 ;
00223     for (int i=0; i<nzer; i++) {
00224         // Distance of previous value of omega from the center of the
00225         //  interval [azer(i), bzer(i)]
00226         double dist = fabs( omega - 0.5 * ( (*azer)(i) + (*bzer)(i) ) ) ;
00227 
00228         if (dist < dist_min) {
00229         dist_min = dist ;
00230         i_dist_min = i ;
00231         }
00232     }
00233     omeg_min = (*azer)(i_dist_min) ;
00234     omeg_max = (*bzer)(i_dist_min) ;
00235     }
00236 
00237     delete azer ; // Tbl allocated by zero_list
00238     delete bzer ; //
00239 
00240     cout << "Bin_ns_bh:orbit_omega : "
00241      << "interval selected for the search of the zero : "
00242      << endl << "  [" << omeg_min << ", " << omeg_max << "] = ["
00243      << omeg_min * f_unit << ", " << omeg_max * f_unit << "] rad/s "
00244      << endl ;
00245 
00246     // Computation of the zero in the selected interval by the secant method
00247     // ---------------------------------------------------------------------
00248 
00249     int nitermax = 200 ;
00250     int niter ;
00251     double precis = 1.e-13 ;
00252     omega = zerosec_b(fonc_bin_ns_bh_orbit, parf, omeg_min, omeg_max,
00253               precis, nitermax, niter) ;
00254 
00255     cout << "Bin_ns_bh::orbit_omega : "
00256      << "Number of iterations in zerosec for omega : "
00257      << niter << endl ;
00258 
00259     cout << "Bin_ns_bh::orbit_omega : omega [rad/s] : "
00260      << omega * f_unit << endl ;
00261 
00262 }
00263 
00264 //***********************************************************
00265 //  Function used for search of the orbital angular velocity
00266 //***********************************************************
00267 
00268 double fonc_bin_ns_bh_orbit(double om, const Param& parf) {
00269 
00270     int relat = parf.get_int() ;
00271 
00272     double xc = parf.get_double(0) ;
00273     double dnulg = parf.get_double(1) ;
00274     double asn2 = parf.get_double(2) ;
00275     double dasn2 = parf.get_double(3) ;
00276     double ny = parf.get_double(4) ;
00277     double dny = parf.get_double(5) ;
00278     double npn = parf.get_double(6) ;
00279     double dnpn = parf.get_double(7) ;
00280     double x_axe = parf.get_double(8) ;
00281 
00282     double xx = xc - x_axe ;
00283     double om2 = om*om ;
00284 
00285     double dphi_cent ;
00286 
00287     if (relat == 1) {
00288     double bpb = om2 * xx*xx - 2*om * ny * xx + npn ;
00289 
00290     dphi_cent = ( asn2* ( om* (ny + xx*dny) - om2*xx - 0.5*dnpn )
00291              - 0.5*bpb* dasn2 )
00292                / ( 1 - asn2 * bpb ) ;
00293     }
00294     else {
00295     cout << "Bin_ns_bh::orbit_omega : "
00296          << "It should be the relativistic calculation !" << endl ;
00297     abort() ;
00298     }
00299 
00300     return dnulg + dphi_cent ;
00301 
00302 }

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