bin_bhns_orbit.C

00001 /*
00002  *  Methods of class Bin_bhns to compute the orbital angular velocity
00003  *
00004  *    (see file bin_bhns.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 bin_bhns_orbit_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bin_bhns/bin_bhns_orbit.C,v 1.2 2008/05/15 19:01:28 k_taniguchi Exp $" ;
00029 
00030 /*
00031  * $Id: bin_bhns_orbit.C,v 1.2 2008/05/15 19:01:28 k_taniguchi Exp $
00032  * $Log: bin_bhns_orbit.C,v $
00033  * Revision 1.2  2008/05/15 19:01:28  k_taniguchi
00034  * Change of some parameters.
00035  *
00036  * Revision 1.1  2007/06/22 01:10:20  k_taniguchi
00037  * *** empty log message ***
00038  *
00039  *
00040  * $Header: /cvsroot/Lorene/C++/Source/Bin_bhns/bin_bhns_orbit.C,v 1.2 2008/05/15 19:01:28 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 "bin_bhns.h"
00052 #include "param.h"
00053 #include "utilitaires.h"
00054 #include "unites.h"
00055 
00056 double func_binbhns_orbit_ks(double , const Param& ) ;
00057 double func_binbhns_orbit_is(double , const Param& ) ;
00058 
00059 //**********************************************************************
00060 
00061 void Bin_bhns::orbit_omega(double fact_omeg_min, double fact_omeg_max) {
00062 
00063     using namespace Unites ;
00064 
00065     if (hole.is_kerrschild()) {
00066 
00067       //--------------------------------------------------------------------
00068       // Evaluation of various quantities at the center of the neutron star
00069       //--------------------------------------------------------------------
00070 
00071       double dnulg, p6sl2, dp6sl2 ;
00072       double shiftx, shifty, dshiftx, dshifty, shift2, dshift2 ;
00073       double x_orb, y_orb, y_separ, xbh_orb, mhsr ;
00074 
00075       const Map& mp = star.get_mp() ;
00076 
00077       const Scalar& lapconf = star.get_lapconf_tot() ;
00078       const Scalar& lapconf_auto = star.get_lapconf_auto() ;
00079       const Scalar& confo = star.get_confo_tot() ;
00080       const Scalar& confo_auto = star.get_confo_auto() ;
00081       const Scalar& loggam = star.get_loggam() ;
00082       const Vector& shift = star.get_shift_tot() ;
00083       const Vector& shift_auto = star.get_shift_auto() ;
00084 
00085       const Vector& dlapconf_comp = star.get_d_lapconf_comp() ;
00086       const Vector& dconfo_comp = star.get_d_confo_comp() ;
00087       const Tensor& dshift_comp = star.get_d_shift_comp() ;
00088 
00089       const double& massbh = hole.get_mass_bh() ;
00090       double mass = ggrav * massbh ;
00091 
00092       //----------------------------------------------------------
00093       // Calculation of d/dX( ln(lapconf) - ln(psi) + ln(Gamma) )
00094       //  at the center of NS
00095       //----------------------------------------------------------
00096 
00097       // Factor to translate x --> X
00098       double factx ;
00099       if (fabs(mp.get_rot_phi()) < 1.e-14) {
00100         factx = 1. ;
00101       }
00102       else {
00103         if (fabs(mp.get_rot_phi() - M_PI) < 1.e-14) {
00104       factx = - 1. ;
00105     }
00106     else {
00107       cout << "Bin_bhns::orbit_omega : unknown value of rot_phi !"
00108            << endl ;
00109       abort() ;
00110     }
00111       }
00112 
00113       Scalar tmp1(mp) ;
00114       tmp1 = loggam ;
00115       tmp1.std_spectral_base() ;
00116 
00117       // d/dX tmp1
00118       dnulg = factx * ( ((lapconf_auto.dsdx()).val_grid_point(0,0,0,0)
00119              + dlapconf_comp(1).val_grid_point(0,0,0,0))
00120             / lapconf.val_grid_point(0,0,0,0)
00121             - ((confo_auto.dsdx()).val_grid_point(0,0,0,0)
00122                + dconfo_comp(1).val_grid_point(0,0,0,0))
00123             / confo.val_grid_point(0,0,0,0)
00124             + (tmp1.dsdx()).val_grid_point(0,0,0,0) ) ;
00125 
00126 
00127       //----------------------------------------------------
00128       // Calculation of psi^6/lapconf^2 at the center of NS
00129       //----------------------------------------------------
00130 
00131       double lapconf_c = lapconf.val_grid_point(0,0,0,0) ;
00132       double confo_c = confo.val_grid_point(0,0,0,0) ;
00133 
00134       p6sl2 = pow(confo_c,6.) / lapconf_c / lapconf_c ;
00135 
00136 
00137       //----------------------------------------------------------
00138       // Calculation of d/dX(psi^6/lapconf^2) at the center of NS
00139       //----------------------------------------------------------
00140 
00141       double dlapconf_c = factx *
00142     ( (lapconf_auto.dsdx()).val_grid_point(0,0,0,0)
00143       + dlapconf_comp(1).val_grid_point(0,0,0,0) ) ;
00144 
00145       double dpsi6_c = 6. * factx * pow(confo_c,5.)
00146     * ((confo_auto.dsdx()).val_grid_point(0,0,0,0)
00147        + dconfo_comp(1).val_grid_point(0,0,0,0)) ;
00148 
00149       dp6sl2 = (dpsi6_c - 2.*pow(confo_c,6.)*dlapconf_c/lapconf_c)
00150     / lapconf_c / lapconf_c ;
00151 
00152 
00153       //--------------------------------------------------------
00154       // Calculation of shift^X and shift^Y at the center of NS
00155       //--------------------------------------------------------
00156 
00157       shiftx = shift(1).val_grid_point(0,0,0,0) ;
00158       shifty = shift(2).val_grid_point(0,0,0,0) ;
00159 
00160 
00161       //------------------------------------------------------------------
00162       // Calculation of d shift^X/dX and d shift^Y/dX at the center of NS
00163       //------------------------------------------------------------------
00164 
00165       dshiftx = factx * ((shift_auto(1).dsdx()).val_grid_point(0,0,0,0)
00166              + dshift_comp(1,1).val_grid_point(0,0,0,0)) ;
00167 
00168       dshifty = factx * ((shift_auto(2).dsdx()).val_grid_point(0,0,0,0)
00169              + dshift_comp(1,2).val_grid_point(0,0,0,0)) ;
00170 
00171 
00172       //--------------------------------------------------------
00173       // Calculation of (shift^X)^2 + (shift^Y)^2 + (shift^Z)^2
00174       //  at the center of NS
00175       //--------------------------------------------------------
00176 
00177       Scalar tmp2(mp) ;
00178       tmp2 = shift(1)%shift(1) + shift(2)%shift(2) + shift(3)%shift(3) ;
00179       tmp2.std_spectral_base() ;
00180 
00181       shift2 = tmp2.val_grid_point(0,0,0,0) ;
00182 
00183 
00184       //----------------------------------------------------------------
00185       // Calculation of d/dX( (shift^X)^2 + (shift^Y)^2 + (shift^Z)^2 )
00186       //  at the center of NS
00187       //----------------------------------------------------------------
00188 
00189       dshift2 = 2.*factx*((shift(1).val_grid_point(0,0,0,0))
00190               * ((shift_auto(1).dsdx()).val_grid_point(0,0,0,0)
00191                  + dshift_comp(1,1).val_grid_point(0,0,0,0))
00192               +(shift(2).val_grid_point(0,0,0,0))
00193               * ((shift_auto(2).dsdx()).val_grid_point(0,0,0,0)
00194                  + dshift_comp(1,2).val_grid_point(0,0,0,0))
00195               +(shift(3).val_grid_point(0,0,0,0))
00196               * ((shift_auto(3).dsdx()).val_grid_point(0,0,0,0)
00197                  + dshift_comp(1,3).val_grid_point(0,0,0,0)) ) ;
00198 
00199 
00200       //-------------------------
00201       // Information of position
00202       //-------------------------
00203 
00204       x_orb = (star.get_mp()).get_ori_x() - x_rot ;
00205       y_orb = (star.get_mp()).get_ori_y() - y_rot ;
00206       y_separ = (star.get_mp()).get_ori_y() - (hole.get_mp()).get_ori_y() ;
00207 
00208       xbh_orb = (hole.get_mp()).get_ori_x() - x_rot ;
00209 
00210        //------------------------------
00211       // Calculation of H_BH / r_bh^4
00212       //------------------------------
00213 
00214       mhsr = mass / pow( separ*separ+y_separ*y_separ, 2.5 ) ;
00215 
00216       // Output
00217       // ------
00218 
00219       cout << "Bin_bhns::orbit_omega: central d(log(lap)+log(Gam))/dX : "
00220        << dnulg << endl ;
00221       cout << "Bin_bhns::orbit_omega: central psi^6/lapconf^2 :         "
00222        << p6sl2 << endl ;
00223       cout << "Bin_bhns::orbit_omega: central d(psi^6/lapconf^2)/dX :   "
00224        << dp6sl2 << endl ;
00225       cout << "Bin_bhns::orbit_omega: central shift^X :                 "
00226        << shiftx << endl ;
00227       cout << "Bin_bhns::orbit_omega: central shift^Y :                 "
00228        << shifty << endl ;
00229       cout << "Bin_bhns::orbit_omega: central d(shift^X)/dX :           "
00230        << dshiftx << endl ;
00231       cout << "Bin_bhns::orbit_omega: central d(shift^Y)/dX :           "
00232        << dshifty << endl ;
00233       cout << "Bin_bhns::orbit_omega: central shift^i shift_i :         "
00234        << shift2 << endl ;
00235       cout << "Bin_bhns::orbit_omega: central d(shift^i shift_i)/dX :   "
00236        << dshift2 << endl ;
00237 
00238 
00239       //---------------------------------------------------------------//
00240       //          Calculation of the orbital angular velocity          //
00241       //---------------------------------------------------------------//
00242 
00243       Param parorb ;
00244       parorb.add_double(dnulg, 0) ;
00245       parorb.add_double(p6sl2, 1) ;
00246       parorb.add_double(dp6sl2, 2) ;
00247       parorb.add_double(shiftx, 3) ;
00248       parorb.add_double(shifty, 4) ;
00249       parorb.add_double(dshiftx, 5) ;
00250       parorb.add_double(dshifty, 6) ;
00251       parorb.add_double(shift2, 7) ;
00252       parorb.add_double(dshift2, 8) ;
00253       parorb.add_double(x_orb, 9) ;
00254       parorb.add_double(y_orb, 10) ;
00255       parorb.add_double(separ, 11) ;
00256       parorb.add_double(y_separ, 12) ;
00257       parorb.add_double(xbh_orb, 13) ;
00258       parorb.add_double(mhsr, 14) ;
00259 
00260       double omega1 = fact_omeg_min * omega ;
00261       double omega2 = fact_omeg_max * omega ;
00262 
00263       cout << "Bin_bhns::orbit_omega: omega1,  omega2 [rad/s] : "
00264        << omega1 * f_unit << "  " << omega2 * f_unit << endl ;
00265 
00266       // Search for the various zeros in the interval [omega1,omega2]
00267       // ------------------------------------------------------------
00268 
00269       int nsub = 50 ;  // total number of subdivisions of the interval
00270       Tbl* azer = 0x0 ;
00271       Tbl* bzer = 0x0 ;
00272       zero_list(func_binbhns_orbit_ks, parorb, omega1, omega2, nsub,
00273         azer, bzer) ;
00274 
00275       // Search for the zero closest to the previous value of omega
00276       // ----------------------------------------------------------
00277 
00278       double omeg_min, omeg_max ;
00279       int nzer = azer->get_taille() ; // number of zeros found by zero_list
00280 
00281       cout << "Bin_bhns::orbit_omega: " << nzer
00282        << " zero(s) found in the interval [omega1,  omega2]." << endl ;
00283       cout << "omega, omega1, omega2 : " << omega << "  " << omega1
00284        << "  " << omega2 << endl ;
00285       cout << "azer : " << *azer << endl ;
00286       cout << "bzer : " << *bzer << endl ;
00287 
00288       if (nzer == 0) {
00289         cout <<
00290       "Bin_bhns::orbit_omega: WARNING : no zero detected in the interval"
00291          << endl << "   [" << omega1 * f_unit << ", "
00292          << omega2 * f_unit << "]  rad/s  !" << endl ;
00293     omeg_min = omega1 ;
00294     omeg_max = omega2 ;
00295       }
00296       else {
00297         double dist_min = fabs(omega2 - omega1) ;
00298     int i_dist_min = -1 ;
00299     for (int i=0; i<nzer; i++) {
00300       // Distance of previous value of omega from the center of the
00301       //  interval [azer(i), bzer(i)]
00302 
00303       double dist = fabs( omega - 0.5 * ( (*azer)(i) + (*bzer)(i) ) ) ;
00304 
00305       if (dist < dist_min) {
00306         dist_min = dist ;
00307         i_dist_min = i ;
00308       }
00309     }
00310     omeg_min = (*azer)(i_dist_min) ;
00311     omeg_max = (*bzer)(i_dist_min) ;
00312       }
00313 
00314       delete azer ; // Tbl allocated by zero_list
00315       delete bzer ;
00316 
00317       cout << "Bin_bhns::orbit_omega: interval selected for the search of the zero : "
00318        << endl << "  [" << omeg_min << ", " << omeg_max << "]  =  ["
00319        << omeg_min * f_unit << ", " << omeg_max * f_unit << "] rad/s " << endl ;
00320 
00321       // Computation of the zero in the selected interval by the secant method
00322       // ---------------------------------------------------------------------
00323 
00324       int nitermax = 200 ;
00325       int niter ;
00326       double precis = 1.e-13 ;
00327       omega = zerosec_b(func_binbhns_orbit_ks, parorb, omeg_min, omeg_max,
00328             precis, nitermax, niter) ;
00329 
00330       cout << "Bin_bhns::orbit_omega: Number of iterations in zerosec for omega : "
00331        << niter << endl ;
00332 
00333       cout << "Bin_bhns::orbit_omega: omega [rad/s] : "
00334        << omega * f_unit << endl ;
00335 
00336 
00337     }
00338     else { // Isotropic coordinates with the maximal slicing
00339 
00340       //--------------------------------------------------------------------
00341       // Evaluation of various quantities at the center of the neutron star
00342       //--------------------------------------------------------------------
00343 
00344       double dnulg, p6sl2, dp6sl2 ;
00345       double shiftx, shifty, dshiftx, dshifty, shift2, dshift2 ;
00346       double x_orb, y_orb ;
00347 
00348       const Map& mp = star.get_mp() ;
00349 
00350       const Scalar& lapconf = star.get_lapconf_tot() ;
00351       const Scalar& lapconf_auto = star.get_lapconf_auto() ;
00352       const Scalar& confo = star.get_confo_tot() ;
00353       const Scalar& confo_auto = star.get_confo_auto() ;
00354       const Scalar& loggam = star.get_loggam() ;
00355       const Vector& shift = star.get_shift_tot() ;
00356       const Vector& shift_auto = star.get_shift_auto() ;
00357 
00358       const Vector& dlapconf_comp = star.get_d_lapconf_comp() ;
00359       const Vector& dconfo_comp = star.get_d_confo_comp() ;
00360       const Tensor& dshift_comp = star.get_d_shift_comp() ;
00361 
00362       //----------------------------------------------------------
00363       // Calculation of d/dX( ln(lapconf) - ln(psi) + ln(Gamma) )
00364       //  at the center of NS
00365       //----------------------------------------------------------
00366 
00367       // Factor to translate x --> X
00368       double factx ;
00369       if (fabs(mp.get_rot_phi()) < 1.e-14) {
00370         factx = 1. ;
00371       }
00372       else {
00373         if (fabs(mp.get_rot_phi() - M_PI) < 1.e-14) {
00374       factx = - 1. ;
00375     }
00376     else {
00377       cout << "Bin_bhns::orbit_omega: unknown value of rot_phi !"
00378            << endl ;
00379       abort() ;
00380     }
00381       }
00382 
00383       Scalar tmp1(mp) ;
00384       tmp1 = loggam ;
00385       tmp1.std_spectral_base() ;
00386 
00387       // d/dX tmp1
00388       dnulg = factx * ( ((lapconf_auto.dsdx()).val_grid_point(0,0,0,0)
00389              + dlapconf_comp(1).val_grid_point(0,0,0,0))
00390             / lapconf.val_grid_point(0,0,0,0)
00391             - ((confo_auto.dsdx()).val_grid_point(0,0,0,0)
00392              + dconfo_comp(1).val_grid_point(0,0,0,0))
00393             / confo.val_grid_point(0,0,0,0)
00394             + (tmp1.dsdx()).val_grid_point(0,0,0,0) ) ;
00395 
00396 
00397       //----------------------------------------------------
00398       // Calculation of psi^6/lapconf^2 at the center of NS
00399       //----------------------------------------------------
00400 
00401       double lapconf_c = lapconf.val_grid_point(0,0,0,0) ;
00402       double confo_c = confo.val_grid_point(0,0,0,0) ;
00403 
00404       p6sl2 = pow(confo_c,6.) / lapconf_c / lapconf_c ;
00405 
00406 
00407       //----------------------------------------------------------
00408       // Calculation of d/dX(psi^6/lapconf^2) at the center of NS
00409       //----------------------------------------------------------
00410 
00411       double dlapconf_c = factx *
00412     ( (lapconf_auto.dsdx()).val_grid_point(0,0,0,0)
00413       + dlapconf_comp(1).val_grid_point(0,0,0,0) ) ;
00414 
00415       double dpsi6_c = 6. * factx * pow(confo_c,5.)
00416     * ((confo_auto.dsdx()).val_grid_point(0,0,0,0)
00417        + dconfo_comp(1).val_grid_point(0,0,0,0)) ;
00418 
00419       dp6sl2 = (dpsi6_c - 2.*pow(confo_c,6.)*dlapconf_c/lapconf_c)
00420     / lapconf_c / lapconf_c ;
00421 
00422 
00423       //--------------------------------------------------------
00424       // Calculation of shift^X and shift^Y at the center of NS
00425       //--------------------------------------------------------
00426 
00427       shiftx = shift(1).val_grid_point(0,0,0,0) ;
00428       shifty = shift(2).val_grid_point(0,0,0,0) ;
00429 
00430 
00431       //------------------------------------------------------------------
00432       // Calculation of d shift^X/dX and d shift^Y/dX at the center of NS
00433       //------------------------------------------------------------------
00434 
00435       dshiftx = factx * ( (shift_auto(1).dsdx()).val_grid_point(0,0,0,0)
00436               + dshift_comp(1,1).val_grid_point(0,0,0,0) ) ;
00437 
00438       dshifty = factx * ( (shift_auto(2).dsdx()).val_grid_point(0,0,0,0)
00439               + dshift_comp(1,2).val_grid_point(0,0,0,0) ) ;
00440 
00441 
00442       //--------------------------------------------------------
00443       // Calculation of (shift^X)^2 + (shift^Y)^2 + (shift^Z)^2
00444       //  at the center of NS
00445       //--------------------------------------------------------
00446 
00447       Scalar tmp2(mp) ;
00448       tmp2 = shift(1)%shift(1) + shift(2)%shift(2) + shift(3)%shift(3) ;
00449       tmp2.std_spectral_base() ;
00450 
00451       shift2 = tmp2.val_grid_point(0,0,0,0) ;
00452 
00453 
00454       //----------------------------------------------------------------
00455       // Calculation of d/dX( (shift^X)^2 + (shift^Y)^2 + (shift^Z)^2 )
00456       //  at the center of NS
00457       //----------------------------------------------------------------
00458 
00459       dshift2 = 2.*factx*( (shift(1).val_grid_point(0,0,0,0))
00460                * ((shift_auto(1).dsdx()).val_grid_point(0,0,0,0)
00461                   + dshift_comp(1,1).val_grid_point(0,0,0,0))
00462                +(shift(2).val_grid_point(0,0,0,0))
00463                * ((shift_auto(2).dsdx()).val_grid_point(0,0,0,0)
00464                   + dshift_comp(1,2).val_grid_point(0,0,0,0))
00465                +(shift(3).val_grid_point(0,0,0,0))
00466                * ((shift_auto(3).dsdx()).val_grid_point(0,0,0,0)
00467                   + dshift_comp(1,3).val_grid_point(0,0,0,0)) ) ;
00468 
00469 
00470       //-------------------------
00471       // Information of position
00472       //-------------------------
00473 
00474       x_orb = (star.get_mp()).get_ori_x() - x_rot ;
00475       y_orb = (star.get_mp()).get_ori_y() - y_rot ;
00476 
00477 
00478       // Output
00479       // ------
00480 
00481       cout << "Bin_bhns::orbit_omega: central d(log(lap)+log(Gam))/dX : "
00482        << dnulg << endl ;
00483       cout << "Bin_bhns::orbit_omega: central psi^6/lapconf^2 :         "
00484        << p6sl2 << endl ;
00485       cout << "Bin_bhns::orbit_omega: central d(psi^6/lapconf^2)/dX :   "
00486        << dp6sl2 << endl ;
00487       cout << "Bin_bhns::orbit_omega: central shift^X :                 "
00488        << shiftx << endl ;
00489       cout << "Bin_bhns::orbit_omega: central shift^Y :                 "
00490        << shifty << endl ;
00491       cout << "Bin_bhns::orbit_omega: central d(shift^X)/dX :           "
00492        << dshiftx << endl ;
00493       cout << "Bin_bhns::orbit_omega: central d(shift^Y)/dX :           "
00494        << dshifty << endl ;
00495       cout << "Bin_bhns::orbit_omega: central shift^i shift_i :         "
00496        << shift2 << endl ;
00497       cout << "Bin_bhns::orbit_omega: central d(shift^i shift_i)/dX :   "
00498        << dshift2 << endl ;
00499 
00500 
00501       //---------------------------------------------------------------//
00502       //          Calculation of the orbital angular velocity          //
00503       //---------------------------------------------------------------//
00504 
00505       Param parorb ;
00506       parorb.add_double(dnulg, 0) ;
00507       parorb.add_double(p6sl2, 1) ;
00508       parorb.add_double(dp6sl2, 2) ;
00509       parorb.add_double(shiftx, 3) ;
00510       parorb.add_double(shifty, 4) ;
00511       parorb.add_double(dshiftx, 5) ;
00512       parorb.add_double(dshifty, 6) ;
00513       parorb.add_double(shift2, 7) ;
00514       parorb.add_double(dshift2, 8) ;
00515       parorb.add_double(x_orb, 9) ;
00516       parorb.add_double(y_orb, 10) ;
00517 
00518       double omega1 = fact_omeg_min * omega ;
00519       double omega2 = fact_omeg_max * omega ;
00520 
00521       cout << "Bin_bhns::orbit_omega: omega1,  omega2 [rad/s] : "
00522        << omega1 * f_unit << "  " << omega2 * f_unit << endl ;
00523 
00524       // Search for the various zeros in the interval [omega1,omega2]
00525       // ------------------------------------------------------------
00526 
00527       int nsub = 50 ;  // total number of subdivisions of the interval
00528       Tbl* azer = 0x0 ;
00529       Tbl* bzer = 0x0 ;
00530       zero_list(func_binbhns_orbit_is, parorb, omega1, omega2, nsub,
00531         azer, bzer) ;
00532 
00533       // Search for the zero closest to the previous value of omega
00534       // ----------------------------------------------------------
00535 
00536       double omeg_min, omeg_max ;
00537       int nzer = azer->get_taille() ; // number of zeros found by zero_list
00538 
00539       cout << "Bin_bhns::orbit_omega: " << nzer
00540        << " zero(s) found in the interval [omega1,  omega2]." << endl ;
00541       cout << "omega, omega1, omega2 : " << omega << "  " << omega1
00542        << "  " << omega2 << endl ;
00543       cout << "azer : " << *azer << endl ;
00544       cout << "bzer : " << *bzer << endl ;
00545 
00546       if (nzer == 0) {
00547         cout <<
00548       "Bin_bhns::orbit_omega: WARNING : no zero detected in the interval"
00549          << endl << "   [" << omega1 * f_unit << ", "
00550          << omega2 * f_unit << "]  rad/s  !" << endl ;
00551     omeg_min = omega1 ;
00552     omeg_max = omega2 ;
00553       }
00554       else {
00555         double dist_min = fabs(omega2 - omega1) ;
00556     int i_dist_min = -1 ;
00557     for (int i=0; i<nzer; i++) {
00558       // Distance of previous value of omega from the center of the
00559       //  interval [azer(i), bzer(i)]
00560 
00561       double dist = fabs( omega - 0.5 * ( (*azer)(i) + (*bzer)(i) ) ) ;
00562 
00563       if (dist < dist_min) {
00564         dist_min = dist ;
00565         i_dist_min = i ;
00566       }
00567     }
00568     omeg_min = (*azer)(i_dist_min) ;
00569     omeg_max = (*bzer)(i_dist_min) ;
00570       }
00571 
00572       delete azer ; // Tbl allocated by zero_list
00573       delete bzer ;
00574 
00575       cout << "Bin_bhns::orbit_omega: interval selected for the search of the zero : "
00576        << endl << "  [" << omeg_min << ", " << omeg_max << "]  =  ["
00577        << omeg_min * f_unit << ", " << omeg_max * f_unit << "] rad/s " << endl ;
00578 
00579       // Computation of the zero in the selected interval by the secant method
00580       // ---------------------------------------------------------------------
00581 
00582       int nitermax = 200 ;
00583       int niter ;
00584       double precis = 1.e-13 ;
00585       omega = zerosec_b(func_binbhns_orbit_is, parorb, omeg_min, omeg_max,
00586             precis, nitermax, niter) ;
00587 
00588       cout << "Bin_bhns::orbit_omega: Number of iterations in zerosec for omega : "
00589        << niter << endl ;
00590 
00591       cout << "Bin_bhns::orbit_omega: omega [rad/s] : "
00592        << omega * f_unit << endl ;
00593 
00594     }
00595 }
00596 
00597 //******************************
00598 // Function for searching omega
00599 //******************************
00600 
00601 double func_binbhns_orbit_ks(double om, const Param& parorb) {
00602 
00603     double dnulg = parorb.get_double(0) ;
00604     double p6sl2 = parorb.get_double(1) ;
00605     double dp6sl2 = parorb.get_double(2) ;
00606     double shiftx = parorb.get_double(3) ;
00607     double shifty = parorb.get_double(4) ;
00608     double dshiftx = parorb.get_double(5) ;
00609     double dshifty = parorb.get_double(6) ;
00610     double shift2 = parorb.get_double(7) ;
00611     double dshift2 = parorb.get_double(8) ;
00612     double x_orb = parorb.get_double(9) ;
00613     double y_orb = parorb.get_double(10) ;
00614     double x_separ = parorb.get_double(11) ;
00615     double y_separ = parorb.get_double(12) ;
00616     double xbh_orb = parorb.get_double(13) ;
00617     double mhsr = parorb.get_double(14) ;
00618 
00619     double om2 = om * om ;
00620     /*
00621     double bpb = om2 * (x_orb * x_orb + y_orb * y_orb)
00622       + 2. * om * (shifty * x_orb - shiftx * y_orb) + shift2
00623       + 2. * mhsr * (x_separ*x_separ+y_separ*y_separ)
00624       * (x_separ*shiftx + y_separ*shifty + om * y_separ * xbh_orb)
00625       * (x_separ*shiftx + y_separ*shifty + om * y_separ * xbh_orb) ;
00626       */
00627 
00628     double bpb = om2 * (x_orb * x_orb + y_orb * y_orb
00629             + 2.*mhsr*(x_separ*x_separ+y_separ*y_separ)
00630             * y_separ * y_separ * xbh_orb * xbh_orb)
00631       + 2.*om*(shifty * x_orb - shiftx * y_orb
00632            + 2.*mhsr*(x_separ*x_separ+y_separ*y_separ)
00633            * (shiftx*x_separ + shifty*y_separ) * y_separ * xbh_orb)
00634       + shift2 + 2.*mhsr*(x_separ*x_separ+y_separ*y_separ)
00635       * (shiftx*x_separ + shifty*y_separ)
00636       * (shiftx*x_separ + shifty*y_separ) ;
00637 
00638     double dlngam0 =
00639       ( 0.5 * dp6sl2 * bpb
00640     + p6sl2 * (0.5*dshift2
00641            + om * (shifty - dshiftx*y_orb + dshifty*x_orb)
00642            + om2 * x_orb
00643            - mhsr * x_separ * (x_separ*shiftx+y_separ*shifty)
00644            * (x_separ*shiftx+y_separ*shifty)
00645            + 2. * mhsr * (x_separ*shiftx+y_separ*shifty)
00646            * (y_separ*y_separ*shiftx - x_separ*y_separ*shifty
00647               +(x_separ*x_separ+y_separ*y_separ)*(x_separ*dshiftx
00648                               +y_separ*dshifty))
00649            + 2. * mhsr * om * y_separ * xbh_orb
00650            * ( (y_separ*y_separ-2.*x_separ*x_separ)*shiftx
00651                - 3. * x_separ * y_separ * shifty
00652                +(x_separ*x_separ+y_separ*y_separ)*(x_separ*dshiftx
00653                                +y_separ*dshifty) )
00654            - 3. * mhsr * om2 * x_separ * y_separ * y_separ * xbh_orb
00655            * xbh_orb)
00656     ) / (1 - p6sl2 * bpb) ;
00657 
00658     return dnulg - dlngam0 ;
00659 
00660 }
00661 
00662 double func_binbhns_orbit_is(double om, const Param& parorb) {
00663 
00664     double dnulg = parorb.get_double(0) ;
00665     double p6sl2 = parorb.get_double(1) ;
00666     double dp6sl2 = parorb.get_double(2) ;
00667     double shiftx = parorb.get_double(3) ;
00668     double shifty = parorb.get_double(4) ;
00669     double dshiftx = parorb.get_double(5) ;
00670     double dshifty = parorb.get_double(6) ;
00671     double shift2 = parorb.get_double(7) ;
00672     double dshift2 = parorb.get_double(8) ;
00673     double x_orb = parorb.get_double(9) ;
00674     double y_orb = parorb.get_double(10) ;
00675 
00676     double om2 = om * om ;
00677 
00678     double bpb = om2 * (x_orb * x_orb + y_orb * y_orb)
00679       + 2. * om * (shifty * x_orb - shiftx * y_orb) + shift2 ;
00680 
00681     double dlngam0 = ( 0.5 * dp6sl2 * bpb
00682                + p6sl2 * (0.5*dshift2
00683                   + om *
00684                   (shifty - dshiftx*y_orb + dshifty*x_orb)
00685                   + om2 * x_orb)
00686                ) / (1 - p6sl2 * bpb) ;
00687 
00688     return dnulg - dlngam0 ;
00689 
00690 }

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