![]() |
Equation of State Sub-Library: Version 0.910
|
00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006-2012, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef GEN_SN_EOS_H 00024 #define GEN_SN_EOS_H 00025 00026 #include <cmath> 00027 #include <iostream> 00028 #include <fstream> 00029 #include <o2scl/constants.h> 00030 #include <o2scl/tensor.h> 00031 #include <o2scl/table.h> 00032 #include <o2scl/eff_boson.h> 00033 #include <o2scl/rel_fermion.h> 00034 #include <o2scl/eff_fermion.h> 00035 #include <o2scl/test_mgr.h> 00036 00037 #ifndef DOXYGENP 00038 namespace o2scl { 00039 #endif 00040 00041 /** \brief A base class for the supernova EOSs [abstract] 00042 00043 This class is experimental. 00044 00045 The EOSs are stored in a set of \ref tensor_grid3 objects 00046 on grids with density in \f$ \mathrm{fm}^{-3} 00047 \f$, electron fraction (unitless) and temperature in MeV. 00048 00049 Some tabulated EOSs do not store data for \ref gen_sn_eos::E, 00050 \ref gen_sn_eos::F, \ref gen_sn_eos::S, and \ref gen_sn_eos::P. 00051 In that case, the grid is set for these objects but the data is 00052 set to zero. To compute these from the data after loading the 00053 EOS table, use \ref gen_sn_eos::compute_eg(). 00054 00055 The function \ref load() loads the entire EOS into memory. 00056 Memory allocation is automatically performed by \ref load(), 00057 but not deallocated until free() or the destructor is called. 00058 00059 After loading, you can interpolate the EOS by using 00060 \ref tensor_grid3::interp directly. For example, 00061 the following returns the mass number at an arbitrary 00062 baryon density, electron fraction, and temperature 00063 assuming the table is stored in <tt>skm.dat</tt>: 00064 \verbatim 00065 ls_eos ls; 00066 ls.load("skm.dat"); 00067 double nb=0.01, Ye=0.2, T=10.0; 00068 cout << A.interp(nb,Ye,T) << endl; 00069 \endverbatim 00070 Interpolation for all EOSs is linear by default. 00071 00072 \todo Ensure all chemical potentials are based on the 00073 same rest masses? 00074 \todo Allow logarithmic grids for any of nb, Ye, or T. 00075 \future Create a \ref table object, possibly using 00076 tensor_grid::vector_slice. 00077 \future Show how matrix_slice and vector_slice can be used 00078 with this object. 00079 \future Could this be a child of hadronic_eos_temp and 00080 then directly used in cold_nstar()? 00081 \future Add option to load and store a separate lepton/photon EOS 00082 \future Add muons and/or pions 00083 */ 00084 class gen_sn_eos { 00085 00086 public: 00087 00088 gen_sn_eos(); 00089 00090 virtual ~gen_sn_eos(); 00091 00092 /// \name Load table 00093 //@{ 00094 /// Load table from filename \c fname 00095 virtual void load(std::string fname)=0; 00096 /// If true, a EOS table was successfully loaded (default false) 00097 bool loaded; 00098 /// True if thermodynamics with leptons has been loaded 00099 bool with_leptons_loaded; 00100 /// True if baryon-only thermodynamics has been loaded 00101 bool baryons_only_loaded; 00102 //@} 00103 00104 /// \name Grid and data sizes 00105 //@{ 00106 /// Size of baryon density grid 00107 size_t n_nB; 00108 /// Size of electron fraction grid 00109 size_t n_Ye; 00110 /// Size of temperature grid 00111 size_t n_T; 00112 /// Number of additional data sets 00113 size_t n_oth; 00114 /// Number of base data sets 00115 static const size_t n_base=16; 00116 //@} 00117 00118 /// \name Data 00119 //@{ 00120 /// Total free energy per baryon in MeV 00121 tensor_grid3 F; 00122 /** \brief Free energy per baryon without lepton and photon 00123 contributions in MeV 00124 */ 00125 tensor_grid3 Fint; 00126 /// Total internal energy per baryon in MeV 00127 tensor_grid3 E; 00128 /** \brief Internal energy per baryon without lepton and photon 00129 contributions in MeV 00130 */ 00131 tensor_grid3 Eint; 00132 /// Total pressure in \f$ \mathrm{MeV}/\mathrm{fm}^3 \f$ 00133 tensor_grid3 P; 00134 /** \brief Pressure without lepton and photon contributions 00135 in \f$ \mathrm{MeV}/\mathrm{fm}^3 \f$ 00136 */ 00137 tensor_grid3 Pint; 00138 /// Total entropy per baryon 00139 tensor_grid3 S; 00140 /// Entry per baryon without lepton and photon contributions 00141 tensor_grid3 Sint; 00142 /// Neutron chemical potential in MeV 00143 tensor_grid3 mun; 00144 /// Proton chemical potential in MeV 00145 tensor_grid3 mup; 00146 /// Proton number 00147 tensor_grid3 Z; 00148 /// Mass number 00149 tensor_grid3 A; 00150 /// Neutron fraction 00151 tensor_grid3 Xn; 00152 /// Proton fraction 00153 tensor_grid3 Xp; 00154 /// Alpha particle fraction 00155 tensor_grid3 Xalpha; 00156 /// Fraction of heavy nuclei 00157 tensor_grid3 Xnuclei; 00158 /// Other data sets 00159 tensor_grid3 other[20]; 00160 /// List of pointers to data 00161 tensor_grid3 *arr[n_base+20]; 00162 //@} 00163 00164 /// \name Memory allocation 00165 //@{ 00166 /// Allocate memory 00167 int alloc(); 00168 00169 /// Free allocated memory 00170 int free(); 00171 //@} 00172 00173 /// \name Interpolation 00174 //@{ 00175 /** \brief Set interpolation managers 00176 */ 00177 int set_interp(base_interp_mgr<uvector_base> &bi1, 00178 base_interp_mgr<uvector_const_subvector> &bi2); 00179 00180 /// Default interpolation object 00181 def_interp_mgr<uvector_base,linear_interp> dim1; 00182 00183 /// Default interpolation object 00184 def_interp_mgr<uvector_const_subvector,linear_interp> dim2; 00185 //@} 00186 00187 /// \name Electron and photon contribution 00188 //@{ 00189 /// Photon 00190 boson photon; 00191 /// Electron; 00192 fermion electron; 00193 /// Desc 00194 rel_fermion relf; 00195 /// Desc 00196 eff_boson effb; 00197 /** \brief Compute the electron and photon contribution for the full 00198 grid 00199 00200 This function computes the data for <tt>E, P, S,</tt> and 00201 <tt>F</tt> by adding electrons and photons to the baryon 00202 contributions stored in <tt>Eint, Pint, Sint,</tt> and 00203 <tt>Fint</tt>. 00204 00205 The electron contribution to the internal energy and free 00206 energy computed by this function includes the electron rest 00207 mass. 00208 */ 00209 int compute_eg(); 00210 //@} 00211 00212 /** \brief Test the free energy and store results in \c tm 00213 00214 This checks that the data in \c Fint 00215 is consistent with that in \c Eint and \c Sint . 00216 */ 00217 int check_free_energy(test_mgr &tm); 00218 00219 /** \brief Verbosity parameter (default 1) 00220 */ 00221 int verbose; 00222 00223 /** \brief Compute energy per baryon of matter in beta equilibrium 00224 at zero temperature at a fixed grid point [abstract] 00225 00226 Given an index \c i for the baryon grid, between 0 and \ref 00227 n_nB (inclusive), this computes the properties of matter in 00228 beta equilibrium at zero temperature by finding the electron 00229 fraction grid point which minimizes \ref E. The baryon density 00230 is returned in \c nb, the energy per baryon in \c E_beta, the 00231 pressure in \c P_beta, the electron fraction in \c Ye_beta, 00232 the proton number in \c Z_beta and the mass number of the 00233 nucleus in \c A_beta. 00234 00235 Some tables explicitly contain zero-temperature data which is 00236 used when it is available. Otherwise, linear interpolation is 00237 used to extrapolate down to zero from the lowest temperature 00238 grid points. 00239 */ 00240 virtual void beta_eq_T0(size_t i, double &nb, double &E_beta, 00241 double &P_beta, double &Ye_beta, 00242 double &Z_beta, double &A_beta)=0; 00243 00244 /** \brief Compute properties of matter in beta equilibrium 00245 at s=4 00246 */ 00247 virtual void beta_eq_sfixed(size_t i, double entr, 00248 double &nb, double &E_beta, 00249 double &P_beta, double &Ye_beta, 00250 double &Z_beta, double &A_beta, 00251 double &T_beta) { 00252 00253 if (loaded==false) { 00254 O2SCL_ERR2("No data loaded in ", 00255 "gen_sn_eos::beta_eq_s4().",gsl_einval); 00256 } 00257 if (i>=n_nB) { 00258 O2SCL_ERR2("Too high for baryon grid in ", 00259 "gen_sn_eos::beta_eq_s4().",gsl_einval); 00260 } 00261 if (with_leptons_loaded==false) { 00262 compute_eg(); 00263 } 00264 // Get baryon density from grid 00265 nb=E.get_grid(0,i); 00266 double Emin=0.0; 00267 // The electron fraction grid point corresponding to the minimum 00268 size_t j_found=0; 00269 // The temperature grid point for s=4 at the Ye grid point j_found 00270 size_t k_min_j=0; 00271 for(size_t j=0;j<n_Ye;j++) { 00272 // For each electron fraction, we need to find the 00273 // grid point surrounding s=4 00274 bool found=false; 00275 size_t k_found=0; 00276 for(size_t k=0;k<n_T-1;k++) { 00277 if (S.get(i,j,k)<entr && S.get(i,j,k+1)>=entr) { 00278 k_found=k; 00279 k=n_T; 00280 found=true; 00281 } 00282 } 00283 if (found==false) { 00284 if (S.get(i,j,0)<entr) { 00285 k_found=n_T-2; 00286 } 00287 } 00288 if (j==0) { 00289 Emin=E.get(i,j,k_found)+(E.get(i,j,k_found+1)-E.get(i,j,k_found))* 00290 (entr-S.get(i,j,k_found))/(S.get(i,j,k_found+1)- 00291 S.get(i,j,k_found)); 00292 } else { 00293 double Ethis=E.get(i,j,k_found)+ 00294 (E.get(i,j,k_found+1)-E.get(i,j,k_found))* 00295 (entr-S.get(i,j,k_found))/(S.get(i,j,k_found+1)- 00296 S.get(i,j,k_found)); 00297 if (Ethis<Emin) { 00298 j_found=j; 00299 Emin=Ethis; 00300 k_min_j=k_found; 00301 } 00302 } 00303 } 00304 // Interpolate final results 00305 E_beta=Emin; 00306 P_beta=P.get(i,j_found,k_min_j)+ 00307 (P.get(i,j_found,k_min_j+1)-P.get(i,j_found,k_min_j))* 00308 (entr-S.get(i,j_found,k_min_j))/ 00309 (S.get(i,j_found,k_min_j+1)-S.get(i,j_found,k_min_j)); 00310 Z_beta=Z.get(i,j_found,k_min_j)+ 00311 (Z.get(i,j_found,k_min_j+1)-Z.get(i,j_found,k_min_j))* 00312 (entr-S.get(i,j_found,k_min_j))/ 00313 (S.get(i,j_found,k_min_j+1)-S.get(i,j_found,k_min_j)); 00314 A_beta=A.get(i,j_found,k_min_j)+ 00315 (A.get(i,j_found,k_min_j+1)-A.get(i,j_found,k_min_j))* 00316 (entr-S.get(i,j_found,k_min_j))/ 00317 (S.get(i,j_found,k_min_j+1)-S.get(i,j_found,k_min_j)); 00318 Ye_beta=E.get_grid(1,j_found); 00319 T_beta=E.get_grid(2,k_min_j)+ 00320 (E.get_grid(2,k_min_j+1)-E.get_grid(2,k_min_j))* 00321 (entr-S.get(i,j_found,k_min_j))/ 00322 (S.get(i,j_found,k_min_j+1)-S.get(i,j_found,k_min_j)); 00323 00324 return; 00325 } 00326 00327 }; 00328 00329 /** \brief The Lattimer-Swesty supernova EOS 00330 00331 This class is experimental. 00332 00333 \note \o2 Does not contain the Lattimer-Swesty EOS, only 00334 provides some code to manipulate it. This class is designed 00335 to be used with the files <tt>ls.dat, sk1.dat, ska.dat</tt> 00336 and <tt>skm.dat</tt> as provided on Jim Lattimer's website, 00337 http://www.astro.sunysb.edu/lattimer/EOS/main.html . 00338 00339 Note that the tables on this website are different 00340 than what is generated from the LS Fortran code. See 00341 \ref oo_eos to read O'Connor and Ott's tables 00342 generated from the LS Fortran code. 00343 00344 The four models are 00345 - LS (K=370, Sv=31) 00346 - SKI' (K=371, Sv=30.4) 00347 - SKa (K=263, Sv=34.5) 00348 - SKM* (K=217, Sv=31.4) 00349 00350 \note In the original table, the full internal energy per baryon 00351 (data section 4 of 26) is apparently based on a rest mass of 00352 \f$ Y_e m_p + (1-Y_e) m_n \f$, while the baryon part of the 00353 internal energy per baryon (data section 13 of 26) is based 00354 on a rest mass of \f$ m_n \f$. This means that 00355 \f[ 00356 E - E_{\mathrm{int}} = E_{\mathrm{eg}} - Y_e (m_n - m_p) 00357 \f] 00358 where \f$ E_{\mathrm{eg}} \f$ is the energy per baryon of 00359 electrons and photons. In order to keep things consistent with 00360 the other EOS tables, when the EOS table is loaded, \ref 00361 gen_sn_eos::Eint is rescaled to a rest mass of \f$ Y_e m_p + 00362 (1-Y_e) m_n \f$ . 00363 00364 See also the documentation at \ref gen_sn_eos. 00365 00366 See \ref Lattimer91 and \ref Lattimer85. 00367 00368 \todo There are still a few points for which the electron/photon 00369 EOS seems to be off. 00370 */ 00371 class ls_eos : public gen_sn_eos { 00372 00373 public: 00374 00375 /// Filling factor for nuclei 00376 tensor_grid3 &fill; 00377 /// Baryon number density inside nuclei in \f$ \mathrm{fm}^{-3} \f$ 00378 tensor_grid3 &nb_in; 00379 /// Derivative of pressure with respect to baryon density 00380 tensor_grid3 &dPdn; 00381 /// Derivative of pressure with respect to temperature 00382 tensor_grid3 &dPdT; 00383 /// Derivative of pressure with respect to electron fraction 00384 tensor_grid3 &dPdY; 00385 /// Derivative of entropy with respect to temperature 00386 tensor_grid3 &dsdT; 00387 /// Derivative of entropy with respect to electron fraction 00388 tensor_grid3 &dsdY; 00389 /// Number of neutrons in skin 00390 tensor_grid3 &Nskin; 00391 /// Baryon density outside nuclei in \f$ \mathrm{fm}^{-3} \f$ 00392 tensor_grid3 &nb_out; 00393 /// Proton fraction outside nuclei 00394 tensor_grid3 &x_out; 00395 /** \brief Out of whackness parameter, 00396 \f$ \mu_n-\mu_p-\mu_e+1.293~\mathrm{MeV} \f$, in MeV 00397 */ 00398 tensor_grid3 μ 00399 00400 ls_eos() : 00401 fill(other[0]), 00402 nb_in(other[1]), 00403 dPdn(other[2]), 00404 dPdT(other[3]), 00405 dPdY(other[4]), 00406 dsdT(other[5]), 00407 dsdY(other[6]), 00408 Nskin(other[7]), 00409 nb_out(other[8]), 00410 x_out(other[9]), 00411 mu(other[10]) { 00412 } 00413 00414 /// Load table from filename \c fname 00415 virtual void load(std::string fname); 00416 00417 /** \brief Check electrons and photons 00418 00419 This checks that the electron and photon thermodynamics 00420 generated by \o2 is consistent with the data in 00421 \c E, \c Eint, \c F, \c Fint, \c P, \c Pint, \c S, 00422 and \c Sint. 00423 */ 00424 int check_eg(test_mgr &tm); 00425 00426 /** \brief Compute properties of matter in beta equilibrium 00427 at zero temperature at a baryon density grid point 00428 00429 This EOS table doesn't have T=0 results, so we extrapolate 00430 from the two low-temperature grid points. 00431 */ 00432 virtual void beta_eq_T0(size_t i, double &nb, double &E_beta, 00433 double &P_beta, double &Ye_beta, 00434 double &Z_beta, double &A_beta) { 00435 00436 if (loaded==false || with_leptons_loaded==false) { 00437 O2SCL_ERR2("No data loaded in ", 00438 "ls_eos::beta_eq_T0().",gsl_einval); 00439 } 00440 if (i>=n_nB) { 00441 O2SCL_ERR2("Too high for baryon grid in ", 00442 "ls_eos::beta_eq_T0().",gsl_einval); 00443 } 00444 // Get baryon density from grid 00445 nb=E.get_grid(0,i); 00446 // Find the minima for the two low-temperature grid points 00447 double ET1=E.get(i,0,0), PT1=P.get(i,0,0); 00448 double ET2=E.get(i,0,1), PT2=P.get(i,0,1); 00449 double ZT1=Z.get(i,0,0), ZT2=Z.get(i,0,1); 00450 double AT1=A.get(i,0,0), AT2=A.get(i,0,1); 00451 Ye_beta=E.get_grid(1,0); 00452 for(size_t j=1;j<n_Ye;j++) { 00453 double Enew=E.get(i,j,0); 00454 if (Enew<ET1) { 00455 ET1=Enew; 00456 ET2=E.get(i,j,1); 00457 PT1=P.get(i,j,0); 00458 PT2=P.get(i,j,1); 00459 ZT1=Z.get(i,j,0); 00460 ZT2=Z.get(i,j,1); 00461 AT1=A.get(i,j,0); 00462 AT2=A.get(i,j,1); 00463 Ye_beta=E.get_grid(1,j); 00464 } 00465 } 00466 // Now extrapolate to T=0 00467 E_beta=ET1-(ET2-ET1)*E.get_grid(2,0)/ 00468 (E.get_grid(2,1)-E.get_grid(2,0)); 00469 P_beta=PT1-(PT2-PT1)*P.get_grid(2,0)/ 00470 (P.get_grid(2,1)-P.get_grid(2,0)); 00471 Z_beta=ZT1-(ZT2-ZT1)*Z.get_grid(2,0)/ 00472 (Z.get_grid(2,1)-Z.get_grid(2,0)); 00473 A_beta=AT1-(AT2-AT1)*A.get_grid(2,0)/ 00474 (A.get_grid(2,1)-A.get_grid(2,0)); 00475 00476 return; 00477 } 00478 00479 }; 00480 00481 /** \brief The EOS tables from O'Connor and Ott 00482 00483 This class reads the HDF5 EOS tables generated by E. O'Connor 00484 and C. Ott in \ref OConnor10. The tables are available from 00485 00486 http://stellarcollapse.org/equationofstate 00487 00488 and are available under a creative commons 00489 attribution-noncommercial-share alike license. This \o2 code to 00490 read those tables is licensed (along with all \o2 code) under 00491 the GPLv3 license (with permission from Evan O'Connor). 00492 00493 The original README file from O'Connor and Ott's EOSdriver 00494 code is available in the \o2 00495 distribution in <tt>doc/o2scl/eos/extras/scollapse_README</tt> 00496 and is reproduced below 00497 00498 \verbinclude scollapse_README 00499 */ 00500 class oo_eos : public gen_sn_eos { 00501 00502 public: 00503 00504 oo_eos() : 00505 cs2(other[0]), 00506 dedt(other[1]), 00507 dpderho(other[2]), 00508 dpdrhoe(other[3]), 00509 gamma(other[4]), 00510 mu_e(other[5]), 00511 muhat(other[6]), 00512 munu(other[7]) { 00513 } 00514 00515 /// Speed of sound in cm^2/s^2 00516 tensor_grid &cs2; 00517 /// C_V in erg/g/K 00518 tensor_grid &dedt; 00519 /// dpderho in dyn*g/cm^2/erg 00520 tensor_grid &dpderho; 00521 /// dpdrhoe in dyn cm^3/cm^2/g 00522 tensor_grid &dpdrhoe; 00523 /// Gamma 00524 tensor_grid γ 00525 /// Electron chemical potential per baryon including rest mass 00526 tensor_grid &mu_e; 00527 /// mun - mup 00528 tensor_grid &muhat; 00529 /// mue - mun + mup 00530 tensor_grid &munu; 00531 /// The original mass density grid from the table in g/cm^3 00532 ovector rho; 00533 /// Energy shift for table storage in erg/g 00534 double energy_shift; 00535 00536 /// Use the J. Lattimer et al. method for handling the chemical potentials 00537 static const size_t ls_mode=0; 00538 /// Use the H. Shen et al. method for handling the chemical potentials 00539 static const size_t stos_mode=1; 00540 00541 /// Load table from filename \c fname 00542 virtual void load(std::string fname, size_t mode); 00543 00544 /// Load table from filename \c fname 00545 virtual void load(std::string fname) { 00546 load(fname,0); 00547 return; 00548 } 00549 00550 /** \brief Compute properties of matter in beta equilibrium 00551 at zero temperature at a baryon density grid point 00552 00553 This EOS table doesn't have T=0 results, so we extrapolate 00554 from the two low-temperature grid points. 00555 */ 00556 virtual void beta_eq_T0(size_t i, double &nb, double &E_beta, 00557 double &P_beta, double &Ye_beta, 00558 double &Z_beta, double &A_beta) { 00559 00560 if (loaded==false || baryons_only_loaded==false) { 00561 O2SCL_ERR2("No data loaded in ", 00562 "oo_eos::beta_eq_T0().",gsl_einval); 00563 } 00564 if (i>=n_nB) { 00565 O2SCL_ERR2("Too high for baryon grid in ", 00566 "oo_eos::beta_eq_T0().",gsl_einval); 00567 } 00568 if (with_leptons_loaded==false) { 00569 compute_eg(); 00570 } 00571 // Get baryon density from grid 00572 nb=E.get_grid(0,i); 00573 // Find the minima for the two low-temperature grid points 00574 double ET1=E.get(i,0,0), PT1=P.get(i,0,0); 00575 double ET2=E.get(i,0,1), PT2=P.get(i,0,1); 00576 double ZT1=Z.get(i,0,0), ZT2=Z.get(i,0,1); 00577 double AT1=A.get(i,0,0), AT2=A.get(i,0,1); 00578 Ye_beta=E.get_grid(1,0); 00579 for(size_t j=1;j<n_Ye;j++) { 00580 double Enew=E.get(i,j,0); 00581 if (Enew<ET1) { 00582 ET1=Enew; 00583 ET2=E.get(i,j,1); 00584 PT1=P.get(i,j,0); 00585 PT2=P.get(i,j,1); 00586 ZT1=Z.get(i,j,0); 00587 ZT2=Z.get(i,j,1); 00588 AT1=A.get(i,j,0); 00589 AT2=A.get(i,j,1); 00590 Ye_beta=E.get_grid(1,j); 00591 } 00592 } 00593 // Now extrapolate to T=0 00594 E_beta=ET1-(ET2-ET1)*E.get_grid(2,0)/ 00595 (E.get_grid(2,1)-E.get_grid(2,0)); 00596 P_beta=PT1-(PT2-PT1)*P.get_grid(2,0)/ 00597 (P.get_grid(2,1)-P.get_grid(2,0)); 00598 Z_beta=ZT1-(ZT2-ZT1)*Z.get_grid(2,0)/ 00599 (Z.get_grid(2,1)-Z.get_grid(2,0)); 00600 A_beta=AT1-(AT2-AT1)*A.get_grid(2,0)/ 00601 (A.get_grid(2,1)-A.get_grid(2,0)); 00602 00603 return; 00604 } 00605 00606 }; 00607 00608 /** \brief The Shen et al. supernova EOS 00609 00610 This class is experimental. 00611 00612 \note \o2 Does not contain the EOS, only provides some code to 00613 manipulate it. This class is designed to be used with the file 00614 which was originally called <tt>eos.tab</tt> and now referred to 00615 as <tt>eos1.tab</tt> and stored e.g. at 00616 http://user.numazu-ct.ac.jp/~sumi/eos/. 00617 00618 In order to force the EOS to a uniform grid, linear 00619 interpolation is used to recast the variation in baryon density, 00620 choosing the grid in baryon density to be the same as the 00621 section in the table with T=0.1 MeV and \f$ Y_p = 0.1 \f$ for 00622 all temperature and proton fraction points. 00623 00624 Also, the original EOS is tabulated for constant proton 00625 fraction, and this \o2 interface assumes that the electron 00626 fraction is equal to the proton fraction. Currently, this is a 00627 problem only at higher densities where muons might appear. 00628 00629 The data for \ref gen_sn_eos::E, \ref gen_sn_eos::F, \ref 00630 gen_sn_eos::S, and \ref gen_sn_eos::P is not stored in the table 00631 but can be computed with \ref gen_sn_eos::compute_eg(). 00632 00633 See also the documentation at \ref gen_sn_eos. 00634 00635 See \ref Shen98 and \ref Shen98b . 00636 00637 \note Thanks to Matthias Hempel for providing the correct 00638 temperature grid. 00639 00640 \future Add the T=0 and Ye=0 data to this class 00641 */ 00642 class stos_eos : public gen_sn_eos { 00643 00644 public: 00645 00646 /** \brief Logarithm of baryon number density in 00647 \f$ \mathrm{g}/\mathrm{cm}^3 \f$ 00648 */ 00649 tensor_grid3 &log_rho; 00650 /// Baryon number density in \f$ \mathrm{fm}^{-3} \f$ 00651 tensor_grid3 &nB; 00652 /// Logarithm of proton fraction 00653 tensor_grid3 &log_Y; 00654 /// Proton fraction 00655 tensor_grid3 &Yp; 00656 /// Nucleon effective mass in MeV 00657 tensor_grid3 &M_star; 00658 /// Fraction of quark matter 00659 tensor_grid3 &quark_frac; 00660 00661 stos_eos() : 00662 log_rho(other[0]), 00663 nB(other[1]), 00664 log_Y(other[2]), 00665 Yp(other[3]), 00666 M_star(other[4]), 00667 quark_frac(other[5]) { 00668 } 00669 00670 static const size_t orig_mode=0; 00671 static const size_t quark_mode=1; 00672 00673 /// Load table from filename \c fname 00674 virtual void load(std::string fname) { 00675 load(fname,0); 00676 } 00677 00678 /// Load table from filename \c fname 00679 virtual void load(std::string fname, size_t mode); 00680 00681 /** \brief Compute properties of matter in beta equilibrium 00682 at zero temperature at a baryon density grid point 00683 00684 This EOS table doesn't have T=0 results, so we extrapolate 00685 from the two low-temperature grid points. 00686 */ 00687 virtual void beta_eq_T0(size_t i, double &nb, double &E_beta, 00688 double &P_beta, double &Ye_beta, 00689 double &Z_beta, double &A_beta) { 00690 if (loaded==false || baryons_only_loaded==false) { 00691 O2SCL_ERR2("No data loaded in ", 00692 "stos_eos::beta_eq_T0().",gsl_einval); 00693 } 00694 if (i>=n_nB) { 00695 O2SCL_ERR2("Too high for baryon grid in ", 00696 "stos_eos::beta_eq_T0().",gsl_einval); 00697 } 00698 if (with_leptons_loaded==false) { 00699 compute_eg(); 00700 } 00701 // Get baryon density from grid 00702 nb=E.get_grid(0,i); 00703 // Find the minima for the two low-temperature grid points 00704 double ET1=E.get(i,0,0), PT1=P.get(i,0,0); 00705 double ET2=E.get(i,0,1), PT2=P.get(i,0,1); 00706 double ZT1=Z.get(i,0,0), ZT2=Z.get(i,0,1); 00707 double AT1=A.get(i,0,0), AT2=A.get(i,0,1); 00708 Ye_beta=E.get_grid(1,0); 00709 for(size_t j=1;j<n_Ye;j++) { 00710 double Enew=E.get(i,j,0); 00711 if (Enew<ET1) { 00712 ET1=Enew; 00713 ET2=E.get(i,j,1); 00714 PT1=P.get(i,j,0); 00715 PT2=P.get(i,j,1); 00716 ZT1=Z.get(i,j,0); 00717 ZT2=Z.get(i,j,1); 00718 AT1=A.get(i,j,0); 00719 AT2=A.get(i,j,1); 00720 Ye_beta=E.get_grid(1,j); 00721 } 00722 } 00723 // Now extrapolate to T=0 00724 E_beta=ET1-(ET2-ET1)*E.get_grid(2,0)/ 00725 (E.get_grid(2,1)-E.get_grid(2,0)); 00726 P_beta=PT1-(PT2-PT1)*P.get_grid(2,0)/ 00727 (P.get_grid(2,1)-P.get_grid(2,0)); 00728 Z_beta=ZT1-(ZT2-ZT1)*Z.get_grid(2,0)/ 00729 (Z.get_grid(2,1)-Z.get_grid(2,0)); 00730 A_beta=AT1-(AT2-AT1)*A.get_grid(2,0)/ 00731 (A.get_grid(2,1)-A.get_grid(2,0)); 00732 00733 return; 00734 } 00735 00736 }; 00737 00738 /** \brief A class to manipulate the G. Shen et al. EOS 00739 00740 This class is experimental. 00741 00742 \note \o2 Does not contain the EOS, only provides some code to 00743 manipulate it. This class was designed to be used with the FSU 00744 models given at 00745 http://cecelia.physics.indiana.edu/gang_shen_eos/FSU/fsu.html 00746 and stored in the files named <tt>FSU1.7eos1.01.dat</tt>, 00747 <tt>FSU2.1eos1.01.dat</tt>, <tt>FSU1.7eosb1.01.dat</tt>, and 00748 <tt>FSU2.1eosb1.01.dat</tt> corresponding to 00749 \ref mode_17, \ref mode_21, \ref mode_17b and \ref mode_21b 00750 respectively. 00751 00752 See also the documentation at \ref gen_sn_eos. 00753 00754 \todo More documentation 00755 */ 00756 class sht_eos : public gen_sn_eos { 00757 00758 public: 00759 00760 /// 1.7 solar masses with leptons and photons 00761 static const size_t mode_17=0; 00762 /// 2.1 solar masses with leptons and photons 00763 static const size_t mode_21=1; 00764 /// 1.7 solar masses without leptons and photons 00765 static const size_t mode_17b=2; 00766 /// 2.1 solar masses without leptons and photons 00767 static const size_t mode_21b=3; 00768 /// NL3 model with leptons and photons 00769 static const size_t mode_NL3=4; 00770 /// NL3 model with leptons and photons 00771 static const size_t mode_NL3b=5; 00772 00773 /// Temperature in MeV 00774 tensor_grid3 &T; 00775 /// Proton fraction 00776 tensor_grid3 &Yp; 00777 /// Baryon number density in \f$ 1/\mathrm{fm}^3 \f$ 00778 tensor_grid3 &nB; 00779 /// Electron chemical potential in MeV 00780 tensor_grid3 &mue; 00781 /// Nucleon effective mass in MeV 00782 tensor_grid3 &M_star; 00783 00784 sht_eos() : 00785 T(other[0]), 00786 Yp(other[1]), 00787 nB(other[2]), 00788 mue(other[3]), 00789 M_star(other[4]) { 00790 } 00791 00792 /// Load table from filename \c fname 00793 virtual void load(std::string fname, size_t mode); 00794 00795 /// Load table from filename \c fname 00796 virtual void load(std::string fname) { 00797 load(fname,0); 00798 } 00799 00800 /** \brief Compute properties of matter in beta equilibrium 00801 at zero temperature at a baryon density grid point 00802 */ 00803 virtual void beta_eq_T0(size_t i, double &nb, double &E_beta, 00804 double &P_beta, double &Ye_beta, 00805 double &Z_beta, double &A_beta) { 00806 if (loaded==false || (with_leptons_loaded==false && 00807 baryons_only_loaded==false)) { 00808 O2SCL_ERR2("No data loaded in ", 00809 "sht_eos::beta_eq_T0().",gsl_einval); 00810 } 00811 if (i>=n_nB) { 00812 O2SCL_ERR2("Too high for baryon grid in ", 00813 "sht_eos::beta_eq_T0().",gsl_einval); 00814 } 00815 if (with_leptons_loaded==false) { 00816 compute_eg(); 00817 } 00818 // Look up baryon density from grid 00819 nb=E.get_grid(0,i); 00820 // Minimize over all electron fractions 00821 E_beta=E.get(i,0,0); 00822 P_beta=P.get(i,0,0); 00823 Ye_beta=E.get_grid(1,0); 00824 Z_beta=Z.get(i,0,0); 00825 A_beta=A.get(i,0,0); 00826 for(size_t j=1;j<n_Ye;j++) { 00827 double Enew=E.get(i,j,0); 00828 if (Enew<E_beta) { 00829 E_beta=Enew; 00830 P_beta=P.get(i,j,0); 00831 Ye_beta=E.get_grid(1,j); 00832 Z_beta=Z.get(i,j,0); 00833 A_beta=A.get(i,j,0); 00834 } 00835 } 00836 return; 00837 } 00838 00839 }; 00840 00841 /** \brief The Hempel et al. supernova EOSs 00842 00843 This class is experimental. 00844 00845 \note \o2 Does not contain the EOS, only provides some code to 00846 manipulate it. This class was designed to be used with the files 00847 <tt>dd2_frdm_eos_shen98format_v1.02.tab</tt>, 00848 <tt>fsg_roca_eos_shen98format_v1.0.tab</tt>, and 00849 <tt>nl3_lala_eos_shen98format_v1.0.tab</tt> as obtained from 00850 http://phys-merger.physik.unibas.ch/~hempel/eos.html. 00851 00852 See also the documentation at \ref gen_sn_eos. 00853 00854 See \ref Hempel10 and \ref Hempel11. 00855 */ 00856 class hfsl_eos : public gen_sn_eos { 00857 00858 public: 00859 00860 /** \brief Logarithm of baryon number density in 00861 \f$ \mathrm{g}/\mathrm{cm}^3 \f$ 00862 */ 00863 tensor_grid3 &log_rho; 00864 /// Baryon number density in \f$ 1/\mathrm{fm}^3 \f$ 00865 tensor_grid3 &nB; 00866 /// Logarithm of proton fraction 00867 tensor_grid3 &log_Y; 00868 /// Proton fraction 00869 tensor_grid3 &Yp; 00870 /// Nucleon effective mass in MeV 00871 tensor_grid3 &M_star; 00872 /// Mass number of light fragments 00873 tensor_grid3 &A_light; 00874 /// Proton number of light fragments 00875 tensor_grid3 &Z_light; 00876 00877 hfsl_eos() : 00878 log_rho(other[0]), 00879 nB(other[1]), 00880 log_Y(other[2]), 00881 Yp(other[3]), 00882 M_star(other[4]), 00883 A_light(other[5]), 00884 Z_light(other[6]) { 00885 } 00886 00887 /// Load table from filename \c fname 00888 virtual void load(std::string fname); 00889 00890 /** \brief Compute properties of matter in beta equilibrium 00891 at zero temperature at a baryon density grid point 00892 */ 00893 virtual void beta_eq_T0(size_t i, double &nb, double &E_beta, 00894 double &P_beta, double &Ye_beta, 00895 double &Z_beta, double &A_beta) { 00896 if (loaded==false || baryons_only_loaded==false) { 00897 O2SCL_ERR2("No data loaded in ", 00898 "hfsl_eos::beta_eq_T0().",gsl_einval); 00899 } 00900 if (i>=n_nB) { 00901 O2SCL_ERR2("Too high for baryon grid in ", 00902 "hfsl_eos::beta_eq_T0().",gsl_einval); 00903 } 00904 if (with_leptons_loaded==false) { 00905 compute_eg(); 00906 } 00907 // Get baryon density from grid 00908 nb=E.get_grid(0,i); 00909 // Find the minima for the two low-temperature grid points 00910 double ET1=E.get(i,0,0), PT1=P.get(i,0,0); 00911 double ET2=E.get(i,0,1), PT2=P.get(i,0,1); 00912 double ZT1=Z.get(i,0,0), ZT2=Z.get(i,0,1); 00913 double AT1=A.get(i,0,0), AT2=A.get(i,0,1); 00914 Ye_beta=E.get_grid(1,0); 00915 for(size_t j=1;j<n_Ye;j++) { 00916 double Enew=E.get(i,j,0); 00917 if (Enew<ET1) { 00918 ET1=Enew; 00919 ET2=E.get(i,j,1); 00920 PT1=P.get(i,j,0); 00921 PT2=P.get(i,j,1); 00922 ZT1=Z.get(i,j,0); 00923 ZT2=Z.get(i,j,1); 00924 AT1=A.get(i,j,0); 00925 AT2=A.get(i,j,1); 00926 Ye_beta=E.get_grid(1,j); 00927 } 00928 } 00929 // Now extrapolate to T=0 00930 E_beta=ET1-(ET2-ET1)*E.get_grid(2,0)/ 00931 (E.get_grid(2,1)-E.get_grid(2,0)); 00932 P_beta=PT1-(PT2-PT1)*P.get_grid(2,0)/ 00933 (P.get_grid(2,1)-P.get_grid(2,0)); 00934 Z_beta=ZT1-(ZT2-ZT1)*Z.get_grid(2,0)/ 00935 (Z.get_grid(2,1)-Z.get_grid(2,0)); 00936 A_beta=AT1-(AT2-AT1)*A.get_grid(2,0)/ 00937 (A.get_grid(2,1)-A.get_grid(2,0)); 00938 00939 return; 00940 } 00941 00942 }; 00943 00944 #ifndef DOXYGENP 00945 } 00946 #endif 00947 00948 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).