etoile.C

00001 /*
00002  * Methods of class Etoile
00003  *
00004  * (see file etoile.h for documentation)
00005  */
00006 
00007 /*
00008  *   Copyright (c) 2000-2001 Eric Gourgoulhon
00009  *   Copyright (c) 2000-2001 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 as published by
00015  *   the Free Software Foundation; either version 2 of the License, or
00016  *   (at your option) any later version.
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 
00030 char etoile_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/etoile.C,v 1.8 2005/01/18 22:36:50 k_taniguchi Exp $" ;
00031 
00032 /*
00033  * $Id: etoile.C,v 1.8 2005/01/18 22:36:50 k_taniguchi Exp $
00034  * $Log: etoile.C,v $
00035  * Revision 1.8  2005/01/18 22:36:50  k_taniguchi
00036  * Delete a pointer for ray_eq(int kk).
00037  *
00038  * Revision 1.7  2005/01/18 20:35:05  k_taniguchi
00039  * Addition of ray_eq(int kk).
00040  *
00041  * Revision 1.6  2005/01/17 20:40:25  k_taniguchi
00042  * Addition of ray_eq_3pis2().
00043  *
00044  * Revision 1.5  2004/03/25 10:29:06  j_novak
00045  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00046  *
00047  * Revision 1.4  2003/10/13 15:23:56  f_limousin
00048  * *** empty log message ***
00049  *
00050  * Revision 1.3  2002/04/09 14:32:15  e_gourgoulhon
00051  * 1/ Added extra parameters in EOS computational functions (argument par)
00052  * 2/ New class MEos for multi-domain EOS
00053  *
00054  * Revision 1.2  2001/12/04 21:27:53  e_gourgoulhon
00055  *
00056  * All writing/reading to a binary file are now performed according to
00057  * the big endian convention, whatever the system is big endian or
00058  * small endian, thanks to the functions fwrite_be and fread_be
00059  *
00060  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00061  * LORENE
00062  *
00063  * Revision 2.14  2000/11/24  13:27:44  eric
00064  * Dans eqution_of_state(): changement leger de ent dans le cas ou l'on a
00065  * deux domaine avant d'appeler l'EOS.
00066  *
00067  * Revision 2.13  2000/09/25  12:22:02  keisuke
00068  * *** empty log message ***
00069  *
00070  * Revision 2.12  2000/09/22  15:50:58  keisuke
00071  * Ajout du membre d_logn_auto_div.
00072  *
00073  * Revision 2.11  2000/09/07  14:34:09  keisuke
00074  * Ajout du membre logn_auto_regu.
00075  *
00076  * Revision 2.10  2000/08/31  15:36:54  eric
00077  * Bases spectrales standards pour nnn, a_car et gam_euler dans le
00078  * constructeur (initialisation a la metrique plate).
00079  *
00080  * Revision 2.9  2000/08/29  11:37:49  eric
00081  * Ajout des membres k_div et logn_auto_div.
00082  *
00083  * Revision 2.8  2000/07/21  12:01:11  eric
00084  * Modif dans Etoile::del_deriv() :
00085  *   appel de Etoile::set_der_0x0() et non de la fonction virtuelle set_der_0x0().
00086  *
00087  * Revision 2.7  2000/03/21  12:39:34  eric
00088  * Le constructeur standard teste la compatibilite de l'EOS avec le
00089  * caractere relativiste de l'etoile.
00090  *
00091  * Revision 2.6  2000/02/21  14:32:40  eric
00092  * gam_euler est initialise a 1 dans le constructeur standard.
00093  * Suppression de l'appel a del_hydro_euler dans equation_of_state().
00094  *
00095  * Revision 2.5  2000/02/09  19:30:47  eric
00096  * La triade de decomposition doit desormais figurer en argument des
00097  *  constructeurs de Tenseur.
00098  *
00099  * Revision 2.4  2000/02/02  09:23:34  eric
00100  * Affichage de la masse.
00101  *
00102  * Revision 2.3  2000/01/28  17:18:10  eric
00103  * Modifs noms des quantites globales.
00104  * Affichage.
00105  *
00106  * Revision 2.2  2000/01/24  17:13:36  eric
00107  * Le mapping mp n'est plus constant.
00108  *
00109  * Revision 2.1  2000/01/24  13:37:22  eric
00110  * *** empty log message ***
00111  *
00112  * Revision 2.0  2000/01/20  17:04:45  eric
00113  * *** empty log message ***
00114  *
00115  *
00116  * $Header: /cvsroot/Lorene/C++/Source/Etoile/etoile.C,v 1.8 2005/01/18 22:36:50 k_taniguchi Exp $
00117  *
00118  */
00119 
00120 // Headers C
00121 #include "math.h"
00122 
00123 // Headers Lorene
00124 #include "etoile.h"
00125 #include "eos.h"
00126 #include "utilitaires.h"
00127 #include "param.h"
00128 #include "unites.h"
00129 
00130                 //--------------//
00131                 // Constructors //
00132                 //--------------//
00133 
00134 // Standard constructor
00135 // --------------------
00136 Etoile::Etoile(Map& mpi, int nzet_i, bool relat, const Eos& eos_i)
00137          : mp(mpi), 
00138            nzet(nzet_i), 
00139            relativistic(relat), 
00140            k_div(0), 
00141            eos(eos_i), 
00142            ent(mpi), 
00143            nbar(mpi), 
00144            ener(mpi), 
00145            press(mpi),  
00146            ener_euler(mpi), 
00147            s_euler(mpi), 
00148            gam_euler(mpi), 
00149            u_euler(mpi, 1, CON, mp.get_bvect_cart()), 
00150            logn_auto(mpi), 
00151            logn_auto_regu(mpi), 
00152            logn_auto_div(mpi),
00153            d_logn_auto_div(mpi, 1, COV, mp.get_bvect_spher()), 
00154            beta_auto(mpi), 
00155            nnn(mpi), 
00156            shift(mpi, 1, CON, mp.get_bvect_cart()),
00157            a_car(mpi) {
00158 
00159     // Check of the EOS
00160     const Eos_poly* p_eos_poly = dynamic_cast<const Eos_poly*>( &eos ) ;      
00161     const Eos_poly_newt* p_eos_poly_newt = 
00162                 dynamic_cast<const Eos_poly_newt*>( &eos ) ;      
00163     const Eos_incomp* p_eos_incomp = dynamic_cast<const Eos_incomp*>( &eos ) ;    
00164     const Eos_incomp_newt* p_eos_incomp_newt = 
00165                 dynamic_cast<const Eos_incomp_newt*>( &eos ) ;    
00166 
00167     if (relativistic) {
00168 
00169     if (p_eos_poly_newt != 0x0) {
00170         cout << 
00171         "Etoile::Etoile : the EOS Eos_poly_newt must not be employed"
00172         << " for a relativistic star ! " << endl ; 
00173         cout << "(Use Eos_poly instead)" << endl ; 
00174         abort() ; 
00175     }
00176     if (p_eos_incomp_newt != 0x0) {
00177         cout << 
00178         "Etoile::Etoile : the EOS Eos_incomp_newt must not be employed"
00179         << " for a relativistic star ! " << endl ; 
00180         cout << "(Use Eos_incomp instead)" << endl ; 
00181         abort() ; 
00182     }
00183 
00184     }
00185     else{
00186 
00187     if ( (p_eos_poly != 0x0) && (p_eos_poly_newt == 0x0) ) {
00188         cout << 
00189         "Etoile::Etoile : the EOS Eos_poly must not be employed"
00190         << " for a Newtonian star ! " << endl ; 
00191         cout << "(Use Eos_poly_newt instead)" << endl ; 
00192         abort() ; 
00193     }
00194     if ( (p_eos_incomp != 0x0) && (p_eos_incomp_newt == 0x0) ) {
00195         cout << 
00196         "Etoile::Etoile : the EOS Eos_incomp must not be employed"
00197         << " for a relativistic star ! " << endl ; 
00198         cout << "(Use Eos_incomp_newt instead)" << endl ; 
00199         abort() ; 
00200     }
00201     
00202     }
00203 
00204 
00205     // Parameter 1/c^2      
00206     unsurc2 = relativistic ? double(1) : double(0) ; 
00207 
00208     // Pointers of derived quantities initialized to zero : 
00209     set_der_0x0() ;
00210 
00211     // All the matter quantities are initialized to zero :
00212     nbar = 0 ; 
00213     ener = 0 ; 
00214     press = 0 ; 
00215     ent = 0 ; 
00216     ener_euler = 0 ; 
00217     s_euler = 0 ; 
00218     gam_euler = 1 ; 
00219     gam_euler.set_std_base() ; 
00220     u_euler = 0 ; 
00221 
00222     // The metric is initialized to the flat one : 
00223     logn_auto = 0 ; 
00224     logn_auto_regu = 0 ; 
00225     logn_auto_div = 0 ; 
00226     d_logn_auto_div = 0 ; 
00227     beta_auto = 0 ; 
00228     nnn = 1 ; 
00229     nnn.set_std_base() ; 
00230     shift = 0 ; 
00231     a_car = 1 ; 
00232     a_car.set_std_base() ; 
00233     
00234 }
00235 
00236 // Copy constructor
00237 // ----------------
00238 Etoile::Etoile(const Etoile& et) 
00239          : mp(et.mp), 
00240            nzet(et.nzet), 
00241            relativistic(et.relativistic), 
00242            unsurc2(et.unsurc2), 
00243            k_div(et.k_div), 
00244            eos(et.eos), 
00245            ent(et.ent), 
00246            nbar(et.nbar), 
00247            ener(et.ener), 
00248            press(et.press),  
00249            ener_euler(et.ener_euler), 
00250            s_euler(et.s_euler), 
00251            gam_euler(et.gam_euler), 
00252            u_euler(et.u_euler), 
00253            logn_auto(et.logn_auto), 
00254            logn_auto_regu(et.logn_auto_regu),
00255            logn_auto_div(et.logn_auto_div), 
00256            d_logn_auto_div(et.d_logn_auto_div), 
00257            beta_auto(et.beta_auto), 
00258            nnn(et.nnn), 
00259            shift(et.shift), 
00260            a_car(et.a_car) {
00261            
00262     set_der_0x0() ;
00263 
00264 }
00265 
00266 // Constructor from a file
00267 // -----------------------
00268 Etoile::Etoile(Map& mpi, const Eos& eos_i, FILE* fich)
00269          : mp(mpi), 
00270            eos(eos_i), 
00271            ent(mpi), 
00272            nbar(mpi), 
00273            ener(mpi), 
00274            press(mpi),  
00275            ener_euler(mpi), 
00276            s_euler(mpi), 
00277            gam_euler(mpi), 
00278            u_euler(mpi, 1, CON, mp.get_bvect_cart()), 
00279            logn_auto(mpi), 
00280            logn_auto_regu(mpi), 
00281            logn_auto_div(mpi),
00282            d_logn_auto_div(mpi, 1, COV, mp.get_bvect_spher()), 
00283            beta_auto(mpi), 
00284            nnn(mpi), 
00285            shift(mpi, 1, CON, mp.get_bvect_cart()),
00286            a_car(mpi) {
00287 
00288     // Etoile parameters
00289     // -----------------
00290 
00291     // nzet and relativistic are read in the file:     
00292     int xx ; 
00293     fread_be(&xx, sizeof(int), 1, fich) ;   
00294     k_div = xx / 1000 ;     // integer part
00295     nzet = xx - k_div * 1000  ;
00296             
00297     fread(&relativistic, sizeof(bool), 1, fich) ;       
00298           
00299     // Parameter 1/c^2 is deduced from relativistic:
00300     unsurc2 = relativistic ? double(1) : double(0) ; 
00301 
00302 
00303     // Equation of state
00304     // -----------------
00305     
00306     // Read of the saved EOS
00307     Eos* p_eos_file = Eos::eos_from_file(fich) ; 
00308     
00309     // Comparison with the assigned EOS:
00310     if (eos != *p_eos_file) {
00311     cout << 
00312     "Etoile::Etoile(const Map&, const Eos&, FILE*) : the EOS given in "
00313     << endl << 
00314     " argument and that read in the file are different !" << endl ; 
00315     abort() ;  
00316     }
00317     
00318     // p_eos_file is no longer required (it was used only for checking the
00319     //  EOS compatibility)
00320     delete p_eos_file ;
00321     
00322     // Read of the saved fields:
00323     // ------------------------
00324     Tenseur ent_file(mp, fich) ; 
00325     ent = ent_file ; 
00326     
00327     Tenseur logn_auto_file(mp, fich) ; 
00328     logn_auto = logn_auto_file ; 
00329     
00330     Tenseur beta_auto_file(mp, fich) ; 
00331     beta_auto = beta_auto_file ; 
00332 
00333     if (k_div == 0) {
00334     logn_auto_div = 0 ; 
00335     d_logn_auto_div = 0 ; 
00336     }
00337     else {
00338 
00339     Tenseur logn_auto_div_file(mp, fich) ; 
00340     logn_auto_div = logn_auto_div_file ; 
00341 
00342     Tenseur d_logn_auto_div_file(mp, mp.get_bvect_spher(), fich) ; 
00343     d_logn_auto_div = d_logn_auto_div_file ; 
00344     }
00345 
00346     logn_auto_regu = logn_auto - logn_auto_div ;
00347 
00348      shift = 0 ;    
00349 
00350     // Pointers of derived quantities initialized to zero 
00351     // --------------------------------------------------
00352     set_der_0x0() ;
00353     
00354 }
00355 
00356                 //------------//
00357                 // Destructor //
00358                 //------------//
00359 
00360 Etoile::~Etoile(){
00361 
00362     del_deriv() ; 
00363 
00364 }
00365 
00366 
00367             //----------------------------------//
00368             // Management of derived quantities //
00369             //----------------------------------//
00370 
00371 void Etoile::del_deriv() const {
00372 
00373     if (p_mass_b != 0x0) delete p_mass_b ; 
00374     if (p_mass_g != 0x0) delete p_mass_g ; 
00375     if (p_ray_eq != 0x0) delete p_ray_eq ; 
00376     if (p_ray_eq_pis2 != 0x0) delete p_ray_eq_pis2 ; 
00377     if (p_ray_eq_pi != 0x0) delete p_ray_eq_pi ; 
00378     if (p_ray_eq_3pis2 != 0x0) delete p_ray_eq_3pis2 ; 
00379     if (p_ray_pole != 0x0) delete p_ray_pole ; 
00380     if (p_l_surf != 0x0) delete p_l_surf ; 
00381     if (p_xi_surf != 0x0) delete p_xi_surf ; 
00382 
00383     Etoile::set_der_0x0() ; 
00384 }               
00385 
00386 
00387 
00388 
00389 void Etoile::set_der_0x0() const {
00390 
00391     p_mass_b = 0x0 ; 
00392     p_mass_g = 0x0 ; 
00393     p_ray_eq = 0x0 ; 
00394     p_ray_eq_pis2 = 0x0 ; 
00395     p_ray_eq_pi = 0x0 ; 
00396     p_ray_eq_3pis2 = 0x0 ; 
00397     p_ray_pole = 0x0 ; 
00398     p_l_surf = 0x0 ; 
00399     p_xi_surf = 0x0 ; 
00400 
00401 }               
00402 
00403 void Etoile::del_hydro_euler() {
00404 
00405     ener_euler.set_etat_nondef() ; 
00406     s_euler.set_etat_nondef() ; 
00407     gam_euler.set_etat_nondef() ; 
00408     u_euler.set_etat_nondef() ; 
00409 
00410     del_deriv() ; 
00411 
00412 }               
00413 
00414 
00415 
00416 
00417                 //--------------//
00418                 //  Assignment  //
00419                 //--------------//
00420 
00421 // Assignment to another Etoile
00422 // ----------------------------
00423 void Etoile::operator=(const Etoile& et) {
00424 
00425     assert( &(et.mp) == &mp ) ;         // Same mapping
00426     assert( &(et.eos) == &eos ) ;       // Same EOS
00427     
00428     nzet = et.nzet ; 
00429     relativistic = et.relativistic ; 
00430     k_div = et.k_div ; 
00431     unsurc2 = et.unsurc2 ; 
00432     
00433     ent = et.ent ;
00434     nbar = et.nbar ; 
00435     ener = et.ener ;
00436     press = et.press ;
00437     ener_euler = et.ener_euler ;
00438     s_euler = et.s_euler ;
00439     gam_euler = et.gam_euler ;
00440     u_euler = et.u_euler ;
00441     logn_auto = et.logn_auto ;
00442     logn_auto_regu = et.logn_auto_regu ;
00443     logn_auto_div = et.logn_auto_div ;
00444     d_logn_auto_div = et.d_logn_auto_div ;
00445     beta_auto = et.beta_auto ;
00446     nnn = et.nnn ;
00447     shift = et.shift ;
00448     a_car = et.a_car ;
00449     
00450     
00451     del_deriv() ;  // Deletes all derived quantities
00452 
00453 }   
00454 
00455 // Assignment of the enthalpy field
00456 // --------------------------------
00457 
00458 void Etoile::set_enthalpy(const Cmp& ent_i) {
00459     
00460     ent = ent_i ; 
00461     
00462     // Update of (nbar, ener, press) :
00463     equation_of_state() ; 
00464     
00465     // The derived quantities are obsolete:
00466     del_deriv() ; 
00467     
00468 }
00469 
00470                 //--------------//
00471                 //    Outputs   //
00472                 //--------------//
00473 
00474 // Save in a file
00475 // --------------
00476 void Etoile::sauve(FILE* fich) const {
00477     
00478     int xx = nzet + k_div * 1000 ;     
00479     fwrite_be(&xx, sizeof(int), 1, fich) ;          
00480 
00481     fwrite(&relativistic, sizeof(bool), 1, fich) ;      
00482 
00483     eos.sauve(fich) ; 
00484     
00485     ent.sauve(fich) ;     
00486     logn_auto.sauve(fich) ;     
00487     beta_auto.sauve(fich) ;  
00488     
00489     if (k_div != 0) {
00490     logn_auto_div.sauve(fich) ; 
00491     d_logn_auto_div.sauve(fich) ; 
00492     }   
00493     
00494 }
00495 
00496 // Printing
00497 // --------
00498 
00499 ostream& operator<<(ostream& ost, const Etoile& et)  {
00500     et >> ost ;
00501     return ost ;
00502 }
00503     
00504 ostream& Etoile::operator>>(ostream& ost) const {
00505     
00506   using namespace Unites ;
00507 
00508     ost << endl ; 
00509     if (relativistic) {
00510     ost << "Relativistic star" << endl ; 
00511     ost << "-----------------" << endl ; 
00512     }
00513     else {
00514     ost << "Newtonian star" << endl ; 
00515     ost << "--------------" << endl ; 
00516     }
00517     
00518     ost << "Number of domains occupied by the star : " << nzet << endl ; 
00519     
00520     ost << "Equation of state : " << endl ; 
00521     ost << eos << endl ; 
00522     
00523     ost << endl << "Central enthalpy : " << ent()(0,0,0,0) << " c^2" << endl ; 
00524     ost << "Central proper baryon density : " << nbar()(0,0,0,0) 
00525     << " x 0.1 fm^-3" << endl ; 
00526     ost << "Central proper energy density : " << ener()(0,0,0,0) 
00527     << " rho_nuc c^2" << endl ; 
00528     ost << "Central pressure : " << press()(0,0,0,0) 
00529     << " rho_nuc c^2" << endl ; 
00530 
00531     ost << endl 
00532     << "Regularization index of the gravitational potential : k_div = "
00533     << k_div << endl ; 
00534     ost << "Central lapse N :      " << nnn()(0,0,0,0) <<  endl ; 
00535     ost <<     "Central value of A^2 : " << a_car()(0,0,0,0) <<  endl ; 
00536 
00537     ost << endl 
00538     << "Coordinate equatorial radius (phi=0) a1 =    " 
00539     << ray_eq()/km << " km" << endl ;  
00540     ost << "Coordinate equatorial radius (phi=pi/2) a2 = " 
00541     << ray_eq_pis2()/km << " km" << endl ;  
00542     ost << "Coordinate equatorial radius (phi=pi):       " 
00543     << ray_eq_pi()/km << " km" << endl ;  
00544     ost << "Coordinate polar radius a3 =                 " 
00545     << ray_pole()/km << " km" << endl ;  
00546     ost << "Axis ratio a2/a1 = " << ray_eq_pis2() / ray_eq() 
00547     << "  a3/a1 = " << ray_pole() / ray_eq() << endl ;  
00548 
00549     ost << endl << "Baryon mass :        " << mass_b() / msol << " M_sol" << endl ; 
00550     ost << "Gravitational mass : " << mass_g() / msol << " M_sol" << endl ; 
00551     
00552     return ost ; 
00553 }
00554 
00555         //-----------------------------------------//
00556         //  Computation of hydro quantities    //
00557         //-----------------------------------------//
00558 
00559 void Etoile::equation_of_state() {
00560 
00561     Cmp ent_eos = ent() ;
00562 
00563 
00564     // Slight rescale of the enthalpy field in case of 2 domains inside the
00565     //  star
00566 
00567 
00568         double epsilon = 1.e-12 ;
00569 
00570     const Mg3d* mg = mp.get_mg() ;
00571         int nz = mg->get_nzone() ;
00572 
00573         Mtbl xi(mg) ;
00574         xi.set_etat_qcq() ;
00575         for (int l=0; l<nz; l++) {
00576             xi.t[l]->set_etat_qcq() ;
00577             for (int k=0; k<mg->get_np(l); k++) {
00578                 for (int j=0; j<mg->get_nt(l); j++) {
00579                     for (int i=0; i<mg->get_nr(l); i++) {
00580                         xi.set(l,k,j,i) =
00581                             mg->get_grille3d(l)->x[i] ;
00582                     }
00583                 }
00584             }
00585 
00586         }
00587 
00588         Cmp fact_ent(mp) ;
00589         fact_ent.allocate_all() ;
00590         
00591         fact_ent.set(0) = 1 + epsilon * xi(0) * xi(0) ;
00592         fact_ent.set(1) = 1 - 0.25 * epsilon * (xi(1) - 1) * (xi(1) - 1) ;
00593         
00594         for (int l=nzet; l<nz; l++) {
00595             fact_ent.set(l) = 1 ;
00596         }
00597 
00598     if (nzet > 1) {
00599 
00600         if (nzet > 2) {
00601         
00602             cout << "Etoile::equation_of_state: not ready yet for nzet > 2 !"
00603                  << endl ;      
00604         }
00605 
00606         ent_eos = fact_ent * ent_eos ;
00607         ent_eos.std_base_scal() ;
00608     }
00609 
00610 
00611 
00612 
00613 
00614     // Call to the EOS (the EOS is called domain by domain in order to
00615     //          allow for the use of MEos)
00616 
00617     Cmp tempo(mp) ;
00618 
00619     nbar.set_etat_qcq() ;
00620     nbar.set() = 0 ;
00621     for (int l=0; l<nzet; l++) {
00622 
00623         Param par ;       // Paramater for multi-domain equation of state
00624         par.add_int(l) ;
00625 
00626         tempo =  eos.nbar_ent(ent_eos, 1, l, &par) ;
00627 
00628         nbar = nbar() + tempo ;
00629 
00630     }
00631 
00632     ener.set_etat_qcq() ;
00633     ener.set() = 0 ;
00634     for (int l=0; l<nzet; l++) {
00635 
00636         Param par ;    // Paramater for multi-domain equation of state
00637         par.add_int(l) ;
00638 
00639         tempo =  eos.ener_ent(ent_eos, 1, l, &par) ;
00640 
00641         ener = ener() + tempo ;
00642 
00643     }
00644 
00645     press.set_etat_qcq() ;
00646     press.set() = 0 ;
00647     for (int l=0; l<nzet; l++) {
00648 
00649         Param par ;     // Paramater for multi-domain equation of state
00650         par.add_int(l) ;
00651 
00652         tempo =  eos.press_ent(ent_eos, 1, l, &par) ;
00653 
00654         press = press() + tempo ;
00655 
00656     }
00657 
00658 
00659     // Set the bases for spectral expansion
00660     nbar.set_std_base() ; 
00661     ener.set_std_base() ; 
00662     press.set_std_base() ; 
00663 
00664     // The Eulerian quantities are obsolete
00665     //## del_hydro_euler() ; 
00666     
00667     // The derived quantities are obsolete
00668     del_deriv() ; 
00669     
00670 }
00671 
00672 void Etoile::hydro_euler() {
00673     
00674     cout << 
00675     "Etoile::hydro_euler : hydro_euler must be called via a derived class"
00676     << endl << " of Etoile !" << endl ; 
00677     
00678     abort() ;        
00679     
00680 }

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