LORENE
eos_compose_fit.C
1 /*
2  * Methods of class Eos_compose_fit
3  *
4  * (see file eos_compose_fit.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2022 Jerome Novak
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License as published by
15  * the Free Software Foundation; either version 2 of the License, or
16  * (at your option) any later version.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 
30 
31 
32 /*
33  * $Id: eos_compose_fit.C,v 1.8 2023/07/12 15:39:30 j_novak Exp $
34  * $Log: eos_compose_fit.C,v $
35  * Revision 1.8 2023/07/12 15:39:30 j_novak
36  * Reversed order of matching the EoS fit: it now starts from the crust (polytrope) and adjusts the kappa of this crust so as to have the right pressure at the highest density in the original table.
37  *
38  * Revision 1.5 2023/03/17 08:36:59 j_novak
39  * Output of "middle" enthalpy (boundary between domains).
40  *
41  * Revision 1.4 2023/01/27 16:10:35 j_novak
42  * A polytrope (class Eos_poly) is used for low and high enthalpies.
43  *
44  * Revision 1.3 2022/07/21 12:34:18 j_novak
45  * Corrected units
46  *
47  * Revision 1.2 2022/07/20 12:59:43 j_novak
48  * Added methods for saving into files and constructor from a file
49  *
50  *
51  * $Header: /cvsroot/Lorene/C++/Source/Eos/eos_compose_fit.C,v 1.8 2023/07/12 15:39:30 j_novak Exp $
52  *
53  */
54 
55 // Headers Lorene
56 #include "eos_compose_fit.h"
57 #include "scalar.h"
58 #include "utilitaires.h"
59 #include "unites.h"
60 
61 
62 namespace Lorene {
63 
64 
65  //----------------------------//
66  // Constructors //
67  //----------------------------//
68 
69 // Standard constructor
70 // --------------------
71  Eos_compose_fit::Eos_compose_fit(const string& par_file_name) :
72  Eos("EoS fitted from CompOSE data"), p_eos_low(nullptr), p_eos_high(nullptr),
73  mg(nullptr), mp(nullptr), log_p(nullptr), log_e(nullptr), log_nb(nullptr),
74  log_cs2(nullptr)
75  {
76  ifstream para_file(par_file_name) ;
77  para_file.ignore(1000, '\n') ;
78  para_file.ignore(1000, '\n') ;
79  read_and_compute(para_file) ;
80 
81  }
82 
83 
84 // Constructor from binary file
85 // ----------------------------
87  Eos(fich), p_eos_low(nullptr), p_eos_high(nullptr), mg(nullptr),
88  mp(nullptr), log_p(nullptr), log_e(nullptr), log_nb(nullptr), log_cs2(nullptr)
89  {
90 
91  char tmp_string[160] ;
92  fread(tmp_string, sizeof(char), 160, fich) ;
93  tablename = tmp_string ;
94 
95  fread_be(&n_coefs, sizeof(int), 1, fich) ;
96  fread_be(&nb_min, sizeof(double), 1, fich) ;
97  fread_be(&nb_mid, sizeof(double), 1, fich) ;
98  fread_be(&nb_max, sizeof(double), 1, fich) ;
99  fread_be(&hmin, sizeof(double), 1, fich) ;
100  fread_be(&hmax, sizeof(double), 1, fich) ;
101  double gamma_low, gamma_high, kappa_low, kappa_high;
102  double m0_low, m0_high, mu0_low, mu0_high ;
103  fread_be(&gamma_low, sizeof(double), 1, fich) ;
104  fread_be(&gamma_high, sizeof(double), 1, fich) ;
105  fread_be(&kappa_low, sizeof(double), 1, fich) ;
106  fread_be(&kappa_high, sizeof(double), 1, fich) ;
107  fread_be(&m0_low, sizeof(double), 1, fich) ;
108  fread_be(&m0_high, sizeof(double), 1, fich) ;
109  fread_be(&mu0_low, sizeof(double), 1, fich) ;
110  fread_be(&mu0_high, sizeof(double), 1, fich) ;
111 
112  p_eos_low = new Eos_poly(gamma_low, kappa_low, m0_low, mu0_low) ;
113  p_eos_high = new Eos_poly(gamma_high, kappa_high, m0_high, mu0_high) ;
114  mg = new Mg3d(fich) ;
115  mp = new Map_af(*mg, fich) ;
116 
117  log_p = new Scalar(*mp, *mg, fich) ;
118  log_e = new Scalar(*mp, *mg, fich) ;
119  log_nb = new Scalar(*mp, *mg, fich) ;
120  log_cs2 = new Scalar(*mp, *mg, fich) ;
121 
122  }
123 
124 
125 
126 // Constructor from a formatted file
127 // ---------------------------------
128 
129  Eos_compose_fit::Eos_compose_fit(ifstream& para_file) :
130  Eos(para_file), p_eos_low(nullptr), p_eos_high(nullptr),
131  mg(nullptr), mp(nullptr), log_p(nullptr), log_e(nullptr), log_nb(nullptr),
132  log_cs2(nullptr)
133  {
134  read_and_compute(para_file) ;
135  }
136 
137 
138  //--------------//
139  // Destructor //
140  //--------------//
141 
143  if (p_eos_low != nullptr) delete p_eos_low ;
144  if (p_eos_high != nullptr) delete p_eos_high ;
145  if (mg != nullptr) delete mg ;
146  if (mp != nullptr) delete mp ;
147  if (log_p != nullptr) delete log_p ;
148  if (log_e != nullptr) delete log_e ;
149  if (log_nb != nullptr) delete log_nb ;
150  if (log_cs2 != nullptr) delete log_cs2 ;
151 }
152 
153  //------------------------//
154  // Comparison operators //
155  //------------------------//
156 
157 
158  bool Eos_compose_fit::operator==(const Eos& eos_i) const {
159 
160  bool resu = true ;
161 
162  if ( eos_i.identify() != identify() ) {
163  cout << "The second EOS is not of type Eos_compose_fit !" << endl ;
164  resu = false ;
165  }
166  return resu ;
167  }
168 
169  bool Eos_compose_fit::operator!=(const Eos& eos_i) const {
170 
171  return !(operator==(eos_i)) ;
172  }
173 
174  //------------//
175  // Outputs //
176  //------------//
177 
178 void Eos_compose_fit::sauve(FILE* fich) const {
179 
180  Eos::sauve(fich) ;
181 
182  char tmp_string[160] ;
183  strcpy(tmp_string, tablename.c_str()) ;
184  fwrite(tmp_string, sizeof(char), 160, fich) ;
185  fwrite_be(&n_coefs, sizeof(int), 1, fich) ;
186  fwrite_be(&nb_min, sizeof(double), 1, fich) ;
187  fwrite_be(&nb_mid, sizeof(double), 1, fich) ;
188  fwrite_be(&nb_max, sizeof(double), 1, fich) ;
189  fwrite_be(&hmin, sizeof(double), 1, fich) ;
190  fwrite_be(&hmax, sizeof(double), 1, fich) ;
191  double gamma_low = p_eos_low->get_gam() ;
192  double kappa_low = p_eos_low->get_kap() ;
193  double m0_low = p_eos_low->get_m_0() ;
194  double mu0_low = p_eos_low->get_mu_0() ;
195  double gamma_high = p_eos_high->get_gam() ;
196  double kappa_high = p_eos_high->get_kap() ;
197  double m0_high = p_eos_high->get_m_0() ;
198  double mu0_high = p_eos_high->get_mu_0() ;
199  fwrite_be(&gamma_low, sizeof(double), 1, fich) ;
200  fwrite_be(&gamma_high, sizeof(double), 1, fich) ;
201  fwrite_be(&kappa_low, sizeof(double), 1, fich) ;
202  fwrite_be(&kappa_high, sizeof(double), 1, fich) ;
203  fwrite_be(&m0_low, sizeof(double), 1, fich) ;
204  fwrite_be(&m0_high, sizeof(double), 1, fich) ;
205  fwrite_be(&mu0_low, sizeof(double), 1, fich) ;
206  fwrite_be(&mu0_high, sizeof(double), 1, fich) ;
207 
208  mg->sauve(fich) ;
209  mp->sauve(fich) ;
210 
211  log_p->sauve(fich) ;
212  log_e->sauve(fich) ;
213  log_nb->sauve(fich) ;
214  log_cs2->sauve(fich) ;
215 
216 }
217 
218  ostream& Eos_compose_fit::operator>>(ostream & ost) const {
219 
220  ost << "EOS of class Eos_compose_fit." << endl ;
221  ost << "Built from files in directory " << tablename << endl ;
222  ost << "Number of coefficients used for the polynomial fit : "
223  << n_coefs << endl ;
224  ost << "nb_min = " << nb_min << ", nb_mid = " << nb_mid
225  << ", nb_max = " << nb_max << " [fm^-3]" << endl ;
226  ost << "hmin = " << hmin << ", hmid = " << exp(mp->val_r_jk(0, 1., 0, 0))
227  << ", hmax = " << hmax << endl ;
228  ost << "EoS for low density part: " << *p_eos_low ;
229  ost << "EoS for high density part: " << *p_eos_high ;
230  ost << *mg << endl ;
231  return ost ;
232 
233 }
234 
235  //----------------------------------------------//
236  // Reading of the table and computing the fit //
237  //----------------------------------------------//
238 
239 void Eos_compose_fit::read_and_compute(ifstream& para_file) {
240 
241  cout << para_file.good() << endl ;
242  para_file >> tablename ;
243  para_file >> nb_mid >> nb_min ;
244  para_file.ignore(1000, '\n') ;
245  assert(nb_min < nb_mid) ;
246  para_file >> n_coefs ;
247  para_file.ignore(1000, '\n') ;
248  int nr ; para_file >> nr ;
249  para_file.ignore(1000, '\n') ;
250 
251  Tbl* logh_read = nullptr ;
252  Tbl* logp_read = nullptr ;
253  Tbl* loge_read = nullptr ;
254  Tbl* lognb_read = nullptr ;
255  Tbl* dlpsdlnb = nullptr ;
256  int nbp = 0 ;
257 
258  read_compose_data(nbp, logh_read, logp_read, loge_read, lognb_read, dlpsdlnb) ;
259 
260  int i_min = 0 ; int i_mid = 0 ;
261 
262  Tbl nb_read = exp((*lognb_read)) ;
263 
264  for (int i=0; i<nbp; i++) {
265  if (nb_read(i) < 10.*nb_min) i_min = i ; // nb_min & nb_mid are in fm^-3
266  if (nb_read(i) < 10.*nb_mid) i_mid = i ; // Lorene units are 0.1 fm^-3
267  }
268 
269  int i_max = nbp - 1 ; //## to improve, depending on the sound speed values
270 
271  int nz = 2 ;
272  double x_min = (*logh_read)(i_min) ;
273  double x_mid = (*logh_read)(i_mid) ;
274  double x_max = (*logh_read)(i_max) ;
275 
276  hmin = exp(x_min) ;
277  hmax = exp(x_max) ;
278 
279  // Boundaries of spectral grid
280  //------------------------------
281  Tbl xtab(nz+1) ; xtab.set_etat_qcq() ;
282  xtab.set(0) = x_min ;
283  xtab.set(1) = x_mid ;
284  xtab.set(nz) = x_max ;
285 
286  // Grid and xx = log(h) coordinate
287  mg = new Mg3d(nz, nr, 1, 1, SYM, SYM) ;
288  mp = new Map_af(*mg, xtab) ;
289 
290  adiabatic_index_fit(i_min, i_mid, *logh_read, *logp_read, *loge_read, *lognb_read,
291  *dlpsdlnb) ;
292 
293  // Cleaning
294  if (logh_read != nullptr) delete logh_read ;
295  if (logp_read != nullptr) delete logp_read ;
296  if (loge_read != nullptr) delete loge_read ;
297  if (lognb_read != nullptr) delete lognb_read ;
298  if (dlpsdlnb != nullptr) delete dlpsdlnb ;
299 }
300 
301 void Eos_compose_fit::write_lorene_table(const string& name_file, int nlines)
302  const {
303 
304  cout << "Writing the Eos_compose_fit object into a Lorene-format file ("
305  << name_file << ") ..." << flush ;
306  if (nlines < 10) {
307  cerr << "Eos_compose_fit::write_lorene_table" << endl ;
308  cerr << "The number of lines to be outputted is too small!" << endl ;
309  cerr << " nlines = " << nlines << endl ;
310  cerr << "Aborting..." << endl ;
311  abort() ;
312  }
313  double rhonuc_cgs = Unites::rhonuc_si * 1.e-3 ;
314  double c2_cgs = Unites::c_si * Unites::c_si * 1.e4 ;
315 
316  ofstream lorene_file(name_file) ;
317  Scalar press = exp(*log_p) * c2_cgs * rhonuc_cgs ;
318  Scalar ener = exp(*log_e) * rhonuc_cgs ;
319  Scalar nbar = exp(*log_nb) * 10. ;
320 
321  lorene_file << "#" << endl ;
322  lorene_file << "# Built from Eos_compose_fit object" << endl ;
323  lorene_file << "# From the table: " << tablename << endl ;
324  lorene_file << "#\n#" << endl ;
325  lorene_file << nlines << "\t Number of lines" << endl ;
326  lorene_file << "#\n#\t n_B [fm^{-3}] rho [g/cm^3] p [dyn/cm^2]" << endl ;
327  lorene_file << "#" << endl ;
328 
329  // const Coord& xx = mp->r ;
330 #ifndef NDEBUG
331  int nz = mg->get_nzone() ;
332  assert (nz > 1) ;
333 #endif
334  double xmin0 = log(1.e-14) ;
335  double xmax0 = log(hmin) ;
336  int nlines0 = nlines / 10 ;
337  double dx = (xmax0 - xmin0)/double(nlines0-1) ;
338  assert(dx>0.) ;
339  double logh = xmin0 ;
340  for (int i=0; i<nlines0; i++) {
341  double ent = exp(logh) ;
342  lorene_file << setprecision(16) ;
343  lorene_file << i << '\t' << nbar_ent_p(ent)/10. << '\t'
344  << ener_ent_p(ent)*rhonuc_cgs / c2_cgs << '\t'
345  << press_ent_p(ent) *c2_cgs * rhonuc_cgs << endl ;
346  logh += dx ;
347  }
348  logh -= dx ;
349  xmin0 = xmax0 ;
350  xmax0 = log(10.*hmax) ; // ## to improve?
351  dx = (xmax0 - xmin0) / double(nlines - nlines0-1) ;
352  for (int i=nlines0; i<nlines; i++) {
353  double ent = exp(logh) ;
354  lorene_file << setprecision(16) ;
355  lorene_file << i << '\t' << nbar_ent_p(ent)/10. << '\t'
356  << ener_ent_p(ent)*rhonuc_cgs / c2_cgs << '\t'
357  << press_ent_p(ent) *c2_cgs * rhonuc_cgs << endl ;
358  logh += dx ;
359  }
360  lorene_file.close() ;
361  cout << " done!" << endl ;
362 }
363 
364 
365 
366 
367  //------------------------------//
368  // Computational routines //
369  //------------------------------//
370 
371 // Baryon density from enthalpy
372 //------------------------------
373 
374 double Eos_compose_fit::nbar_ent_p(double ent, const Param* ) const {
375 
376  if ( ent > hmin ) {
377  if (ent < hmax) {
378  double logh0 = log( ent ) ;
379  return exp(log_nb->val_point(logh0, 0. ,0.)) ;
380  }
381  else {
382  assert(p_eos_high != nullptr) ;
383  return p_eos_high->nbar_ent_p(ent) ;
384  }
385  }
386  else{
387  assert(p_eos_low != nullptr) ;
388  return p_eos_low->nbar_ent_p(ent) ;
389  }
390 }
391 
392 // Energy density from enthalpy
393 //------------------------------
394 
395 double Eos_compose_fit::ener_ent_p(double ent, const Param* ) const {
396 
397  if ( ent > hmin ) {
398  if (ent < hmax) {
399  double logh0 = log( ent ) ;
400  return exp(log_e->val_point(logh0, 0., 0.)) ;
401  }
402  else {
403  assert(p_eos_high != nullptr) ;
404  return p_eos_high->ener_ent_p(ent) ;
405  }
406  }
407  else{
408  assert(p_eos_low != nullptr) ;
409  return p_eos_low->ener_ent_p(ent) ;
410  }
411 }
412 
413 // Pressure from enthalpy
414 //------------------------
415 
416 double Eos_compose_fit::press_ent_p(double ent, const Param* ) const {
417 
418  if ( ent > hmin ) {
419  if (ent < hmax) {
420  double logh0 = log( ent ) ;
421  return exp(log_p->val_point(logh0, 0., 0.)) ;
422  }
423  else {
424  assert(p_eos_high != nullptr) ;
425  return p_eos_high->press_ent_p(ent) ;
426  }
427  }
428  else{
429  assert(p_eos_low != nullptr) ;
430  return p_eos_low->press_ent_p(ent) ;
431  }
432 }
433 
434 // dln(n)/ln(H) from enthalpy
435 //---------------------------
436 
437 double Eos_compose_fit::der_nbar_ent_p(double ent, const Param* ) const {
438 
439  return ent / csound_square_ent_p(ent) ;
440 }
441 
442 
443 // dln(e)/ln(H) from enthalpy
444 //---------------------------
445 
446 double Eos_compose_fit::der_ener_ent_p(double ent, const Param* ) const {
447 
448  double ener = ener_ent_p(ent) ;
449  double press = press_ent_p(ent) ;
450  double cs2 = csound_square_ent_p(ent) ;
451 
452  return (ener + press)*ent / (ener * cs2) ;
453 }
454 
455 
456 // dln(p)/ln(H) from enthalpy
457 //---------------------------
458 
459 double Eos_compose_fit::der_press_ent_p(double ent, const Param* ) const {
460 
461  double ener = ener_ent_p(ent) ;
462  double press = press_ent_p(ent) ;
463 
464  return ent*(ener + press)/press ;
465 
466 }
467 
468 
469 // dln(p)/dln(n) from enthalpy
470 //---------------------------
471 
472 double Eos_compose_fit::der_press_nbar_p(double ent, const Param*) const {
473 
474  double dlpsdlh0 = der_press_ent_p(ent) ;
475  double dlnsdlh0 = der_nbar_ent_p(ent) ;
476 
477  return dlpsdlh0 / dlnsdlh0 ;
478 
479 }
480 
481 double Eos_compose_fit::csound_square_ent_p(double ent, const Param*) const {
482 
483  if ( ent > hmin ) {
484  if (ent < hmax) {
485  double logh0 = log( ent ) ;
486  return exp(log_cs2->val_point(logh0, 0., 0.)) ;
487  }
488  else {
489  assert(p_eos_high != nullptr) ;
490  return p_eos_high->csound_square_ent_p(ent) ;
491  }
492  }
493  else{
494  assert(p_eos_low != nullptr) ;
495  return p_eos_low->csound_square_ent_p(ent) ;
496  }
497 }
498 
499 } // end of namespace Lorene
virtual double ener_ent_p(double ent, const Param *par=0x0) const
Computes the total energy density from the log-enthalpy.
Definition: eos_poly.C:398
virtual double nbar_ent_p(double ent, const Param *par=0x0) const
Higher bound in density.
void read_compose_data(int &nbp, Tbl *&logh, Tbl *&logp, Tbl *&loge, Tbl *&lognb, Tbl *&gam1)
Reads Compose data and stores the values of thermodynamic quantities.
virtual void sauve(FILE *) const
Save in a file.
Definition: map_af.C:616
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
Eos_compose_fit(const string &param_file)
Constructor from a parameter file.
virtual double der_press_ent_p(double ent, const Param *par=0x0) const
Computes the logarithmic derivative from the log-enthalpy.
virtual void sauve(FILE *) const
Save in a file.
const Mg3d * mg
Multi-grid defining the number of domains and spectral coefficients.
Lorene prototypes.
Definition: app_hor.h:67
Equation of state base class.
Definition: eos.h:209
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
virtual int identify() const =0
Returns a number to identify the sub-classe of Eos the object belongs to.
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
virtual double der_nbar_ent_p(double ent, const Param *par=0x0) const
Computes the logarithmic derivative from the log-enthalpy.
virtual bool operator!=(const Eos &) const
Comparison operator (difference)
string tablename
Name of the file containing the tabulated data.
virtual double press_ent_p(double ent, const Param *par=0x0) const
Computes the pressure from the log-enthalpy.
virtual double der_press_nbar_p(double ent, const Param *par=0x0) const
Computes the logarithmic derivative from the log-enthalpy.
virtual double csound_square_ent_p(double, const Param *par=0x0) const
Computes the sound speed squared from the enthapy with extra parameters (virtual function implemente...
virtual double val_r_jk(int l, double xi, int j, int k) const
Returns the value of the radial coordinate r for a given and a given collocation point in in a give...
double get_gam() const
Returns the adiabatic index (cf. Eq. (3))
Definition: eos_poly.C:271
void read_and_compute(ifstream &)
Reads the Compose data and makes the fit.
virtual void sauve(FILE *) const
Save in a file.
Definition: eos.C:189
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
virtual double ener_ent_p(double ent, const Param *par=0x0) const
Computes the total energy density from the log-enthalpy.
void adiabatic_index_fit(int i_min, int i_max, const Tbl &logh_read, const Tbl &logp_read, const Tbl &loge_read, const Tbl &lognb_read, const Tbl &gam1_read)
From the read values, makes the fit on the adiabatic index and deduces the other quantities from ther...
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
Scalar * log_cs2
Table of .
void sauve(FILE *fd, bool save_base=false) const
Saves into a file.
Definition: mg3d.C:500
virtual ostream & operator>>(ostream &) const
Operator >>
double get_kap() const
Returns the pressure coefficient (cf.
Definition: eos_poly.C:275
double get_m_0() const
Return the individual particule mass (cf.
Definition: eos_poly.C:279
Scalar * log_nb
Table of .
Scalar * log_e
Table of .
int n_coefs
Number of coeficients for polynomial regression.
const Map_af * mp
Mapping in .
double nb_max
Higher bound on the density, above which the EoS is assumed to be a polytrope.
virtual double nbar_ent_p(double ent, const Param *par=0x0) const
Computes the baryon density from the log-enthalpy.
Definition: eos_poly.C:383
Parameter storage.
Definition: param.h:125
Polytropic equation of state (relativistic case).
Definition: eos.h:812
const Eos_poly * p_eos_low
Pointer on a polytropic EoS for the description of low densities (nb<nb_min)
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
virtual double csound_square_ent_p(double ent, const Param *par=0x0) const
Computes the sound speed squared from the enthapy with extra parameters (virtual function implemente...
Definition: eos_poly.C:469
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition: fwrite_be.C:73
virtual double press_ent_p(double ent, const Param *par=0x0) const
Computes the pressure from the log-enthalpy.
Definition: eos_poly.C:410
double get_mu_0() const
Return the relativistic chemical potential at zero pressure [unit: , with ].
Definition: eos_poly.C:283
const Eos_poly * p_eos_high
Pointer on a polytropic EoS for the description of high densities (nb>nb_max)
virtual double der_ener_ent_p(double ent, const Param *par=0x0) const
Computes the logarithmic derivative from the log-enthalpy.
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
Definition: scalar.C:896
double nb_min
Lower bound in baryon density, below which the EoS is assumed to be a polytrope.
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition: fread_be.C:72
Multi-domain grid.
Definition: grilles.h:279
Affine radial mapping.
Definition: map.h:2042
double hmin
Values of enthalpy corresponding to nb_min and nb_max.
double nb_mid
Middle bound in baryon density, above which which the EoS is determined from the polynomial fit to th...
Scalar * log_p
Table of .
void write_lorene_table(const string &, int nlines=200) const
Save into a table in Lorene format.
virtual bool operator==(const Eos &) const
Comparison operator (egality)
Basic array class.
Definition: tbl.h:164
virtual int identify() const
Returns a number to identify the sub-classe of Eos the object belongs to.
virtual ~Eos_compose_fit()
Destructor.