eos_tabul.C

00001 /*
00002  *  Methods of class Eos_tabul
00003  *
00004  *  (see file eos.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2000-2001 Eric Gourgoulhon
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 eos_tabul_C[] = "$Header: /cvsroot/Lorene/C++/Source/Eos/eos_tabul.C,v 1.10 2010/02/16 11:14:50 j_novak Exp $" ;
00031 
00032 /*
00033  * $Id: eos_tabul.C,v 1.10 2010/02/16 11:14:50 j_novak Exp $
00034  * $Log: eos_tabul.C,v $
00035  * Revision 1.10  2010/02/16 11:14:50  j_novak
00036  * More verbose opeining of the file.
00037  *
00038  * Revision 1.9  2010/02/02 13:22:16  j_novak
00039  * New class Eos_Compstar.
00040  *
00041  * Revision 1.8  2004/03/25 10:29:02  j_novak
00042  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
00043  *
00044  * Revision 1.7  2003/11/25 13:42:50  m_bejger
00045  * read_table written in more ordered way
00046  *
00047  * Revision 1.6  2003/11/21 16:11:16  m_bejger
00048  * added the computation dlnp/dlnn_b, dlnn/dlnH
00049  *
00050  * Revision 1.5  2003/05/30 07:50:06  e_gourgoulhon
00051  *
00052  * Reformulate the computation of der_nbar_ent
00053  * Added computation of der_press_ent_p.
00054  *
00055  * Revision 1.4  2003/05/15 09:42:47  e_gourgoulhon
00056  *
00057  * Now computes d ln / dln H.
00058  *
00059  * Revision 1.3  2002/10/16 14:36:35  j_novak
00060  * Reorganization of #include instructions of standard C++, in order to
00061  * use experimental version 3 of gcc.
00062  *
00063  * Revision 1.2  2002/04/09 14:32:15  e_gourgoulhon
00064  * 1/ Added extra parameters in EOS computational functions (argument par)
00065  * 2/ New class MEos for multi-domain EOS
00066  *
00067  * Revision 1.1.1.1  2001/11/20 15:19:27  e_gourgoulhon
00068  * LORENE
00069  *
00070  * Revision 2.3  2001/09/13  13:35:49  eric
00071  * Suppression des affichages dans read_table().
00072  *
00073  * Revision 2.2  2001/02/07  09:48:05  eric
00074  * Suppression de la fonction derent_ent_p.
00075  * Ajout des fonctions donnant les derivees de l'EOS:
00076  *      der_nbar_ent_p
00077  *      der_ener_ent_p
00078  *      der_press_ent_p
00079  *
00080  * Revision 2.1  2000/11/23  00:10:20  eric
00081  * Enthalpie minimale fixee a 1e-14.
00082  *
00083  * Revision 2.0  2000/11/22  19:31:30  eric
00084  * *** empty log message ***
00085  *
00086  *
00087  * $Header: /cvsroot/Lorene/C++/Source/Eos/eos_tabul.C,v 1.10 2010/02/16 11:14:50 j_novak Exp $
00088  *
00089  */
00090 
00091 // headers C
00092 #include <stdlib.h>
00093 #include <string.h>
00094 #include <math.h>
00095 
00096 // Headers Lorene
00097 #include "eos.h"
00098 #include "tbl.h"
00099 #include "utilitaires.h"
00100 #include "unites.h"
00101 
00102 
00103 void interpol_herm(const Tbl& , const Tbl&, const Tbl&, double, int&,
00104            double&, double& ) ;
00105 
00106 void interpol_linear(const Tbl&, const Tbl&, double, int&, double&) ;
00107 
00108             //----------------------------//
00109             //      Constructors          //
00110             //----------------------------//
00111 
00112 // Standard constructor
00113 // --------------------         
00114 Eos_tabul::Eos_tabul(const char* name_i, const char* table,
00115              const char* path) : Eos(name_i) {  
00116 
00117     strcpy(tablename, path) ;
00118     strcat(tablename, "/") ;
00119     strcat(tablename, table) ;
00120     
00121     read_table() ;  
00122 
00123 }
00124 
00125 // Standard constructor with full filename
00126 // ---------------------------------------
00127 Eos_tabul::Eos_tabul(const char* name_i, const char* file_name) 
00128   : Eos(name_i) {   
00129 
00130     strcpy(tablename, file_name) ;
00131     
00132     read_table() ;  
00133 
00134 }
00135 
00136 
00137 // Constructor from binary file
00138 // ----------------------------
00139 Eos_tabul::Eos_tabul(FILE* fich) : Eos(fich) {
00140 
00141        fread(tablename, sizeof(char), 160, fich) ;      
00142 
00143        read_table() ;
00144 
00145 }
00146 
00147 
00148 
00149 // Constructor from a formatted file
00150 // ---------------------------------
00151 Eos_tabul::Eos_tabul(ifstream& fich, const char* table) : Eos(fich) {
00152 
00153     char path[160] ;    
00154     
00155     fich.getline(path, 160) ;
00156     
00157     strcpy(tablename, path) ;
00158     strcat(tablename, "/") ;
00159     strcat(tablename, table) ;
00160     
00161     read_table() ;  
00162 
00163 }
00164 
00165 Eos_tabul::Eos_tabul(ifstream& fich) : Eos(fich) {
00166 
00167     char file_name[160] ;   
00168     
00169     fich.getline(file_name, 160) ;
00170     
00171     strcpy(tablename, file_name) ;
00172     
00173     read_table() ;  
00174 
00175 }
00176 
00177 
00178             //--------------//
00179             //  Destructor  //
00180             //--------------//
00181 
00182 Eos_tabul::~Eos_tabul(){
00183     delete logh ;
00184     delete logp ;
00185     delete dlpsdlh ;
00186         delete lognb ;
00187         delete dlpsdlnb ;
00188 }
00189 
00190             //------------//
00191             //  Outputs   //
00192             //------------//
00193 
00194 void Eos_tabul::sauve(FILE* fich) const {
00195 
00196     Eos::sauve(fich) ;
00197 
00198         fwrite(tablename, sizeof(char), 160, fich) ;        
00199 
00200 }
00201 
00202             //------------------------//
00203             //  Reading of the table  //
00204             //------------------------//
00205             
00206 void Eos_tabul::read_table() {
00207 
00208   using namespace Unites ;
00209         
00210     ifstream fich(tablename) ;
00211 
00212     if (!fich) {
00213       cout << "Eos_tabul::read_table(): " << endl ;
00214       cout << "Problem in opening the EOS file!" << endl ;
00215       cout << "Aborting..." << endl ;
00216       abort() ;
00217     }
00218 
00219     for (int i=0; i<5; i++) {       //  jump over the file
00220       fich.ignore(1000, '\n') ;             // header
00221         }                                       //
00222 
00223         int nbp ;
00224         fich >> nbp ; fich.ignore(1000, '\n') ;   // number of data
00225     if (nbp<=0) {
00226       cout << "Eos_tabul::read_table(): " << endl ;
00227       cout << "Wrong value for the number of lines!" << endl ;
00228       cout << "nbp = " << nbp << endl ;
00229       cout << "Aborting..." << endl ;
00230       abort() ;
00231     }
00232 
00233     for (int i=0; i<3; i++) {       //  jump over the table
00234       fich.ignore(1000, '\n') ;
00235         }                                      
00236 
00237         press = new double[nbp] ;
00238         nb    = new double[nbp] ;
00239         ro    = new double[nbp] ; 
00240  
00241         logh = new Tbl(nbp) ;
00242         logp = new Tbl(nbp) ;
00243         dlpsdlh = new Tbl(nbp) ;
00244     lognb = new Tbl(nbp) ;
00245         dlpsdlnb = new Tbl(nbp) ;
00246         
00247         logh->set_etat_qcq() ;
00248         logp->set_etat_qcq() ;
00249         dlpsdlh->set_etat_qcq() ;
00250         lognb->set_etat_qcq() ;
00251         dlpsdlnb->set_etat_qcq() ;      
00252         
00253     double rhonuc_cgs = rhonuc_si * 1e-3 ;
00254     double c2_cgs = c_si * c_si * 1e4 ;     
00255         
00256         int no ;
00257         double nb_fm3, rho_cgs, p_cgs ;
00258         
00259         double ww = 0 ;
00260         
00261         for (int i=0; i<nbp; i++) {
00262         
00263             fich >> no ;
00264             fich >> nb_fm3 ;
00265             fich >> rho_cgs ;
00266             fich >> p_cgs ; fich.ignore(1000,'\n') ;            
00267         if ( (nb_fm3<0) || (rho_cgs<0) || (p_cgs < 0) ){
00268           cout << "Eos_tabul::read_table(): " << endl ;
00269           cout << "Negative value in table!" << endl ;
00270           cout << "nb = " << nb_fm3 << ", rho = " << rho_cgs <<
00271             ", p = " << p_cgs << endl ;
00272           cout << "Aborting..." << endl ;
00273           abort() ;
00274         }
00275             double psc2_cgs = p_cgs / c2_cgs ;
00276             double h = log( (rho_cgs + psc2_cgs) /
00277                                 (10 * nb_fm3 * rhonuc_cgs) ) ;
00278             
00279                 press[i] = psc2_cgs ; 
00280                 nb[i]    = nb_fm3 ;
00281                 ro[i]    = rho_cgs ; 
00282 
00283             if (i==0) { ww = h ; }
00284             
00285             h = h - ww + 1.e-14 ;           
00286             
00287 //##
00288 //  cout.precision(8) ;
00289 //  cout << i << "  " << rho_cgs << "  " << p_cgs << "  " << h
00290 //        << "  "  << nb_fm3 << endl ;
00291 //##        
00292             
00293             logh->set(i) = log10( h ) ;
00294             logp->set(i) = log10( psc2_cgs / rhonuc_cgs ) ;
00295             dlpsdlh->set(i) = h * (rho_cgs + psc2_cgs) / psc2_cgs ;
00296         lognb->set(i) = log10(nb_fm3) ;
00297 
00298         }
00299 
00300         double p0, p1, p2, n0, n1, n2, dpdnb; 
00301 
00302     // special case: i=0
00303 
00304           p0 = log(press[0]);
00305           p1 = log(press[1]);
00306           p2 = log(press[2]);
00307 
00308           n0 = log(nb[0]);
00309           n1 = log(nb[1]);
00310           n2 = log(nb[2]);
00311 
00312           dpdnb = p0*(2*n0-n1-n2)/(n0-n1)/(n0-n2) +
00313                      p1*(n0-n2)/(n1-n0)/(n1-n2) +
00314                          p2*(n0-n1)/(n2-n0)/(n2-n1) ;
00315 
00316           dlpsdlnb->set(0) = dpdnb ; 
00317 
00318 /*
00319       cout << exp(n0) << "  " << ro[0] << "  " 
00320                << c2_cgs*exp(p0) << "  " << dpdnb << endl ;
00321 */
00322     for(int i=1;i<nbp-1;i++) { 
00323 
00324           p0 = log(press[i-1]);
00325           p1 = log(press[i]);
00326           p2 = log(press[i+1]);
00327 
00328           n0 = log(nb[i-1]);
00329           n1 = log(nb[i]);
00330           n2 = log(nb[i+1]);
00331 
00332           dpdnb = p0*(n1-n2)/(n0-n1)/(n0-n2) +
00333                      p1*(2*n1-n0-n2)/(n1-n0)/(n1-n2) +
00334                          p2*(n1-n0)/(n2-n0)/(n2-n1) ;
00335 
00336       dlpsdlnb->set(i) = dpdnb ;
00337 
00338 /*
00339       cout << exp(n1) << "  " << ro[i] << "  " 
00340                << c2_cgs*exp(p1) << "  " << dpdnb << endl ;
00341 */
00342 
00343         } 
00344         
00345     // special case: i=nbp-1
00346 
00347           p0 = log(press[nbp-3]);
00348           p1 = log(press[nbp-2]);
00349           p2 = log(press[nbp-1]);
00350 
00351           n0 = log(nb[nbp-3]);
00352           n1 = log(nb[nbp-2]);
00353           n2 = log(nb[nbp-1]);
00354 
00355           dpdnb = p0*(n2-n1)/(n0-n1)/(n0-n2) +
00356                      p1*(n2-n0)/(n1-n0)/(n1-n2) +
00357                          p2*(2*n2-n0-n1)/(n2-n0)/(n2-n1) ;
00358 
00359       dlpsdlnb->set(nbp-1) = dpdnb ;
00360 
00361 /*
00362       cout << exp(n2) << "  " << ro[nbp-1] << "  " 
00363                << c2_cgs*exp(p2) << "  " << dpdnb << endl ;
00364 */
00365                 
00366         hmin = pow( double(10), (*logh)(0) ) ;
00367         hmax = pow( double(10), (*logh)(nbp-1) ) ;
00368         
00369 //##            
00370 //      cout << "logh : " << endl ;
00371 //      cout << *logh << endl ;
00372 //      
00373 //
00374 //      cout << "logp : " << endl ;
00375 //      cout << *logp << endl ;
00376 //      
00377 //      cout << "dlpsdlh : " << endl ;
00378 //      cout << *dlpsdlh << endl ;
00379 //
00380 //      cout << "dlpsdlnb : " << endl ;
00381 //      cout << *dlpsdlnb << endl ;
00382 //
00383 //      
00384 //      cout << "hmin, hmax : " << hmin << "  " << hmax << endl ;
00385 //##        
00386         
00387         fich.close();
00388  
00389         delete [] press ; 
00390         delete [] nb ; 
00391         delete [] ro ; 
00392     
00393 }
00394 
00395 
00396             //------------------------------//
00397             //    Computational routines    //
00398             //------------------------------//
00399 
00400 // Baryon density from enthalpy
00401 //------------------------------
00402 
00403 double Eos_tabul::nbar_ent_p(double ent, const Param* ) const {
00404 
00405     static int i_near = logh->get_taille() / 2 ;
00406 
00407     if ( ent > hmin ) {
00408            if (ent > hmax) {
00409             cout << "Eos_tabul::nbar_ent_p : ent > hmax !" << endl ;
00410             abort() ;
00411            }
00412            double logh0 = log10( ent ) ;
00413            double logp0 ;
00414            double dlpsdlh0 ;
00415            interpol_herm(*logh, *logp, *dlpsdlh, logh0, i_near, logp0,
00416                  dlpsdlh0) ;
00417 
00418            double pp = pow(double(10), logp0) ;
00419 
00420            return pp / ent * dlpsdlh0 * exp(-ent) ;
00421     }
00422     else{
00423     return 0 ;
00424     }
00425 }
00426 
00427 // Energy density from enthalpy
00428 //------------------------------
00429 
00430 double Eos_tabul::ener_ent_p(double ent, const Param* ) const {
00431 
00432     static int i_near = logh->get_taille() / 2 ;
00433 
00434     if ( ent > hmin ) {
00435            if (ent > hmax) {
00436             cout << "Eos_tabul::ener_ent_p : ent > hmax !" << endl ;
00437             abort() ;
00438            }
00439            double logh0 = log10( ent ) ;
00440            double logp0 ;
00441            double dlpsdlh0 ;
00442            interpol_herm(*logh, *logp, *dlpsdlh, logh0, i_near, logp0,
00443                  dlpsdlh0) ;
00444 
00445            double pp = pow(double(10), logp0) ;
00446 
00447            return pp / ent * dlpsdlh0 - pp ;
00448     }
00449     else{
00450     return 0 ;
00451     }
00452 }
00453 
00454 // Pressure from enthalpy
00455 //------------------------
00456 
00457 double Eos_tabul::press_ent_p(double ent, const Param* ) const {
00458 
00459     static int i_near = logh->get_taille() / 2 ;
00460 
00461     if ( ent > hmin ) {
00462            if (ent > hmax) {
00463             cout << "Eos_tabul::press_ent_p : ent > hmax !" << endl ;
00464             abort() ;
00465            }
00466 
00467            double logh0 = log10( ent ) ;
00468            double logp0 ;
00469            double dlpsdlh0 ;
00470            interpol_herm(*logh, *logp, *dlpsdlh, logh0, i_near, logp0,
00471                  dlpsdlh0) ;
00472 
00473            return pow(double(10), logp0) ;
00474     }
00475     else{
00476     return 0 ;
00477     }
00478 }
00479 
00480 // dln(n)/ln(H) from enthalpy 
00481 //---------------------------
00482 
00483 double Eos_tabul::der_nbar_ent_p(double ent, const Param* ) const {
00484 
00485     if ( ent > hmin ) {
00486            if (ent > hmax) {
00487             cout << "Eos_tabul::der_nbar_ent_p : ent > hmax !" << endl ;
00488             abort() ;
00489            }
00490 
00491        double zeta = der_press_ent_p(ent) / der_press_nbar_p(ent) ; 
00492 
00493        return zeta ; 
00494       
00495     }
00496     else 
00497 
00498     return 1./(der_press_nbar_p(hmin)-1.) ;  // to ensure continuity
00499 
00500 }
00501 
00502 
00503 // dln(e)/ln(H) from enthalpy 
00504 //---------------------------
00505 
00506 double Eos_tabul::der_ener_ent_p(double ent, const Param* ) const {
00507 
00508     if ( ent > hmin ) {
00509            if (ent > hmax) {
00510             cout << "Eos_tabul::der_ener_ent_p : ent > hmax !" << endl ;
00511             abort() ;
00512            }
00513 
00514            cout << "Eos_tabul::der_ener_ent_p : not ready yet !" << endl ;
00515            abort() ;
00516        return 0 ;
00517     }
00518     else{
00519     return 0 ;
00520     }
00521 }
00522 
00523 
00524 // dln(p)/ln(H) from enthalpy 
00525 //---------------------------
00526 
00527 double Eos_tabul::der_press_ent_p(double ent, const Param* ) const {
00528 
00529     static int i_near = logh->get_taille() / 2 ;
00530 
00531     if ( ent > hmin ) {
00532            if (ent > hmax) {
00533             cout << "Eos_tabul::der_press_ent_p : ent > hmax !" << endl ;
00534             abort() ;
00535            }
00536 
00537            double logh0 = log10( ent ) ;
00538            double logp0 ;
00539            double dlpsdlh0 ;
00540            interpol_herm(*logh, *logp, *dlpsdlh, logh0, i_near, logp0,
00541                  dlpsdlh0) ;
00542 
00543            return dlpsdlh0 ;
00544 
00545     }
00546     else{
00547         
00548         return 0 ; 
00549     // return der_press_ent_p(hmin) ; // to ensure continuity
00550     }
00551 }
00552 
00553 
00554 // dln(p)/dln(n) from enthalpy 
00555 //---------------------------
00556 
00557 double Eos_tabul::der_press_nbar_p(double ent, const Param*) const {
00558 
00559     static int i_near = logh->get_taille() / 2 ;
00560 
00561     if ( ent > hmin ) {
00562            if (ent > hmax) {
00563             cout << "Eos_tabul::der_press_nbar_p : ent > hmax !" << endl ;
00564             abort() ;
00565            }
00566 
00567            double logh0 = log10( ent ) ;
00568            double dlpsdlnb0 ;
00569 
00570            interpol_linear(*logh, *dlpsdlnb, logh0, i_near, dlpsdlnb0) ;
00571  
00572 //     cout << "gamma: " << dlpsdlnb0 << " ent: " << ent << endl; 
00573 
00574           return dlpsdlnb0 ;
00575 
00576     }
00577     else{
00578         
00579          return (*dlpsdlnb)(0) ; 
00580      // return der_press_nbar_p(hmin) ; // to ensure continuity
00581 
00582     }
00583 
00584           
00585 }

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