star.C

00001 /*
00002  * Methods of class Star
00003  *
00004  * (see file star.h for documentation)
00005  */
00006 
00007 /*
00008  *   Copyright (c) 2004 Francois Limousin
00009  *  
00010  *   Copyright (c) 2000-2001 Eric Gourgoulhon (for preceding class Etoile)
00011  *   Copyright (c) 2000-2001 Keisuke Taniguchi (for preceding class Etoile)
00012  *
00013  *   This file is part of LORENE.
00014  *
00015  *   LORENE is free software; you can redistribute it and/or modify
00016  *   it under the terms of the GNU General Public License as published by
00017  *   the Free Software Foundation; either version 2 of the License, or
00018  *   (at your option) any later version.
00019  *
00020  *   LORENE is distributed in the hope that it will be useful,
00021  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00022  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00023  *   GNU General Public License for more details.
00024  *
00025  *   You should have received a copy of the GNU General Public License
00026  *   along with LORENE; if not, write to the Free Software
00027  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00028  *
00029  */
00030 
00031 
00032 char star_C[] = "$Header: /cvsroot/Lorene/C++/Source/Star/star.C,v 1.19 2010/02/02 12:45:16 e_gourgoulhon Exp $" ;
00033 
00034 /*
00035  * $Id: star.C,v 1.19 2010/02/02 12:45:16 e_gourgoulhon Exp $
00036  * $Log: star.C,v $
00037  * Revision 1.19  2010/02/02 12:45:16  e_gourgoulhon
00038  * Improved the display (operator>>)
00039  *
00040  * Revision 1.18  2010/01/26 16:49:03  e_gourgoulhon
00041  * Commented the test on the relativistic character of the EOS: the
00042  * relativity parameter is not defined (yet !) in the base class Star.
00043  *
00044  * Revision 1.17  2007/11/06 16:22:03  j_novak
00045  * The data member stress_euler is now a Sym_tensor instead of a Tensor.
00046  *
00047  * Revision 1.16  2007/06/21 19:53:47  k_taniguchi
00048  * Addition of p_ray_eq_3pis2
00049  *
00050  * Revision 1.15  2005/09/14 12:30:52  f_limousin
00051  * Saving of fields lnq and logn in class Star.
00052  *
00053  * Revision 1.14  2005/09/13 19:38:31  f_limousin
00054  * Reintroduction of the resolution of the equations in cartesian coordinates.
00055  *
00056  * Revision 1.13  2005/02/17 17:29:04  f_limousin
00057  * Change the name of some quantities to be consistent with other classes
00058  * (for instance nnn is changed to nn, shift to beta, beta to lnq...)
00059  *
00060  * Revision 1.12  2005/01/05 17:43:03  f_limousin
00061  * u_euler is now constructed in the spherical triad to be consistent
00062  * with all the others vectors ans tensors.
00063  *
00064  * Revision 1.11  2004/11/11 16:29:49  j_novak
00065  * set_der_0x0 is no longer virtual (to be coherent with Tensor/Scalar classes).
00066  *
00067  * Revision 1.10  2004/06/22 12:48:08  f_limousin
00068  * Change qq, qq_auto and qq_comp to beta, beta_auto and beta_comp.
00069  *
00070  * Revision 1.9  2004/06/07 16:21:35  f_limousin
00071  * Add outputs
00072  *
00073  * Revision 1.8  2004/04/08 16:32:10  f_limousin
00074  * The new variable is ln(Q) instead of Q=psi^2*N. It improves the
00075  * convergence of the code.
00076  *
00077  * Revision 1.7  2004/03/25 10:29:26  j_novak
00078  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00079  *
00080  * Revision 1.6  2004/03/08 11:48:00  f_limousin
00081  * Error in del_deriv() and set_der_0x0() : p_mass_b and p_mass_g were
00082  * missing. And so they were never recomputed.
00083  *
00084  * Revision 1.5  2004/02/27 09:36:46  f_limousin
00085  * u_euler is now constructed on a cartesian basis instead
00086  * of a spherical one.
00087  *
00088  * Revision 1.4  2004/02/21 17:05:13  e_gourgoulhon
00089  * Method Scalar::point renamed Scalar::val_grid_point.
00090  * Method Scalar::set_point renamed Scalar::set_grid_point.
00091  *
00092  * Revision 1.3  2004/01/20 15:16:58  f_limousin
00093  * First version
00094  *
00095  *
00096  * $Header: /cvsroot/Lorene/C++/Source/Star/star.C,v 1.19 2010/02/02 12:45:16 e_gourgoulhon Exp $
00097  *
00098  */
00099 
00100 // Headers C
00101 #include "math.h"
00102 
00103 // Headers Lorene
00104 #include "star.h"
00105 #include "eos.h"
00106 #include "utilitaires.h"
00107 #include "param.h"
00108 #include "unites.h"
00109 
00110 
00111                 //--------------//
00112                 // Constructors //
00113                 //--------------//
00114 
00115 // Standard constructor
00116 // --------------------
00117 Star::Star(Map& mpi, int nzet_i, const Eos& eos_i)
00118          : mp(mpi), 
00119            nzet(nzet_i), 
00120            eos(eos_i), 
00121            ent(mpi), 
00122            nbar(mpi), 
00123            ener(mpi), 
00124            press(mpi),  
00125            ener_euler(mpi), 
00126            s_euler(mpi), 
00127            gam_euler(mpi), 
00128            u_euler(mpi, CON, mp.get_bvect_spher()), 
00129            stress_euler(mpi, CON, mp.get_bvect_spher()),
00130            logn(mpi), 
00131            nn(mpi), 
00132            beta(mpi, CON, mp.get_bvect_spher()),
00133            lnq(mpi),
00134            gamma(mp.flat_met_spher()){
00135     
00136  
00137     // Check of the EOS
00138 //     const Eos_poly_newt* p_eos_poly_newt = 
00139 //              dynamic_cast<const Eos_poly_newt*>( &eos ) ; 
00140 //    
00141 //     const Eos_incomp_newt* p_eos_incomp_newt = 
00142 //              dynamic_cast<const Eos_incomp_newt*>( &eos ) ; 
00143 //    
00144 // 
00145 //     
00146 //     if (p_eos_poly_newt != 0x0) {
00147 //  cout << 
00148 //      "Star::Star : the EOS Eos_poly_newt must not be employed"
00149 //       << " for a relativistic star ! " << endl ; 
00150 //  cout << "(Use Eos_poly instead)" << endl ; 
00151 //  abort() ; 
00152 //  }
00153 //     if (p_eos_incomp_newt != 0x0) {
00154 //  cout << 
00155 //      "Star::Star : the EOS Eos_incomp_newt must not be employed"
00156 //       << " for a relativistic star ! " << endl ; 
00157 //  cout << "(Use Eos_incomp instead)" << endl ; 
00158 //  abort() ; 
00159 //     }
00160 //  
00161     // Pointers of derived quantities initialized to zero : 
00162     set_der_0x0() ;
00163 
00164     // All the matter quantities are initialized to zero :
00165     nbar = 0 ; 
00166     ener = 0 ; 
00167     press = 0 ; 
00168     ent = 0 ; 
00169     ener_euler = 0 ; 
00170     s_euler = 0 ; 
00171     gam_euler = 1. ; 
00172     gam_euler.std_spectral_base() ; 
00173     u_euler.set_etat_zero() ; 
00174     stress_euler.set_etat_zero() ;
00175 
00176     // The metric is initialized to the flat one : 
00177     Metric flat(mp.flat_met_spher()) ;
00178     flat.cov() ;
00179     gamma = flat ;
00180 
00181     logn = 0 ; 
00182     nn = 1. ; 
00183     nn.std_spectral_base() ; 
00184     beta.set_etat_zero() ; 
00185     lnq = 0 ;
00186 
00187 }
00188 
00189 // Copy constructor
00190 // ----------------
00191 Star::Star(const Star& et) 
00192          : mp(et.mp), 
00193            nzet(et.nzet), 
00194            eos(et.eos), 
00195            ent(et.ent), 
00196            nbar(et.nbar), 
00197            ener(et.ener), 
00198            press(et.press),  
00199            ener_euler(et.ener_euler), 
00200            s_euler(et.s_euler), 
00201            gam_euler(et.gam_euler), 
00202            u_euler(et.u_euler), 
00203            stress_euler(et.stress_euler),
00204            logn(et.logn), 
00205            nn(et.nn), 
00206            beta(et.beta),
00207            lnq(et.lnq),
00208            gamma(et.gamma){
00209            
00210     set_der_0x0() ;
00211 
00212 }
00213 
00214 // Constructor from a file
00215 // -----------------------
00216 Star::Star(Map& mpi, const Eos& eos_i, FILE* fich)
00217          : mp(mpi), 
00218            eos(eos_i), 
00219            ent(mpi), 
00220            nbar(mpi), 
00221            ener(mpi), 
00222            press(mpi),  
00223            ener_euler(mpi), 
00224            s_euler(mpi), 
00225            gam_euler(mpi), 
00226            u_euler(mpi, CON, mp.get_bvect_spher()), 
00227            stress_euler(mpi, CON, mp.get_bvect_spher()), 
00228            logn(mpi, *(mpi.get_mg()), fich), 
00229            nn(mpi), 
00230            beta(mpi, CON, mp.get_bvect_spher()),
00231            lnq(mpi, *(mpi.get_mg()), fich),
00232            gamma(mpi.flat_met_spher()){
00233 
00234     // Star parameters
00235     // -----------------
00236 
00237     // nzet is read in the file:     
00238     int xx ; 
00239     fread_be(&xx, sizeof(int), 1, fich) ;   
00240     nzet = xx ;
00241             
00242     // Equation of state
00243     // -----------------
00244     
00245     // Read of the saved EOS
00246     Eos* p_eos_file = Eos::eos_from_file(fich) ; 
00247     
00248     // Comparison with the assigned EOS:
00249     if (eos != *p_eos_file) {
00250     cout << 
00251     "Star::Star(const Map&, const Eos&, FILE*) : the EOS given in "
00252     << endl << 
00253     " argument and that read in the file are different !" << endl ; 
00254     abort() ;  
00255     }
00256     
00257     // p_eos_file is no longer required (it was used only for checking the
00258     //  EOS compatibility)
00259     delete p_eos_file ;
00260     
00261     // Read of the saved fields:
00262     // ------------------------
00263     Scalar ent_file(mp, *(mp.get_mg()), fich) ; 
00264     ent = ent_file ; 
00265     u_euler.set_etat_zero() ; 
00266     stress_euler.set_etat_zero() ;
00267     nn = 1 ;
00268     beta.set_etat_zero() ;
00269 
00270     // Pointers of derived quantities initialized to zero 
00271     // --------------------------------------------------
00272     set_der_0x0() ;
00273     
00274 }
00275 
00276                 //------------//
00277                 // Destructor //
00278                 //------------//
00279 
00280 Star::~Star(){
00281 
00282     del_deriv() ; 
00283 
00284 }
00285 
00286 
00287             //----------------------------------//
00288             // Management of derived quantities //
00289             //----------------------------------//
00290 
00291 void Star::del_deriv() const {
00292 
00293     if (p_mass_b != 0x0) delete p_mass_b ; 
00294     if (p_mass_g != 0x0) delete p_mass_g ; 
00295     if (p_ray_eq != 0x0) delete p_ray_eq ; 
00296     if (p_ray_eq_pis2 != 0x0) delete p_ray_eq_pis2 ; 
00297     if (p_ray_eq_pi != 0x0) delete p_ray_eq_pi ; 
00298     if (p_ray_eq_3pis2 != 0x0) delete p_ray_eq_3pis2 ;
00299     if (p_ray_pole != 0x0) delete p_ray_pole ; 
00300     if (p_l_surf != 0x0) delete p_l_surf ; 
00301     if (p_xi_surf != 0x0) delete p_xi_surf ; 
00302 
00303     Star::set_der_0x0() ; 
00304 }               
00305 
00306 
00307 
00308 
00309 void Star::set_der_0x0() const {
00310 
00311     p_mass_b = 0x0 ; 
00312     p_mass_g = 0x0 ; 
00313     p_ray_eq = 0x0 ; 
00314     p_ray_eq_pis2 = 0x0 ; 
00315     p_ray_eq_pi = 0x0 ; 
00316     p_ray_eq_3pis2 = 0x0 ;
00317     p_ray_pole = 0x0 ; 
00318     p_l_surf = 0x0 ; 
00319     p_xi_surf = 0x0 ; 
00320 
00321 }               
00322 
00323 void Star::del_hydro_euler() {
00324 
00325     ener_euler.set_etat_nondef() ; 
00326     s_euler.set_etat_nondef() ; 
00327     gam_euler.set_etat_nondef() ; 
00328     u_euler.set_etat_nondef() ; 
00329     stress_euler.set_etat_nondef() ;
00330 
00331     del_deriv() ; 
00332 
00333 }               
00334 
00335 
00336 
00337 
00338                 //--------------//
00339                 //  Assignment  //
00340                 //--------------//
00341 
00342 // Assignment to another Star
00343 // ----------------------------
00344 void Star::operator=(const Star& et) {
00345 
00346     assert( &(et.mp) == &mp ) ;         // Same mapping
00347     assert( &(et.eos) == &eos ) ;       // Same EOS
00348     
00349     nzet = et.nzet ; 
00350     ent = et.ent ;
00351     nbar = et.nbar ; 
00352     ener = et.ener ;
00353     press = et.press ;
00354     ener_euler = et.ener_euler ;
00355     s_euler = et.s_euler ;
00356     gam_euler = et.gam_euler ;
00357     u_euler = et.u_euler ;
00358     stress_euler = et.stress_euler ;
00359     logn = et.logn ;
00360     nn = et.nn ;
00361     beta = et.beta ;
00362     lnq = et.lnq ;
00363     gamma = et.gamma ;
00364 
00365     del_deriv() ;  // Deletes all derived quantities
00366 
00367 }   
00368 
00369 // Assignment of the enthalpy field
00370 // --------------------------------
00371 
00372 void Star::set_enthalpy(const Scalar& ent_i) {
00373     
00374     ent = ent_i ; 
00375     
00376     // Update of (nbar, ener, press) :
00377     equation_of_state() ; 
00378     
00379     // The derived quantities are obsolete:
00380     del_deriv() ; 
00381     
00382 }
00383 
00384                 //--------------//
00385                 //    Outputs   //
00386                 //--------------//
00387 
00388 // Save in a file
00389 // --------------
00390 void Star::sauve(FILE* fich) const {
00391 
00392     logn.sauve(fich) ;
00393     lnq.sauve(fich) ;
00394   
00395     int xx = nzet ;     
00396     fwrite_be(&xx, sizeof(int), 1, fich) ;          
00397   
00398     eos.sauve(fich) ; 
00399     ent.sauve(fich) ;     
00400 }
00401 
00402 // Printing
00403 // --------
00404 
00405 ostream& operator<<(ostream& ost, const Star& et)  {
00406     et >> ost ;
00407     return ost ;
00408 }
00409     
00410 ostream& Star::operator>>(ostream& ost) const {
00411     
00412   using namespace Unites ;
00413 
00414     ost << endl ; 
00415       
00416     ost << "Number of domains occupied by the star : " << nzet << endl ; 
00417     
00418     ost << "Equation of state : " << endl ; 
00419     ost << eos << endl ; 
00420     
00421     ost << endl << "Central enthalpy : " << ent.val_grid_point(0,0,0,0) << " c^2" << endl ; 
00422     ost << "Central proper baryon density : " << nbar.val_grid_point(0,0,0,0) 
00423     << " x 0.1 fm^-3" << endl ; 
00424     ost << "Central proper energy density : " << ener.val_grid_point(0,0,0,0) 
00425     << " rho_nuc c^2" << endl ; 
00426     ost << "Central pressure : " << press.val_grid_point(0,0,0,0) 
00427     << " rho_nuc c^2" << endl ; 
00428     
00429     ost << endl ;
00430     ost << "Central lapse N :      " << nn.val_grid_point(0,0,0,0) <<  endl ; 
00431 //    ost << "Central value of lnq : " << lnq.val_grid_point(0,0,0,0) <<  endl ; 
00432   
00433     ost << endl 
00434     << "Coordinate equatorial radius (phi=0) a1 =    " 
00435     << ray_eq()/km << " km" << endl ;  
00436     ost << "Coordinate equatorial radius (phi=pi/2) a2 = " 
00437     << ray_eq_pis2()/km << " km" << endl ;  
00438     ost << "Coordinate equatorial radius (phi=pi):       " 
00439     << ray_eq_pi()/km << " km" << endl ;  
00440     ost << "Coordinate polar radius a3 =                 " 
00441     << ray_pole()/km << " km" << endl ;  
00442     ost << "Axis ratio a2/a1 = " << ray_eq_pis2() / ray_eq() 
00443     << "  a3/a1 = " << ray_pole() / ray_eq() << endl ;  
00444     ost << endl << "Baryon mass :        " << mass_b() / msol << " M_sol" << endl ; 
00445     ost << "Gravitational mass : " << mass_g() / msol << " M_sol" << endl ; 
00446     
00447 
00448     return ost ; 
00449 }
00450 
00451         //-----------------------------------------//
00452         //  Computation of hydro quantities    //
00453         //-----------------------------------------//
00454 
00455 void Star::equation_of_state() {
00456 
00457     Scalar ent_eos = ent ;
00458 
00459 
00460     // Slight rescale of the enthalpy field in case of 2 domains inside the
00461     //  star
00462 
00463 
00464         double epsilon = 1.e-12 ;
00465 
00466     const Mg3d* mg = mp.get_mg() ;
00467         int nz = mg->get_nzone() ;
00468 
00469         Mtbl xi(mg) ;
00470         xi.set_etat_qcq() ;
00471         for (int l=0; l<nz; l++) {
00472             xi.t[l]->set_etat_qcq() ;
00473             for (int k=0; k<mg->get_np(l); k++) {
00474                 for (int j=0; j<mg->get_nt(l); j++) {
00475                     for (int i=0; i<mg->get_nr(l); i++) {
00476                         xi.set(l,k,j,i) =
00477                             mg->get_grille3d(l)->x[i] ;
00478                     }
00479                 }
00480             }
00481 
00482         }
00483 
00484         Scalar fact_ent(mp) ;
00485         fact_ent.allocate_all() ;
00486         
00487         fact_ent.set_domain(0) = 1 + epsilon * xi(0) * xi(0) ;
00488         fact_ent.set_domain(1) = 1 - 0.25 * epsilon * (xi(1) - 1) * (xi(1) - 1) ;
00489         
00490         for (int l=nzet; l<nz; l++) {
00491             fact_ent.set_domain(l) = 1 ;
00492         }
00493 
00494     if (nzet > 1) {
00495 
00496         if (nzet > 2) {
00497         
00498             cout << "Star::equation_of_state: not ready yet for nzet > 2 !"
00499                  << endl ;      
00500         }
00501 
00502         ent_eos = fact_ent * ent_eos ;
00503         ent_eos.std_spectral_base() ;
00504     }
00505 
00506 
00507 
00508 
00509 
00510     // Call to the EOS (the EOS is called domain by domain in order to
00511     //          allow for the use of MEos)
00512 
00513     Scalar tempo(mp) ;
00514 
00515     nbar.set_etat_qcq() ;
00516     nbar = 0 ;
00517     for (int l=0; l<nzet; l++) {
00518 
00519         Param par ;       // Paramater for multi-domain equation of state
00520         par.add_int(l) ;
00521 
00522         tempo =  eos.nbar_ent(ent_eos, 1, l, &par) ;
00523 
00524         nbar = nbar + tempo ;
00525 
00526     }
00527 
00528     ener.set_etat_qcq() ;
00529     ener = 0 ;
00530     for (int l=0; l<nzet; l++) {
00531 
00532         Param par ;    // Paramater for multi-domain equation of state
00533         par.add_int(l) ;
00534 
00535         tempo =  eos.ener_ent(ent_eos, 1, l, &par) ;
00536 
00537         ener = ener + tempo ;
00538 
00539     }
00540 
00541     press.set_etat_qcq() ;
00542     press = 0 ;
00543     for (int l=0; l<nzet; l++) {
00544 
00545         Param par ;     // Paramater for multi-domain equation of state
00546         par.add_int(l) ;
00547 
00548         tempo =  eos.press_ent(ent_eos, 1, l, &par) ;
00549 
00550         press = press + tempo ;
00551 
00552     }
00553 
00554 
00555     // Set the bases for spectral expansion
00556     nbar.std_spectral_base() ; 
00557     ener.std_spectral_base() ; 
00558     press.std_spectral_base() ; 
00559 
00560     // The derived quantities are obsolete
00561     del_deriv() ; 
00562     
00563 }
00564 
00565 void Star::hydro_euler() {
00566     
00567     cout << 
00568     "Star::hydro_euler : hydro_euler must be called via a derived class"
00569     << endl << " of Star !" << endl ; 
00570     
00571     abort() ;        
00572     
00573 }

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