binaire.C

00001 /*
00002  * Methods of class Binaire
00003  *
00004  * (see file binaire.h for documentation)
00005  */
00006 
00007 /*
00008  *   Copyright (c) 2000-2003 Eric Gourgoulhon
00009  *
00010  *   This file is part of LORENE.
00011  *
00012  *   LORENE is free software; you can redistribute it and/or modify
00013  *   it under the terms of the GNU General Public License as published by
00014  *   the Free Software Foundation; either version 2 of the License, or
00015  *   (at your option) any later version.
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 
00029 char binaire_bin_C[] = "$Header: /cvsroot/Lorene/C++/Source/Binaire/binaire.C,v 1.7 2004/03/25 10:28:59 j_novak Exp $" ;
00030 
00031 /*
00032  * $Id: binaire.C,v 1.7 2004/03/25 10:28:59 j_novak Exp $
00033  * $Log: binaire.C,v $
00034  * Revision 1.7  2004/03/25 10:28:59  j_novak
00035  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00036  *
00037  * Revision 1.6  2003/09/15 20:24:24  e_gourgoulhon
00038  * Improvements in the function write_global.
00039  *
00040  * Revision 1.5  2003/09/15 15:10:12  e_gourgoulhon
00041  * Added the member function write_global.
00042  *
00043  * Revision 1.4  2002/12/17 21:18:46  e_gourgoulhon
00044  * Suppression of Etoile_bin::set_companion.
00045  *
00046  * Revision 1.3  2001/12/20 13:03:24  k_taniguchi
00047  * Addition of the Komar mass, the virial error by Gourgoulhon and Bonazzola, and the virial error by Friedman, Uryu, and Shibata.
00048  *
00049  * Revision 1.2  2001/12/04 21:27:52  e_gourgoulhon
00050  *
00051  * All writing/reading to a binary file are now performed according to
00052  * the big endian convention, whatever the system is big endian or
00053  * small endian, thanks to the functions fwrite_be and fread_be
00054  *
00055  * Revision 1.1.1.1  2001/11/20 15:19:30  e_gourgoulhon
00056  * LORENE
00057  *
00058  * Revision 2.7  2001/06/25  12:55:46  eric
00059  * Traitement des etoiles compagnons (appel de Etoile_bin::set_companion dans
00060  * les constructeurs).
00061  *
00062  * Revision 2.6  2000/07/07  14:10:17  eric
00063  * AJout de display_poly.
00064  *
00065  * Revision 2.5  2000/03/13  14:25:52  eric
00066  *  Ajout des membres p_ham_constr et p_mom_constr.
00067  *
00068  * Revision 2.4  2000/02/18  14:52:44  eric
00069  * Ajout des membres p_virial et p_total_ener.
00070  * p_masse_adm --> p_mass_adm.
00071  *
00072  * Revision 2.3  2000/02/12  17:36:25  eric
00073  * Ajout de la fonction separation().
00074  *
00075  * Revision 2.2  2000/02/04  17:14:47  eric
00076  * Ajout du membre ref_triad.
00077  *
00078  * Revision 2.1  2000/02/02  10:40:36  eric
00079  * Modif affichage.
00080  *
00081  * Revision 2.0  2000/01/31  15:57:49  eric
00082  * *** empty log message ***
00083  *
00084  *
00085  * $Header: /cvsroot/Lorene/C++/Source/Binaire/binaire.C,v 1.7 2004/03/25 10:28:59 j_novak Exp $
00086  *
00087  */
00088 
00089 // Headers C
00090 #include "math.h"
00091 
00092 // Headers Lorene
00093 #include "binaire.h"
00094 #include "eos.h"
00095 #include "utilitaires.h"
00096 #include "unites.h"     
00097 
00098                 //--------------//
00099                 // Constructors //
00100                 //--------------//
00101 
00102 // Standard constructor
00103 // --------------------
00104 
00105 Binaire::Binaire(Map& mp1, int nzet1, const Eos& eos1, int irrot1, 
00106          Map& mp2, int nzet2, const Eos& eos2, int irrot2,
00107          int relat) 
00108          : ref_triad(0., "Absolute frame Cartesian basis"),  
00109            star1(mp1, nzet1, relat, eos1, irrot1, ref_triad), 
00110            star2(mp2, nzet2, relat, eos2, irrot2, ref_triad)
00111 {
00112 
00113     et[0] = &star1 ; 
00114     et[1] = &star2 ; 
00115     
00116     omega = 0 ; 
00117     x_axe = 0 ; 
00118 
00119     // Pointers of derived quantities initialized to zero :
00120     set_der_0x0() ;
00121 }
00122 
00123 // Copy constructor
00124 // ----------------
00125 Binaire::Binaire(const Binaire& bibi) 
00126         : ref_triad(0., "Absolute frame Cartesian basis"), 
00127           star1(bibi.star1), 
00128           star2(bibi.star2),
00129           omega(bibi.omega), 
00130           x_axe(bibi.x_axe) 
00131 {
00132     et[0] = &star1 ; 
00133     et[1] = &star2 ; 
00134 
00135     // Pointers of derived quantities initialized to zero :
00136     set_der_0x0() ;    
00137 }
00138 
00139 // Constructor from a file
00140 // -----------------------
00141 Binaire::Binaire(Map& mp1, const Eos& eos1, Map& mp2, const Eos& eos2, 
00142          FILE* fich)
00143         : ref_triad(0., "Absolute frame Cartesian basis"), 
00144           star1(mp1, eos1, ref_triad, fich), 
00145           star2(mp2, eos2, ref_triad, fich) 
00146 {
00147     et[0] = &star1 ; 
00148     et[1] = &star2 ; 
00149 
00150     // omega and x_axe are read in the file:
00151     fread_be(&omega, sizeof(double), 1, fich) ;     
00152     fread_be(&x_axe, sizeof(double), 1, fich) ;     
00153 
00154     // Pointers of derived quantities initialized to zero : 
00155     set_der_0x0() ;    
00156     
00157 }           
00158 
00159                 //------------//
00160                 // Destructor //
00161                 //------------//
00162 
00163 Binaire::~Binaire(){
00164 
00165     del_deriv() ; 
00166 
00167 }
00168 
00169             //----------------------------------//
00170             // Management of derived quantities //
00171             //----------------------------------//
00172 
00173 void Binaire::del_deriv() const {
00174 
00175     if (p_mass_adm != 0x0) delete p_mass_adm ; 
00176     if (p_mass_kom != 0x0) delete p_mass_kom ; 
00177     if (p_angu_mom != 0x0) delete p_angu_mom ; 
00178     if (p_total_ener != 0x0) delete p_total_ener ; 
00179     if (p_virial != 0x0) delete p_virial ; 
00180     if (p_virial_gb != 0x0) delete p_virial_gb ; 
00181     if (p_virial_fus != 0x0) delete p_virial_fus ; 
00182     if (p_ham_constr != 0x0) delete p_ham_constr ; 
00183     if (p_mom_constr != 0x0) delete p_mom_constr ; 
00184 
00185     set_der_0x0() ; 
00186 }               
00187 
00188 
00189 
00190 
00191 void Binaire::set_der_0x0() const {
00192 
00193     p_mass_adm = 0x0 ; 
00194     p_mass_kom = 0x0 ; 
00195     p_angu_mom = 0x0 ; 
00196     p_total_ener = 0x0 ; 
00197     p_virial = 0x0 ; 
00198     p_virial_gb = 0x0 ; 
00199     p_virial_fus = 0x0 ; 
00200     p_ham_constr = 0x0 ; 
00201     p_mom_constr = 0x0 ; 
00202 
00203 }               
00204 
00205 
00206                 //--------------//
00207                 //  Assignment  //
00208                 //--------------//
00209 
00210 // Assignment to another Binaire
00211 // -----------------------------
00212 
00213 void Binaire::operator=(const Binaire& bibi) {
00214 
00215     assert( bibi.ref_triad == ref_triad ) ; 
00216     
00217     star1 = bibi.star1 ; 
00218     star2 = bibi.star2 ; 
00219     
00220     omega = bibi.omega ; 
00221     x_axe = bibi.x_axe ; 
00222     
00223     // ref_triad remains unchanged 
00224     
00225     del_deriv() ;  // Deletes all derived quantities
00226     
00227 }
00228 
00229                 //--------------//
00230                 //    Outputs   //
00231                 //--------------//
00232 
00233 // Save in a file
00234 // --------------
00235 void Binaire::sauve(FILE* fich) const {
00236     
00237     star1.sauve(fich) ; 
00238     star2.sauve(fich) ; 
00239     
00240     fwrite_be(&omega, sizeof(double), 1, fich) ;        
00241     fwrite_be(&x_axe, sizeof(double), 1, fich) ;        
00242     
00243 }
00244 
00245 // Printing
00246 // --------
00247 ostream& operator<<(ostream& ost, const Binaire& bibi)  {
00248     bibi >> ost ;
00249     return ost ;
00250 }
00251     
00252 
00253 ostream& Binaire::operator>>(ostream& ost) const {
00254 
00255   using namespace Unites ; 
00256 
00257     ost << endl ; 
00258     ost << "Binary system" << endl ; 
00259     ost << "=============" << endl ; 
00260     ost << endl << 
00261     "Orbital angular velocity : " << omega * f_unit << " rad/s" << endl ; 
00262     ost << endl << 
00263     "Coordinate separation between the two stellar centers : " 
00264     << separation() / km  << " km" << endl ; 
00265     ost << 
00266     "Absolute coordinate X of the rotation axis : " << x_axe / km 
00267         << " km" << endl ; 
00268     ost << endl << "Star 1 : " << endl ; 
00269     ost << "======   " << endl ; 
00270     ost << star1 << endl ; 
00271     ost << "Star 2 : " << endl ; 
00272     ost << "======   " << endl ; 
00273     ost << star2 << endl ; 
00274     return ost ;
00275 }
00276 
00277 // Display in polytropic units
00278 // ---------------------------
00279 
00280 void Binaire::display_poly(ostream& ost) const {
00281 
00282   using namespace Unites ; 
00283 
00284     const Eos* p_eos1 = &( star1.get_eos() ) ; 
00285     const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;    
00286 
00287     if ((p_eos_poly != 0x0) && ( star1.get_eos() == star2.get_eos() )) {
00288 
00289     double kappa = p_eos_poly->get_kap() ; 
00290     double gamma = p_eos_poly->get_gam() ;  ; 
00291     double kap_ns2 = pow( kappa,  0.5 /(gamma-1) ) ; 
00292     
00293     // Polytropic unit of length in terms of r_unit : 
00294     double r_poly = kap_ns2 / sqrt(ggrav) ; 
00295     
00296     // Polytropic unit of time in terms of t_unit :
00297     double t_poly = r_poly ; 
00298 
00299     // Polytropic unit of mass in terms of m_unit :
00300     double m_poly = r_poly / ggrav ; 
00301     
00302     // Polytropic unit of angular momentum in terms of j_unit :
00303     double j_poly = r_poly * r_poly / ggrav ; 
00304     
00305     ost.precision(10) ; 
00306     ost << endl << "Quantities in polytropic units : " << endl ; 
00307     ost  << "==============================" << endl ; 
00308     ost << " ( r_poly = " << r_poly / km << " km )" << endl ; 
00309     ost << "  d_e_max   : " << separation() / r_poly << endl ; 
00310     ost << "  d_G       : " 
00311          << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly 
00312          << endl ; 
00313     ost << "  Omega   : " << omega * t_poly << endl ; 
00314     ost << "  J   : " << angu_mom()(2) / j_poly << endl ; 
00315     ost << "  M_ADM   : " << mass_adm() / m_poly << endl ; 
00316     ost << "  M_Komar : " << mass_kom() / m_poly << endl ; 
00317     ost << "  E   : " << total_ener() / m_poly << endl ; 
00318     ost << "  M_bar(star 1) : " << star1.mass_b() / m_poly << endl ; 
00319     ost << "  M_bar(star 2) : " << star2.mass_b() / m_poly << endl ; 
00320     ost << "  R_0(star 1)   : " << 
00321     0.5 * ( star1.ray_eq() + star1.ray_eq_pi() ) / r_poly << endl ;  
00322     ost << "  R_0(star 2)   : " << 
00323     0.5 * ( star2.ray_eq() + star2.ray_eq_pi() ) / r_poly << endl ;  
00324     
00325     }
00326     
00327 
00328 } 
00329 
00330 void Binaire::write_global(ostream& ost) const {
00331 
00332   using namespace Unites ; 
00333 
00334     const Map& mp1 = star1.get_mp() ;
00335     const Mg3d* mg1 = mp1.get_mg() ;
00336     int nz1 = mg1->get_nzone() ; 
00337 
00338     ost.precision(5) ;
00339     ost << "# Grid 1 : " << nz1 << "x"
00340         << mg1->get_nr(0) << "x" << mg1->get_nt(0) << "x" << mg1->get_np(0) 
00341         << "  R_out(l) [km] : " ;
00342     for (int l=0; l<nz1; l++) {
00343         ost << " " << mp1.val_r(l, 1., M_PI/2, 0) / km ; 
00344     }
00345     ost << endl ; 
00346     
00347     ost << "#     VE(M)    "
00348         << "      VE(GB)   "
00349         << "     VE(FUS)   " << endl ; 
00350         
00351     ost.setf(ios::scientific) ; 
00352     ost.width(14) ; 
00353     ost << virial() ; ost.width(14) ;
00354     ost << virial_gb() ; ost.width(14) ; 
00355     ost << virial_fus() << endl ; 
00356 
00357     ost << "#      d [km]         "  
00358         << "       d_G [km]       "
00359         << "     d/(a1 +a1')      "
00360         << "       f [Hz]         "
00361         << "    M_ADM [M_sol]     "     
00362         << "   J [G M_sol^2/c]    "  << endl ;   
00363 
00364     ost.precision(14) ;
00365     ost.width(20) ; 
00366     ost << separation() / km ; ost.width(22) ;
00367     ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / km ; ost.width(22) ;
00368     ost << separation() / (star1.ray_eq() + star2.ray_eq()) ; ost.width(22) ;
00369     ost << omega / (2*M_PI)* f_unit ; ost.width(22) ;
00370     ost << mass_adm() / msol ; ost.width(22) ;
00371     ost << angu_mom()(2)/ ( qpig / (4* M_PI) * msol*msol) << endl ; 
00372                 
00373     ost << "#     H_c(1)[c^2]     "
00374         << "    e_c(1)[rho_nuc]   " 
00375         << "    M_B(1) [M_sol]    "
00376         << "     r_eq(1) [km]     "
00377         << "        a2/a1(1)      " 
00378         << "        a3/a1(1)      " << endl ; 
00379         
00380     ost.width(20) ; 
00381     ost << star1.get_ent()()(0,0,0,0) ; ost.width(22) ;
00382     ost << star1.get_ener()()(0,0,0,0) ; ost.width(22) ;
00383     ost << star1.mass_b() / msol ; ost.width(22) ;  
00384     ost << star1.ray_eq() / km ; ost.width(22) ; 
00385     ost << star1.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
00386     ost << star1.ray_pole() / star1.ray_eq() << endl ;
00387         
00388     ost << "#     H_c(2)[c^2]     "
00389         << "    e_c(2)[rho_nuc]   " 
00390         << "    M_B(2) [M_sol]    "
00391         << "     r_eq(2) [km]     "
00392         << "        a2/a1(2)      " 
00393         << "        a3/a1(2)      " << endl ; 
00394         
00395     ost.width(20) ; 
00396     ost << star2.get_ent()()(0,0,0,0) ; ost.width(22) ;
00397     ost << star2.get_ener()()(0,0,0,0) ; ost.width(22) ;
00398     ost << star2.mass_b() / msol ; ost.width(22) ;  
00399     ost << star2.ray_eq() / km ; ost.width(22) ; 
00400     ost << star2.ray_eq_pis2() / star1.ray_eq() ; ost.width(22) ;
00401     ost << star2.ray_pole() / star1.ray_eq() << endl ;
00402     
00403     // Quantities in polytropic units if the EOS is a polytropic one
00404     // -------------------------------------------------------------
00405     const Eos* p_eos1 = &( star1.get_eos() ) ; 
00406     const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( p_eos1 ) ;    
00407 
00408     if ((p_eos_poly != 0x0) && ( star1.get_eos() == star2.get_eos() )) {
00409 
00410         double kappa = p_eos_poly->get_kap() ; 
00411         double gamma = p_eos_poly->get_gam() ;  ; 
00412         double kap_ns2 = pow( kappa,  0.5 /(gamma-1.) ) ; 
00413     
00414         // Polytropic unit of length in terms of r_unit : 
00415         double r_poly = kap_ns2 / sqrt(ggrav) ; 
00416     
00417         // Polytropic unit of time in terms of t_unit :
00418         double t_poly = r_poly ; 
00419 
00420         // Polytropic unit of mass in terms of m_unit :
00421         double m_poly = r_poly / ggrav ; 
00422     
00423         // Polytropic unit of angular momentum in terms of j_unit :
00424         double j_poly = r_poly * r_poly / ggrav ; 
00425     
00426         ost << "#      d [poly]       "  
00427             << "       d_G [poly]     "
00428             << "     Omega [poly]     "
00429             << "     M_ADM [poly]     "     
00430             << "       J [poly]       "  
00431             << "    M_B(1) [poly]     "
00432             << "    M_B(2) [poly]     " << endl ; 
00433         
00434         ost.width(20) ; 
00435         ost << separation() / r_poly ; ost.width(22) ;
00436         ost << ( star2.xa_barycenter() - star1.xa_barycenter() ) / r_poly ; ost.width(22) ; 
00437         ost << omega * t_poly ; ost.width(22) ;
00438         ost << mass_adm() / m_poly ; ost.width(22) ;
00439         ost << angu_mom()(2) / j_poly ; ost.width(22) ;
00440         ost << star1.mass_b() / m_poly ; ost.width(22) ;
00441         ost << star2.mass_b() / m_poly << endl ; 
00442 
00443     }
00444 
00445 }
00446 
00447 
00448    
00449             //-------------------------------//
00450             //      Miscellaneous        //
00451             //-------------------------------//
00452 
00453 double Binaire::separation() const {
00454     
00455     double dx = star1.get_mp().get_ori_x() - star2.get_mp().get_ori_x() ; 
00456     double dy = star1.get_mp().get_ori_y() - star2.get_mp().get_ori_y() ; 
00457     double dz = star1.get_mp().get_ori_z() - star2.get_mp().get_ori_z() ; 
00458     
00459     return sqrt( dx*dx + dy*dy + dz*dz ) ; 
00460     
00461 }

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