binary.C

00001 /*
00002  * Methods of class Binary
00003  *
00004  */
00005 
00006 /*
00007  *   Copyright (c) 2004 Francois Limousin
00008  *
00009  *   This file is part of LORENE.
00010  *
00011  *   LORENE is free software; you can redistribute it and/or modify
00012  *   it under the terms of the GNU General Public License as published by
00013  *   the Free Software Foundation; either version 2 of the License, or
00014  *   (at your option) any later version.
00015  *
00016  *   LORENE is distributed in the hope that it will be useful,
00017  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00018  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00019  *   GNU General Public License for more details.
00020  *
00021  *   You should have received a copy of the GNU General Public License
00022  *   along with LORENE; if not, write to the Free Software
00023  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00024  *
00025  */
00026 char Binary_C[] = "$Header: /cvsroot/Lorene/C++/Source/Binary/binary.C,v 1.14 2005/09/15 14:39:14 e_gourgoulhon Exp $" ;
00027 
00028 /*
00029  * $Id: binary.C,v 1.14 2005/09/15 14:39:14 e_gourgoulhon Exp $
00030  * $Log: binary.C,v $
00031  * Revision 1.14  2005/09/15 14:39:14  e_gourgoulhon
00032  * Added printing of angular momentum in display_poly.
00033  *
00034  * Revision 1.13  2005/09/13 19:38:31  f_limousin
00035  * Reintroduction of the resolution of the equations in cartesian coordinates.
00036  *
00037  * Revision 1.12  2005/02/24 17:31:27  f_limousin
00038  * Update of the function decouple().
00039  *
00040  * Revision 1.11  2005/02/18 13:14:06  j_novak
00041  * Changing of malloc/free to new/delete + suppression of some unused variables
00042  * (trying to avoid compilation warnings).
00043  *
00044  * Revision 1.10  2004/07/21 11:47:10  f_limousin
00045  * Add / Delete comments.
00046  *
00047  * Revision 1.9  2004/05/25 14:25:12  f_limousin
00048  * Add the virial theorem for conformally flat configurations.
00049  *
00050  * Revision 1.8  2004/03/25 10:29:01  j_novak
00051  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00052  *
00053  * Revision 1.7  2004/03/23 10:00:47  f_limousin
00054  * Minor changes.
00055  *
00056  * Revision 1.6  2004/02/27 09:59:33  f_limousin
00057  * Modification in the routine decouple().
00058  *
00059  * Revision 1.5  2004/02/21 17:05:12  e_gourgoulhon
00060  * Method Scalar::point renamed Scalar::val_grid_point.
00061  * Method Scalar::set_point renamed Scalar::set_grid_point.
00062  *
00063  * Revision 1.4  2004/01/20 15:21:12  f_limousin
00064  * First version
00065  *
00066  *
00067  * $Header: /cvsroot/Lorene/C++/Source/Binary/binary.C,v 1.14 2005/09/15 14:39:14 e_gourgoulhon Exp $
00068  *
00069  */
00070 
00071 // Headers C
00072 #include <math.h>
00073 
00074 // Headers Lorene
00075 #include "binary.h"
00076 #include "eos.h"
00077 #include "utilitaires.h"
00078 #include "graphique.h"
00079 #include "param.h"
00080 #include "unites.h"     
00081 
00082                 //--------------//
00083                 // Constructors //
00084                 //--------------//
00085 
00086 // Standard constructor
00087 // --------------------
00088 
00089 Binary::Binary(Map& mp1, int nzet1, const Eos& eos1, int irrot1, 
00090            Map& mp2, int nzet2, const Eos& eos2, int irrot2, 
00091            int conf_flat) 
00092                  : star1(mp1, nzet1, eos1, irrot1, conf_flat), 
00093            star2(mp2, nzet2, eos2, irrot2, conf_flat)
00094 {
00095 
00096     et[0] = &star1 ; 
00097     et[1] = &star2 ; 
00098     
00099     omega = 0 ; 
00100     x_axe = 0 ; 
00101 
00102     // Pointers of derived quantities initialized to zero :
00103     set_der_0x0() ;
00104 }
00105 
00106 // Copy constructor
00107 // ----------------
00108 Binary::Binary(const Binary& bibi) 
00109         : star1(bibi.star1), 
00110           star2(bibi.star2),
00111           omega(bibi.omega), 
00112           x_axe(bibi.x_axe) 
00113 {
00114     et[0] = &star1 ; 
00115     et[1] = &star2 ; 
00116 
00117     // Pointers of derived quantities initialized to zero :
00118     set_der_0x0() ;    
00119 }
00120 
00121 // Constructor from a file
00122 // -----------------------
00123 Binary::Binary(Map& mp1, const Eos& eos1, Map& mp2, const Eos& eos2, 
00124            FILE* fich)
00125         : star1(mp1, eos1, fich), 
00126           star2(mp2, eos2, fich) 
00127 {
00128     et[0] = &star1 ; 
00129     et[1] = &star2 ; 
00130 
00131     // omega and x_axe are read in the file:
00132     fread_be(&omega, sizeof(double), 1, fich) ;     
00133     fread_be(&x_axe, sizeof(double), 1, fich) ;     
00134 
00135     // Pointers of derived quantities initialized to zero : 
00136     set_der_0x0() ;    
00137     
00138 }           
00139 
00140                 //------------//
00141                 // Destructor //
00142                 //------------//
00143 
00144 Binary::~Binary(){
00145 
00146     del_deriv() ; 
00147 
00148 }
00149 
00150             //----------------------------------//
00151             // Management of derived quantities //
00152             //----------------------------------//
00153 
00154 void Binary::del_deriv() const {
00155 
00156     if (p_mass_adm != 0x0) delete p_mass_adm ; 
00157     if (p_mass_kom != 0x0) delete p_mass_kom ; 
00158     if (p_angu_mom != 0x0) delete p_angu_mom ; 
00159     if (p_total_ener != 0x0) delete p_total_ener ; 
00160     if (p_virial != 0x0) delete p_virial ; 
00161     if (p_ham_constr != 0x0) delete p_ham_constr ; 
00162     if (p_mom_constr != 0x0) delete p_mom_constr ; 
00163 
00164     set_der_0x0() ; 
00165 }               
00166 
00167 
00168 
00169 
00170 void Binary::set_der_0x0() const {
00171 
00172     p_mass_adm = 0x0 ; 
00173     p_mass_kom = 0x0 ; 
00174     p_angu_mom = 0x0 ; 
00175     p_total_ener = 0x0 ; 
00176     p_virial = 0x0 ; 
00177     p_ham_constr = 0x0 ; 
00178     p_mom_constr = 0x0 ; 
00179 
00180 }               
00181 
00182 
00183                 //--------------//
00184                 //  Assignment  //
00185                 //--------------//
00186 
00187 // Assignment to another Binary
00188 // --------------------------------
00189 
00190 void Binary::operator=(const Binary& bibi) {
00191 
00192     star1 = bibi.star1 ; 
00193     star2 = bibi.star2 ; 
00194     
00195     omega = bibi.omega ; 
00196     x_axe = bibi.x_axe ; 
00197     
00198     del_deriv() ;  // Deletes all derived quantities
00199     
00200 }
00201 
00202                 //--------------//
00203                 //    Outputs   //
00204                 //--------------//
00205 
00206 // Save in a file
00207 // --------------
00208 void Binary::sauve(FILE* fich) const {
00209     
00210     star1.sauve(fich) ; 
00211     star2.sauve(fich) ; 
00212     
00213     fwrite_be(&omega, sizeof(double), 1, fich) ;        
00214     fwrite_be(&x_axe, sizeof(double), 1, fich) ;        
00215     
00216 }
00217 
00218 // Printing
00219 // --------
00220 ostream& operator<<(ostream& ost, const Binary& bibi)  {
00221     bibi >> ost ;
00222     return ost ;
00223 }
00224     
00225 
00226 ostream& Binary::operator>>(ostream& ost) const {
00227 
00228   using namespace Unites ;
00229 
00230     ost << endl ; 
00231     ost << "Binary neutron stars" << endl ; 
00232     ost << "=============" << endl ; 
00233     ost << endl << 
00234     "Orbital angular velocity : " << omega * f_unit << " rad/s" << endl ; 
00235     ost << endl << 
00236     "Coordinate separation between the two stellar centers : " 
00237     << separation() / km  << " km" << endl ; 
00238     ost << 
00239     "Absolute coordinate X of the rotation axis : " << x_axe / km 
00240         << " km" << endl ; 
00241     ost << endl << "Star 1 : " << endl ; 
00242     ost << "======   " << endl ; 
00243     ost << star1 << endl ; 
00244     ost << "Star 2 : " << endl ; 
00245     ost << "======   " << endl ; 
00246     ost << star2 << endl ; 
00247     return ost ;
00248 }
00249 
00250 // Display in polytropic units
00251 // ---------------------------
00252 
00253 void Binary::display_poly(ostream& ost) const {
00254 
00255   using namespace Unites ;
00256 
00257     const Eos* p_eos1 = &( star1.get_eos() ) ; 
00258     const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;    
00259 
00260     if (p_eos_poly != 0x0) {
00261 
00262     assert( star1.get_eos() == star2.get_eos() ) ; 
00263 
00264     double kappa = p_eos_poly->get_kap() ; 
00265     double gamma = p_eos_poly->get_gam() ;  ; 
00266     double kap_ns2 = pow( kappa,  0.5 /(gamma-1) ) ; 
00267     
00268     // Polytropic unit of length in terms of r_unit : 
00269     double r_poly = kap_ns2 / sqrt(ggrav) ; 
00270     
00271     // Polytropic unit of time in terms of t_unit :
00272     double t_poly = r_poly ; 
00273 
00274     // Polytropic unit of mass in terms of m_unit :
00275     double m_poly = r_poly / ggrav ; 
00276     
00277     // Polytropic unit of angular momentum in terms of j_unit :
00278     double j_poly = r_poly * r_poly / ggrav ; 
00279     
00280     ost.precision(10) ; 
00281     ost << endl << "Quantities in polytropic units : " << endl ; 
00282     ost  << "==============================" << endl ; 
00283     ost << " ( r_poly = " << r_poly / km << " km )" << endl ; 
00284     ost << "  d_e_max   : " << separation() / r_poly << endl ; 
00285     ost << "  d_G       : " 
00286          << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly 
00287          << endl ; 
00288     ost << "  Omega   : " << omega * t_poly << endl ; 
00289     ost << "  J   : " << angu_mom()(2) / j_poly << endl ; 
00290     ost << "  M_ADM   : " << mass_adm() / m_poly << endl ;      
00291     ost << "  M_Komar : " << mass_kom() / m_poly << endl ; 
00292     ost << "  M_bar(star 1) : " << star1.mass_b() / m_poly << endl ; 
00293     ost << "  M_bar(star 2) : " << star2.mass_b() / m_poly << endl ; 
00294     ost << "  R_0(star 1)   : " << 
00295     0.5 * ( star1.ray_eq() + star1.ray_eq_pi() ) / r_poly << endl ;  
00296     ost << "  R_0(star 2)   : " << 
00297     0.5 * ( star2.ray_eq() + star2.ray_eq_pi() ) / r_poly << endl ;  
00298     
00299     }
00300     
00301 
00302 } 
00303 
00304 
00305 void Binary::fait_decouple () {
00306 
00307     int nz_un = star1.get_mp().get_mg()->get_nzone() ;
00308     int nz_deux = star2.get_mp().get_mg()->get_nzone() ;
00309     
00310     // We determine R_limite :
00311     double distance = star2.get_mp().get_ori_x() - star1.get_mp().get_ori_x() ;
00312     cout << "distance = " << distance << endl ;
00313     double lim_un = distance/2. ;
00314     double lim_deux = distance/2. ;
00315     double int_un = distance/6. ;
00316     double int_deux = distance/6. ;
00317     
00318     // The functions used.
00319     Scalar fonction_f_un (star1.get_mp()) ;
00320     fonction_f_un = 0.5*pow(
00321     cos((star1.get_mp().r-int_un)*M_PI/2./(lim_un-int_un)), 2.)+0.5 ;
00322     fonction_f_un.std_spectral_base();
00323     
00324     Scalar fonction_g_un (star1.get_mp()) ;
00325     fonction_g_un = 0.5*pow
00326     (sin((star1.get_mp().r-int_un)*M_PI/2./(lim_un-int_un)), 2.) ;
00327     fonction_g_un.std_spectral_base();
00328     
00329     Scalar fonction_f_deux (star2.get_mp()) ;
00330     fonction_f_deux = 0.5*pow(
00331     cos((star2.get_mp().r-int_deux)*M_PI/2./(lim_deux-int_deux)), 2.)+0.5 ;
00332     fonction_f_deux.std_spectral_base();
00333     
00334     Scalar fonction_g_deux (star2.get_mp()) ;
00335     fonction_g_deux = 0.5*pow(
00336     sin((star2.get_mp().r-int_deux)*M_PI/2./(lim_un-int_deux)), 2.) ;
00337     fonction_g_deux.std_spectral_base();
00338     
00339     // The functions total :
00340     Scalar decouple_un (star1.get_mp()) ;
00341     decouple_un.allocate_all() ;
00342     Scalar decouple_deux (star2.get_mp()) ;
00343     decouple_deux.allocate_all() ;
00344     
00345     Mtbl xabs_un (star1.get_mp().xa) ;
00346     Mtbl yabs_un (star1.get_mp().ya) ;
00347     Mtbl zabs_un (star1.get_mp().za) ;
00348         
00349     Mtbl xabs_deux (star2.get_mp().xa) ;
00350     Mtbl yabs_deux (star2.get_mp().ya) ;
00351     Mtbl zabs_deux (star2.get_mp().za) ;
00352         
00353     double xabs, yabs, zabs, air_un, air_deux, theta, phi ;
00354         
00355     for (int l=0 ; l<nz_un ; l++) {
00356     int nr = star1.get_mp().get_mg()->get_nr (l) ;
00357         
00358     if (l==nz_un-1)
00359         nr -- ;
00360         
00361     int np = star1.get_mp().get_mg()->get_np (l) ;
00362     int nt = star1.get_mp().get_mg()->get_nt (l) ;
00363         
00364     for (int k=0 ; k<np ; k++)
00365         for (int j=0 ; j<nt ; j++)
00366         for (int i=0 ; i<nr ; i++) {
00367                 
00368             xabs = xabs_un (l, k, j, i) ;
00369             yabs = yabs_un (l, k, j, i) ;
00370             zabs = zabs_un (l, k, j, i) ;
00371                 
00372             // Coordinates of the point
00373             star1.get_mp().convert_absolute 
00374             (xabs, yabs, zabs, air_un, theta, phi) ;
00375             star2.get_mp().convert_absolute 
00376             (xabs, yabs, zabs, air_deux, theta, phi) ;
00377         
00378             if (air_un <= lim_un)
00379             if (air_un < int_un)
00380                 decouple_un.set_grid_point(l, k, j, i) = 1 ;
00381             else
00382             // Close to star 1 :
00383             decouple_un.set_grid_point(l, k, j, i) = 
00384                 fonction_f_un.val_grid_point(l, k, j, i) ;
00385             else 
00386             if (air_deux <= lim_deux)
00387                 if (air_deux < int_deux)
00388                 decouple_un.set_grid_point(l, k, j, i) = 0 ;
00389                 else
00390             // Close to star 2 :
00391                 decouple_un.set_grid_point(l, k, j, i) = 
00392         fonction_g_deux.val_point (air_deux, theta, phi) ;
00393         
00394             else
00395                 // Far from each star :
00396                 decouple_un.set_grid_point(l, k, j, i) = 0.5 ;
00397         }
00398                 
00399         // Case infinity :
00400         if (l==nz_un-1)
00401             for (int k=0 ; k<np ; k++)
00402             for (int j=0 ; j<nt ; j++)
00403                 decouple_un.set_grid_point(nz_un-1, k, j, nr)=0.5 ;
00404         }
00405     
00406     for (int l=0 ; l<nz_deux ; l++) {
00407     int nr = star2.get_mp().get_mg()->get_nr (l) ;
00408         
00409     if (l==nz_deux-1)
00410         nr -- ;
00411         
00412     int np = star2.get_mp().get_mg()->get_np (l) ;
00413     int nt = star2.get_mp().get_mg()->get_nt (l) ;
00414         
00415     for (int k=0 ; k<np ; k++)
00416         for (int j=0 ; j<nt ; j++)
00417         for (int i=0 ; i<nr ; i++) {
00418                 
00419             xabs = xabs_deux (l, k, j, i) ;
00420             yabs = yabs_deux (l, k, j, i) ;
00421             zabs = zabs_deux (l, k, j, i) ;
00422                 
00423             // coordinates of the point  :
00424             star1.get_mp().convert_absolute 
00425             (xabs, yabs, zabs, air_un, theta, phi) ;
00426             star2.get_mp().convert_absolute 
00427             (xabs, yabs, zabs, air_deux, theta, phi) ;
00428             
00429             if (air_deux <= lim_deux)
00430             if (air_deux < int_deux)
00431                 decouple_deux.set_grid_point(l, k, j, i) = 1 ;
00432             else
00433             // close to star two :
00434             decouple_deux.set_grid_point(l, k, j, i) = 
00435                 fonction_f_deux.val_grid_point(l, k, j, i) ;
00436             else 
00437             if (air_un <= lim_un)
00438                 if (air_un < int_un)
00439                 decouple_deux.set_grid_point(l, k, j, i) = 0 ;
00440                 else
00441             // close to star one :
00442                 decouple_deux.set_grid_point(l, k, j, i) = 
00443              fonction_g_un.val_point (air_un, theta, phi) ;
00444         
00445             else
00446                 // Far from each star :
00447                 decouple_deux.set_grid_point(l, k, j, i) = 0.5 ;
00448         }
00449                 
00450         // Case infinity :
00451         if (l==nz_deux-1)
00452             for (int k=0 ; k<np ; k++)
00453             for (int j=0 ; j<nt ; j++)
00454              decouple_deux.set_grid_point(nz_un-1, k, j, nr)=0.5 ;
00455    }
00456    
00457     decouple_un.std_spectral_base() ;
00458     decouple_deux.std_spectral_base() ;
00459 
00460     int nr = star1.get_mp().get_mg()->get_nr (1) ;
00461     int nt = star1.get_mp().get_mg()->get_nt (1) ;
00462     int np = star1.get_mp().get_mg()->get_np (1) ;
00463     cout << "decouple_un"  << endl << norme(decouple_un/(nr*nt*np)) << endl ;
00464     cout << "decouple_deux"  << endl << norme(decouple_deux/(nr*nt*np)) 
00465      << endl ;
00466 
00467     star1.decouple = decouple_un ;
00468     star2.decouple = decouple_deux ;
00469 }
00470 
00471 void Binary::write_global(ostream& ost) const {
00472 
00473   using namespace Unites ;
00474 
00475   const Map& mp1 = star1.get_mp() ;
00476   const Mg3d* mg1 = mp1.get_mg() ;
00477   int nz1 = mg1->get_nzone() ; 
00478 
00479   ost.precision(5) ;
00480   ost << "# Grid 1 : " << nz1 << "x"
00481       << mg1->get_nr(0) << "x" << mg1->get_nt(0) << "x" << mg1->get_np(0) 
00482       << "  R_out(l) [km] : " ;
00483   for (int l=0; l<nz1; l++) {
00484     ost << " " << mp1.val_r(l, 1., M_PI/2, 0) / km ; 
00485   }
00486   ost << endl ; 
00487 
00488   ost << "#     VE(M)      " << endl ;
00489   
00490   
00491   ost.setf(ios::scientific) ; 
00492   ost.width(14) ; 
00493   ost << virial() << endl ;
00494   
00495   ost << "#      d [km]         "  
00496       << "       d_G [km]       "
00497       << "     d/(a1 +a1')      "
00498       << "       f [Hz]         "
00499       << "    M_ADM [M_sol]     "     
00500       << "    M_ADM_vol [M_sol]     "     
00501       << "    M_Komar [M_sol]     "     
00502       << "    M_Komar_vol [M_sol]     "     
00503       << "   J [G M_sol^2/c]    "  << endl ;   
00504   
00505   ost.precision(14) ;
00506   ost.width(20) ; 
00507   ost << separation() / km ; ost.width(22) ;
00508   ost   << ( star2.xa_barycenter() - star1.xa_barycenter() ) / km ; ost.width(22) ;
00509   ost   << separation() / (star1.ray_eq() + star2.ray_eq()) ; ost.width(22) ;
00510   ost   << omega / (2*M_PI)* f_unit ; ost.width(22) ;
00511   ost   << mass_adm() / msol ; ost.width(22) ; 
00512   ost   << mass_adm_vol() / msol ; ost.width(22) ; 
00513   ost   << mass_kom() / msol ; ost.width(22) ; 
00514   ost   << mass_kom_vol() / msol ; ost.width(22) ; 
00515   ost   << angu_mom()(2)/ ( qpig / (4* M_PI) * msol*msol) << endl ; 
00516   
00517   ost << "#     H_c(1)[c^2]     "
00518       << "    e_c(1)[rho_nuc]   " 
00519       << "    M_B(1) [M_sol]    "
00520       << "     r_eq(1) [km]     "
00521       << "        a2/a1(1)    " 
00522         << "        a3/a1(1)      " << endl ; 
00523   
00524   ost.width(20) ; 
00525   ost << star1.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
00526   ost   << star1.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
00527   ost   << star1.mass_b() / msol ; ost.width(22) ;  
00528   ost << star1.ray_eq() / km ; ost.width(22) ; 
00529   ost   << star1.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
00530   ost   << star1.ray_pole() / star1.ray_eq() << endl ;
00531   
00532   ost << "#     H_c(2)[c^2]     "
00533       << "    e_c(2)[rho_nuc]   " 
00534       << "    M_B(2) [M_sol]    "
00535       << "     r_eq(2) [km]     "
00536       << "        a2/a1(2)    " 
00537       << "        a3/a1(2)    " << endl ; 
00538   
00539   ost.width(20) ; 
00540   ost << star2.get_ent().val_grid_point(0,0,0,0) ; ost.width(22) ;
00541   ost   << star2.get_ener().val_grid_point(0,0,0,0) ; ost.width(22) ;
00542   ost   << star2.mass_b() / msol ; ost.width(22) ;  
00543   ost << star2.ray_eq() / km ; ost.width(22) ; 
00544   ost   << star2.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
00545   ost   << star2.ray_pole() / star1.ray_eq() << endl ;
00546   
00547   // Quantities in polytropic units if the EOS is a polytropic one
00548   // -------------------------------------------------------------
00549   const Eos* p_eos1 = &( star1.get_eos() ) ; 
00550   const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;      
00551   
00552   if ((p_eos_poly != 0x0) && ( star1.get_eos() == star2.get_eos() )) {
00553       
00554       double kappa = p_eos_poly->get_kap() ; 
00555       double gamma = p_eos_poly->get_gam() ;  ; 
00556       double kap_ns2 = pow( kappa,  0.5 /(gamma-1.) ) ; 
00557       
00558       // Polytropic unit of length in terms of r_unit : 
00559       double r_poly = kap_ns2 / sqrt(ggrav) ; 
00560       
00561       // Polytropic unit of time in terms of t_unit :
00562       double t_poly = r_poly ; 
00563       
00564       // Polytropic unit of mass in terms of m_unit :
00565       double m_poly = r_poly / ggrav ; 
00566       
00567       // Polytropic unit of angular momentum in terms of j_unit :
00568       double j_poly = r_poly * r_poly / ggrav ; 
00569       
00570       ost << "#      d [poly]       "  
00571       << "       d_G [poly]     "
00572       << "     Omega [poly]     "
00573       << "     M_ADM [poly]     "     
00574       << "       J [poly]       "  
00575       << "    M_B(1) [poly]     "
00576       << "    M_B(2) [poly]     " << endl ; 
00577       
00578       ost.width(20) ; 
00579       ost << separation() / r_poly ; ost.width(22) ;
00580       ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly ; ost.width(22) ; 
00581       ost << omega * t_poly ; ost.width(22) ;
00582       ost << mass_adm() / m_poly ; ost.width(22) ;
00583       ost << angu_mom()(2) / j_poly ; ost.width(22) ;
00584       ost << star1.mass_b() / m_poly ; ost.width(22) ;
00585       ost << star2.mass_b() / m_poly << endl ; 
00586       
00587   }
00588   
00589 }
00590 
00591 
00592    
00593             //-------------------------------//
00594             //      Miscellaneous        //
00595             //-------------------------------//
00596 
00597 double Binary::separation() const {
00598     
00599     double dx = star1.mp.get_ori_x() - star2.mp.get_ori_x() ; 
00600     double dy = star1.mp.get_ori_y() - star2.mp.get_ori_y() ; 
00601     double dz = star1.mp.get_ori_z() - star2.mp.get_ori_z() ; 
00602     
00603     return sqrt( dx*dx + dy*dy + dz*dz ) ; 
00604     
00605 }

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