53 n_lim(n_lim0), m_n(m_n0) {
55 cerr <<
"Deprecated constructor, please do not use ! Aborting..." << endl; abort() ;
56 set_name(
"Pseudo-polytropic fit of cold, beta-equilibrated EoS") ;
59 coefs =
new Tbl(coefs0) ;
60 n_coefs = coefs->get_taille() ;
62 double x_lim =
log(n_lim*0.1) ;
63 double alpha = (*coefs)(n_coefs-1) ;
64 double sum_poly = 0., sum_der_poly = 0. ;
65 for (
int i=0; i<n_coefs-1; i++){
66 sum_poly += (*coefs)(i)*
pow(x_lim,
double(i)) ;
67 sum_der_poly += (i == 0) ? 0. : (*coefs)(i)*
double(i)*
pow(x_lim,
double(i-1)) ;
69 ent_lim = log1p(mu_0 +
exp(alpha*x_lim)*((alpha+1)*sum_poly + sum_der_poly)) ;
77 n_lim(eosi.n_lim), ent_lim(eosi.ent_lim), m_n(eosi.m_n), mu_0(eosi.mu_0), n_coefs(eosi.n_coefs) {}
84 fread_be(&n_coefs,
sizeof(
int), 1, fich) ;
85 coefs =
new Tbl(fich) ;
86 fread_be(&kappa_low,
sizeof(
double), 1, fich) ;
87 fread_be(&gamma_low,
sizeof(
double), 1, fich) ;
88 fread_be(&n_lim,
sizeof(
double), 1, fich) ;
89 fread_be(&ent_lim,
sizeof(
double), 1, fich) ;
90 fread_be(&m_n,
sizeof(
double), 1, fich) ;
91 fread_be(&mu_0,
sizeof(
double), 1, fich) ;
92 fread_be(&Lambda,
sizeof(
double), 1, fich) ;
94 eos_low =
new Eos_poly(gamma_low, kappa_low) ;
101 fich >> n_coefs ; fich.ignore(1000,
'\n') ;
104 for (
int i=0; i < n_coefs; i++)
105 fich >> coefs->
set(i) ;
106 fich.ignore(1000,
'\n') ;
108 fich >> n_lim ; fich.ignore(1000,
'\n') ;
110 fich >> kappa_low ; fich.ignore(1000,
'\n') ;
112 fich >> gamma_low ; fich.ignore(1000,
'\n') ;
114 fich >> m_n ; fich.ignore(1000,
'\n') ;
116 fich >> mu_0 ; fich.ignore(1000,
'\n') ;
118 fich >> Lambda ; fich.ignore(1000,
'\n') ;
120 double x_lim =
log(n_lim*0.1) ;
121 double alpha = (*coefs)(n_coefs-1) ;
122 double sum_poly = 0., sum_der_poly = 0. ;
123 for (
int i=0; i<n_coefs-1; i++){
124 sum_poly += (*coefs)(i)*
pow(x_lim,
double(i)) ;
125 sum_der_poly += (i == 0) ? 0. : (*coefs)(i)*
double(i)*
pow(x_lim,
double(i-1)) ;
127 ent_lim = log1p(mu_0 +
exp(alpha*x_lim)*((alpha+1)*sum_poly + sum_der_poly)) ;
130 eos_low =
new Eos_poly(gamma_low, kappa_low) ;
141 if (coefs != 0x0)
delete coefs ;
142 if (eos_low != 0x0)
delete eos_low ;
155 n_coefs = eosi.n_coefs ;
156 gamma_low = eosi.gamma_low ;
157 kappa_low = eosi.kappa_low ;
159 ent_lim = eosi.ent_lim ;
161 Lambda = eosi.Lambda ;
162 eos_low = eosi.eos_low ;
174 cout <<
"The second EOS is not of type Pseudo_polytrope_1D !" << endl ;
181 for (
int i=0; i < n_coefs; i++)
182 if ((*eos.coefs)(i) != (*coefs)(i)) {
184 <<
"The two Pseudo_polytrope_1D have different coefficients: " << (*coefs)(i) <<
" <-> " 185 << (*eos.coefs)(i) << endl ;
189 if (eos.n_lim != n_lim) {
191 <<
"The two Pseudo_polytrope_1D have different limiting densities n_lim: " << n_lim <<
" <-> " 192 << eos.n_lim << endl ;
196 if (eos.ent_lim != ent_lim) {
198 <<
"The two Pseudo_polytrope_1D have different limiting enthalpies ent_lim: " << ent_lim <<
" <-> " 199 << eos.ent_lim << endl ;
203 if (eos.n_coefs != n_coefs) {
205 <<
"The two Pseudo_polytrope_1D have different number of coefficients: " << n_coefs <<
" <-> " 206 << eos.n_coefs << endl ;
210 if (eos.kappa_low != kappa_low) {
212 <<
"The two Pseudo_polytrope_1D have different kappa in low polytrope: " << kappa_low <<
" <-> " 213 << eos.kappa_low << endl ;
217 if (eos.gamma_low != gamma_low) {
219 <<
"The two Pseudo_polytrope_1D have different kappa in low polytrope: " << gamma_low <<
" <-> " 220 << eos.gamma_low << endl ;
224 if (eos.mu_0!= mu_0) {
226 <<
"The two Pseudo_polytrope_1D have different mu_0 in low polytrope: " << mu_0 <<
" <-> " 227 << eos.mu_0 << endl ;
252 fwrite_be(&n_coefs,
sizeof(
int), 1, fich) ;
254 fwrite_be(&kappa_low,
sizeof(
double), 1, fich) ;
255 fwrite_be(&gamma_low,
sizeof(
double), 1, fich) ;
256 fwrite_be(&n_lim,
sizeof(
double), 1, fich) ;
257 fwrite_be(&ent_lim,
sizeof(
double), 1, fich) ;
258 fwrite_be(&m_n,
sizeof(
double), 1, fich) ;
259 fwrite_be(&mu_0,
sizeof(
double), 1, fich) ;
260 fwrite_be(&Lambda,
sizeof(
double), 1, fich) ;
267 ost << setprecision(16) <<
"EOS of class Pseudo_polytrope_1D (analytical fit of cold EoS): " <<
'\n' ;
268 ost <<
" Coefficients: " << *coefs <<
'\n' ;
269 ost << setprecision(16) <<
" Baryon mass: " << m_n <<
" [MeV]" <<
'\n' ;
270 ost << setprecision(16) <<
" mu_0-1: " << mu_0 <<
" and Lambda: " << Lambda <<
'\n' ;
271 ost <<
"Low densities polytrope :" <<
'\n' ;
272 cout << *eos_low <<
'\n' ;
273 ost << setprecision(16) <<
" Limiting density n_lim: " << 0.1*n_lim <<
" [fm^-3]" <<
'\n' ;
274 ost << setprecision(16) <<
" Limiting enthalpy h_lim: " << ent_lim <<
" [c^2]" << endl;
290 if (ent > ent_lim ) {
292 auto ent_nb_p = [&](
double nbar) {
293 double xx =
log(nbar*0.1) ;
294 double alpha = (*coefs)(n_coefs-1) ;
295 double sum_poly = 0.;
296 for (
int i=0; i < n_coefs-1; i++){
297 double cc1 = (alpha+1.)*(*coefs)(i) ;
298 double cc2 = (i == n_coefs - 2) ? 0. :
double(i+1)*(*coefs)(i+1) ;
299 sum_poly += (cc1 + cc2)*
pow(xx,
double(i));
301 double arg =
exp(alpha*xx)*sum_poly + mu_0 ;
302 return log1p(arg) - ent ;
305 double a = n_lim/2., b =
max(100.*ent, 50.*n_lim) ;
306 double f0 = ent_nb_p(a), f1 = ent_nb_p(b) ;
307 double c=1., c_old=2., f2 ;
309 assert( f0 * f1 < 0.) ;
310 while(fabs((c-c_old)/c)>eps) {
312 c = (b*f0 - a*f1)/(f0 - f1) ;
326 else if ( (0 < ent) && (ent < ent_lim)){
339 if (ent > ent_lim ) {
342 double xx =
log(nbar*0.1) ;
344 double alpha = (*coefs)(n_coefs-1) ;
345 double sum_poly = 0.;
346 for (
int i=0; i < n_coefs-1; i++)
347 sum_poly += (*coefs)(i)*
pow(xx,
double(i));
348 return m_n*Unites::mevpfm3*(
exp(xx)*(1.+mu_0) +
exp((alpha+1) * xx)*sum_poly) - Lambda ;
350 else if ( (0 < ent) && (ent < ent_lim)){
363 if (ent > ent_lim ) {
367 double xx =
log(nbar*0.1) ;
369 double alpha = (*coefs)(n_coefs-1) ;
370 double sum_poly = 0.;
371 for (
int i=0; i < n_coefs-1; i++){
372 double cc1 = alpha*(*coefs)(i) ;
373 double cc2 = (i == n_coefs - 2) ? 0. :
double(i+1)*(*coefs)(i+1) ;
374 sum_poly += (cc1 + cc2)*
pow(xx,
double(i));
377 return m_n*Unites::mevpfm3*(
exp((alpha+1.)*xx) * sum_poly) + Lambda;
379 else if ( (0 < ent) && (ent < ent_lim)){
412 if (ent > ent_lim ) {
414 double xx =
log(nbar*0.1) ;
416 double alpha = (*coefs)(n_coefs-1) ;
417 double sum_poly = 0., sum_der_poly = 0., sum_der2_poly = 0.;
418 for (
int i=0; i < n_coefs-1; i++){
419 sum_poly += (*coefs)(i)*
pow(xx,
double(i));
420 sum_der_poly += (i == 0) ? 0. :
double(i) * (*coefs)(i) *
pow(xx,
double(i-1)) ;
421 sum_der2_poly += (i < 2) ? 0. :
double(i) * double(i-1) * (*coefs)(i) *
pow(xx,
double(i-2)) ;
423 double num = (alpha*(alpha+1.) * sum_poly + (2.*alpha+1.) * sum_der_poly + sum_der2_poly)*
exp(alpha*xx) ;
424 double denom = 1. + mu_0 + ((alpha+1.) * sum_poly + sum_der_poly)*
exp(alpha*xx) ;
428 else if ( (0 < ent) && (ent < ent_lim)){
virtual double ener_ent_p(double ent, const Param *par=0x0) const
Computes the total energy density from the log-enthalpy.
Cmp log(const Cmp &)
Neperian logarithm.
Cmp exp(const Cmp &)
Exponential.
virtual int identify() const
Returns a number to identify the sub-classe of Eos the object belongs to.
Equation of state base class.
double & set(int i)
Read/write of a particular element (index i) (1D case)
virtual int identify() const =0
Returns a number to identify the sub-classe of Eos the object belongs to.
virtual double nbar_ent_p(double ent, const Param *par=0x0) const
Computes the baryon density from the specific enthalpy.
virtual double csound_square_ent_p(double, const Param *) const
Computes the sound speed squared from the enthapy with extra parameters (virtual function implemente...
virtual void sauve(FILE *) const
Save in a file.
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
int get_ndim() const
Gives the number of dimensions (ie dim.ndim)
virtual bool operator==(const Eos &) const
Comparison operator (egality)
virtual double ener_ent_p(double ent, const Param *par=0x0) const
Computes the total energy density from the specific enthalpy.
Pseudo_polytrope_1D(const Tbl &, double, double)
Standard constructor.
virtual double nbar_ent_p(double ent, const Param *par=0x0) const
Computes the baryon density from the log-enthalpy.
virtual double press_ent_p(double ent, const Param *par=0x0) const
Computes the pressure from the specific enthalpy.
Polytropic equation of state (relativistic case).
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
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...
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.
virtual double press_ent_p(double ent, const Param *par=0x0) const
Computes the pressure from the log-enthalpy.
Cmp pow(const Cmp &, int)
Power .
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
void operator=(const Pseudo_polytrope_1D &)
Assignment to another Pseudo_polytrope_1D.
virtual ostream & operator>>(ostream &) const
Operator >>
virtual double der_press_ent_p(double ent, const Param *par=0x0) const
Computes the logarithmic derivative from the specific enthalpy.
void c_est_pas_fait(const char *)
Helpful function to say something is not implemented yet.
void set_name(const char *name_i)
Sets the EOS name.
void sauve(FILE *) const
Save in a file.
virtual ~Pseudo_polytrope_1D()
Destructor.
virtual double der_nbar_ent_p(double ent, const Param *par=0x0) const
Computes the logarithmic derivative from the specific enthalpy.
virtual double der_ener_ent_p(double ent, const Param *par=0x0) const
Computes the logarithmic derivative from the specific enthalpy.
virtual bool operator!=(const Eos &) const
Comparison operator (difference)
virtual void sauve(FILE *) const
Save in a file.