36 m_v.resize(m_TPS->Particles().size());
38 m_TAG =
"ThermalModelEVDiagonal";
51 if (
m_v.size() != m_TPS->Particles().size())
52 m_v.resize(m_TPS->Particles().size());
53 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
60 if (ri.size() != m_TPS->Particles().size()) {
61 std::cerr <<
"**WARNING** " << m_TAG <<
"::FillVirial(const std::vector<double> & ri): size of ri does not match number of hadrons in the list" << std::endl;
64 m_v.resize(m_TPS->Particles().size());
65 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i)
71 if (vi.size() != m_TPS->Particles().size()) {
72 std::cerr <<
"**WARNING** " << m_TAG <<
"::FillVirialEV(const std::vector<double> & vi): size of vi does not match number of hadrons in the list" << std::endl;
80 m_v = std::vector<double>(m_TPS->Particles().size(), 0.);
82 ifstream fin(filename.c_str());
85 fin.getline(cc, 2000);
86 string tmp = string(cc);
90 istringstream iss(elems[0]);
93 if (iss >> pdgid >> b) {
94 int ind = m_TPS->PdgToId(pdgid);
104 ofstream fout(filename.c_str());
105 fout <<
"# List of eigenvolume parameters to be used in the Diagonal excluded-volume HRG model"
107 fout <<
"# Only particles with a non-zero eigenvolume parameter are listed here"
109 fout <<
"#" << std::setw(14) <<
"pdg"
110 << std::setw(15) <<
"v_i[fm^3]"
112 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
114 fout << std::setw(15) << m_TPS->Particle(i).PdgId();
115 fout << std::setw(15) <<
m_v[i];
124 if (i < 0 || i >=
static_cast<int>(
m_v.size()))
151 return m_TPS->Particles()[i].ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i] + dMu);
157 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
158 double dMu = -
m_v[i] * P;
166 BroydenEquationsDEV eqs(
this);
167 BroydenJacobianDEV jac(
this);
169 BroydenSolutionCriteriumDEV crit(
this, 1.0E-8);
173 if (m_Parameters.T == 0.0)
178 std::vector<double> x(1, x0);
180 x = broydn.
Solve(x, &crit);
184 double PressureNew = mnc * exp(x[0]);
198 m_LastCalculationSuccessFlag =
false;
199 else m_LastCalculationSuccessFlag =
true;
205 m_FluctuationsCalculated =
false;
213 double densityid = 0., suppression = 0.;
217#pragma omp parallel for reduction(+:densityid) reduction(+:suppression) if(useOpenMP)
218 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
232 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
235 m_wnSum += m_densities[i] * m_TPS->Particles()[i].ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i] -
m_v[i] *
m_Pressure);
245 int NN = m_densities.size();
246 vector<double> tN(NN);
249 vector<double> chi2id(NN);
250 for (
int i = 0; i < NN; ++i) {
253 chi2id[i] *= m_densities[i] / tN[i];
256 m_PrimCorrel.resize(NN);
257 for (
int i = 0; i < NN; ++i) m_PrimCorrel[i].resize(NN);
258 m_dmusdmu = m_PrimCorrel;
259 m_TotalCorrel = m_PrimCorrel;
261 for (
int i = 0; i < NN; ++i) {
262 for (
int j = 0; j < NN; ++j) {
263 m_PrimCorrel[i][j] = 0.;
264 if (i == j) m_PrimCorrel[i][j] += chi2id[i];
265 m_PrimCorrel[i][j] += -
m_v[i] * m_densities[j] * chi2id[i];
266 m_PrimCorrel[i][j] += -
m_v[j] * m_densities[i] * chi2id[j];
268 for (
size_t k = 0; k < m_densities.size(); ++k)
269 tmp +=
m_v[k] *
m_v[k] * chi2id[k];
270 m_PrimCorrel[i][j] += m_densities[i] * m_densities[j] * tmp;
272 m_dmusdmu[i][j] = (i == j ? 1. : 0.) -
m_v[i] * m_densities[j];
276 for (
int i = 0; i < NN; ++i) {
277 m_wprim[i] = m_PrimCorrel[i][i];
278 if (m_densities[i] > 0.) m_wprim[i] *= m_Parameters.T / m_densities[i];
279 else m_wprim[i] = 1.;
285 if (!IsTemperatureDerivativesCalculated())
291 double eid = 0., dTbn = 0., bn = 0., dTeidterm = 0., dmueidterm = 0.;
292 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
295 dTbn +=
m_v[i] * m_dndT[i];
296 bn +=
m_v[i] * m_densities[i];
301 ret = -dTbn * eid + (1. - bn) * (dTeidterm + dmueidterm);
307 int N = m_TPS->ComponentsNumber();
308 m_dndT = vector<double>(N, 0.);
309 m_dmusdT = vector<double>(N, 0.);
310 m_PrimChi2sdT = vector<vector<double>>(N, vector<double>(N, 0.));
312 double T = m_Parameters.T;
314 vector<double> nids(N, 0.), dniddTs(N, 0.), chi2ids(N, 0.), dchi2idsdT(N, 0.), chi3ids(N, 0.);
315 for (
int i = 0; i < N; ++i) {
323 double s = EntropyDensity();
324 double sum_chi2idb2 = 0., sum_bns= 0., sum_bn = 0., sum_dTbns = 0., sum_dchi2dTidb2 = 0., sum_chi3idb3 = 0.;
325 for (
int i = 0; i < N; ++i) {
326 m_dmusdT[i] = -
m_v[i] * s;
327 sum_chi2idb2 +=
m_v[i] *
m_v[i] * chi2ids[i];
328 sum_bns +=
m_v[i] * nids[i];
329 sum_bn +=
m_v[i] * m_densities[i];
330 sum_dTbns +=
m_v[i] * dniddTs[i];
331 sum_dchi2dTidb2 +=
m_v[i] *
m_v[i] * dchi2idsdT[i];
332 sum_chi3idb3 +=
m_v[i] *
m_v[i] *
m_v[i] * chi3ids[i];
335 double sum_dbndT = (1. - sum_bn) * (1. - sum_bn) * (sum_dTbns - sum_chi2idb2 * s);
337 for (
int i = 0; i < N; ++i) {
338 m_dndT[i] = -sum_dbndT * nids[i] + (1. - sum_bn) * (dniddTs[i] + chi2ids[i] * m_dmusdT[i]);
341 if (IsSusceptibilitiesCalculated()) {
342 for (
int i = 0; i < N; ++i) {
343 for (
int j = 0; j < N; ++j) {
344 m_PrimChi2sdT[i][j] = 0.;
347 tval += -chi2ids[j] * m_densities[i] *
m_v[j];
348 tval += sum_chi2idb2 * m_densities[j] * m_densities[i];
351 tval += -chi2ids[i] * m_densities[j] *
m_v[i];
352 m_PrimChi2sdT[i][j] += -sum_dbndT * tval;
355 tval += -dchi2idsdT[j] * m_densities[i] *
m_v[j];
356 tval += -chi3ids[j] * m_dmusdT[j] * m_densities[i] *
m_v[j];
357 tval += -chi2ids[j] * m_dndT[i] *
m_v[j];
358 tval += sum_dchi2dTidb2 * m_densities[j] * m_densities[i];
359 tval += -sum_chi3idb3 * s * m_densities[j] * m_densities[i];
360 tval += sum_chi2idb2 * (m_dndT[j] * m_densities[i] + m_densities[j] * m_dndT[i]);
362 tval += dchi2idsdT[i];
363 tval += chi3ids[i] * m_dmusdT[i];
365 tval += -dchi2idsdT[i] * m_densities[j] *
m_v[i];
366 tval += -chi3ids[i] * m_dmusdT[i] * m_densities[j] *
m_v[i];
367 tval += -chi2ids[i] * m_dndT[j] *
m_v[i];
369 m_PrimChi2sdT[i][j] += (1. - sum_bn) * tval;
375 m_TemperatureDerivativesCalculated =
true;
379 m_TemperatureDerivativesCalculated =
true;
390 m_FluctuationsCalculated =
true;
392 for (
size_t i = 0; i < m_wprim.size(); ++i) {
397 for (
size_t i = 0; i < m_wprim.size(); ++i) {
402 if (IsTemperatureDerivativesCalculated()) {
403 m_TemperatureDerivativesCalculated =
false;
411 vector<double> ret(order + 1, 0.);
414 for (
size_t i = 0; i < m_densities.size(); ++i)
415 ret[0] += chgs[i] * m_densities[i];
419 if (order < 2)
return ret;
424 int NN = m_densities.size();
426 vector<double> MuStar(NN, 0.);
427 for (
int i = 0; i < NN; ++i) {
428 MuStar[i] = m_Chem[i] +
MuShift(i);
432 MatrixXd densMatrix(2 * NN, 2 * NN);
433 VectorXd solVector(2 * NN), xVector(2 * NN);
435 vector<double> DensitiesId(m_densities.size()), chi2id(m_densities.size());
436 for (
int i = 0; i < NN; ++i) {
438 chi2id[i] = m_TPS->Particles()[i].chi(2, m_Parameters, m_UseWidth, MuStar[i]);
441 for (
int i = 0; i < NN; ++i)
442 for (
int j = 0; j < NN; ++j) {
443 densMatrix(i, j) =
m_v[j] * DensitiesId[i];
444 if (i == j) densMatrix(i, j) += 1.;
447 for (
int i = 0; i < NN; ++i)
448 for (
int j = 0; j < NN; ++j)
449 densMatrix(i, NN + j) = 0.;
451 for (
int i = 0; i < NN; ++i) {
452 densMatrix(i, NN + i) = 0.;
453 for (
int k = 0; k < NN; ++k) {
454 densMatrix(i, NN + i) +=
m_v[k] * m_densities[k];
456 densMatrix(i, NN + i) = (densMatrix(i, NN + i) - 1.) * chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T;
459 for (
int i = 0; i < NN; ++i)
460 for (
int j = 0; j < NN; ++j) {
461 densMatrix(NN + i, j) = 0.;
462 densMatrix(NN + i, NN + j) =
m_v[i] * DensitiesId[j];
463 if (i == j) densMatrix(NN + i, NN + j) += 1.;
467 PartialPivLU<MatrixXd> decomp(densMatrix);
470 vector<double> dni(NN, 0.), dmus(NN, 0.);
472 for (
int i = 0; i < NN; ++i) {
474 xVector[NN + i] = chgs[i];
477 solVector = decomp.solve(xVector);
479 for (
int i = 0; i < NN; ++i) {
480 dni[i] = solVector[i];
481 dmus[i] = solVector[NN + i];
484 for (
int i = 0; i < NN; ++i)
485 ret[1] += chgs[i] * dni[i];
489 if (order < 3)
return ret;
491 vector<double> d2ni(NN, 0.), d2mus(NN, 0.);
493 vector<double> chi3id(m_densities.size());
494 for (
int i = 0; i < NN; ++i)
495 chi3id[i] = m_TPS->Particles()[i].chi(3, m_Parameters, m_UseWidth, MuStar[i]);
497 for (
int i = 0; i < NN; ++i) {
501 for (
int j = 0; j < NN; ++j) tmp +=
m_v[j] * dni[j];
502 tmp = -2. * tmp * chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[i];
506 for (
int j = 0; j < NN; ++j) tmp +=
m_v[j] * m_densities[j];
507 tmp = -(tmp - 1.) * chi3id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * dmus[i] * dmus[i];
510 for (
int i = 0; i < NN; ++i) {
511 xVector[NN + i] = 0.;
514 for (
int j = 0; j < NN; ++j) tmp += -
m_v[i] * dmus[j] * chi2id[j] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[j];
516 xVector[NN + i] = tmp;
519 solVector = decomp.solve(xVector);
521 for (
int i = 0; i < NN; ++i) {
522 d2ni[i] = solVector[i];
523 d2mus[i] = solVector[NN + i];
526 for (
int i = 0; i < NN; ++i)
527 ret[2] += chgs[i] * d2ni[i];
532 if (order < 4)
return ret;
535 vector<double> d3ni(NN, 0.), d3mus(NN, 0.);
537 vector<double> chi4id(m_densities.size());
538 for (
int i = 0; i < NN; ++i)
539 chi4id[i] = m_TPS->Particles()[i].chi(4, m_Parameters, m_UseWidth, MuStar[i]);
541 vector<double> dnis(NN, 0.);
542 for (
int i = 0; i < NN; ++i) {
543 dnis[i] = chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[i];
546 vector<double> d2nis(NN, 0.);
547 for (
int i = 0; i < NN; ++i) {
548 d2nis[i] = chi3id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * dmus[i] * dmus[i] +
549 chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * d2mus[i];
552 for (
int i = 0; i < NN; ++i) {
556 for (
int j = 0; j < NN; ++j) tmp +=
m_v[j] * dni[j];
557 tmp = -3. * tmp * d2nis[i];
561 for (
int j = 0; j < NN; ++j) tmp +=
m_v[j] * d2ni[j];
562 tmp = -3. * tmp * dnis[i];
566 for (
int j = 0; j < NN; ++j) tmps +=
m_v[j] * m_densities[j];
568 tmp = -(tmps - 1.) * chi3id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * d2mus[i] * 3. * dmus[i];
571 tmp = -(tmps - 1.) * chi4id[i] * pow(
xMath::GeVtoifm(), 3) * dmus[i] * dmus[i] * dmus[i];
574 for (
int i = 0; i < NN; ++i) {
575 xVector[NN + i] = 0.;
578 for (
int j = 0; j < NN; ++j) tmp += -2. *
m_v[i] * d2mus[j] * dnis[j];
579 xVector[NN + i] += tmp;
582 for (
int j = 0; j < NN; ++j) tmp += -
m_v[i] * dmus[j] * d2nis[j];
583 xVector[NN + i] += tmp;
586 solVector = decomp.solve(xVector);
588 for (
int i = 0; i < NN; ++i) {
589 d3ni[i] = solVector[i];
590 d3mus[i] = solVector[NN + i];
593 for (
int i = 0; i < NN; ++i)
594 ret[3] += chgs[i] * d3ni[i];
606 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
617 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
646 if (
id >= 0. &&
id <
static_cast<int>(
m_v.size()))
654 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
655 tEV +=
m_v[i] * m_densities[i] / 4.;
662 if (i >= 0 && i <
static_cast<int>(
m_v.size()))
677 std::vector<double> ThermalModelEVDiagonal::BroydenEquationsDEV::Equations(
const std::vector<double>& x)
679 std::vector<double> ret(1);
680 double pressure = m_mnc * exp(x[0]);
682 ret[0] = pressure - m_THM->
Pressure(pressure);
687 std::vector<double> ThermalModelEVDiagonal::BroydenJacobianDEV::Jacobian(
const std::vector<double>& x)
689 double pressure = m_mnc * exp(x[0]);
692 for (
size_t i = 0; i < m_THM->Densities().size(); ++i)
693 ret += m_THM->VirialCoefficient(i, i) * m_THM->DensityId(i, pressure);
697 return std::vector<double>(1, ret);
700 std::vector<double> ThermalModelEVDiagonal::BroydenEquationsDEVOrig::Equations(
const std::vector<double>& x)
702 std::vector<double> ret(1);
703 ret[0] = x[0] - m_THM->Pressure(x[0]);
707 std::vector<double> ThermalModelEVDiagonal::BroydenJacobianDEVOrig::Jacobian(
const std::vector<double>& x)
709 const double &pressure = x[0];
712 for (
size_t i = 0; i < m_THM->Densities().size(); ++i)
713 ret += m_THM->VirialCoefficient(i, i) * m_THM->DensityId(i, pressure);
716 return std::vector<double>(1, ret);
719 bool ThermalModelEVDiagonal::BroydenSolutionCriteriumDEV::IsSolved(
const std::vector<double>& x,
const std::vector<double>& f,
const std::vector<double>& )
const
722 for (
size_t i = 0; i < x.size(); ++i) {
723 maxdiff = std::max(maxdiff, fabs(f[i]) / m_THM->m_Pressure);
725 return (maxdiff < m_MaximumError);
Contains some functions to deal with excluded volumes.
map< string, double > params
Class implementing the Broyden method to solve a system of non-linear equations.
double MaxDifference() const
virtual std::vector< double > Solve(const std::vector< double > &x0, BroydenSolutionCriterium *solcrit=NULL, int max_iterations=MAX_ITERS)
int MaxIterations() const
virtual void CalculateSusceptibilityMatrix()
Calculates the conserved charges susceptibility matrix.
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
virtual void CalculateParticleChargeCorrelationMatrix()
Calculates the matrix of correlators between primordial (and also final) particle numbers and conserv...
@ DiagonalEV
Diagonal excluded volume model.
virtual void CalculateFeeddown()
Calculates the total densities which include feeddown contributions.
virtual void ValidateCalculation()
Checks whether issues have occured during the calculation of particle densities in the CalculateDensi...
virtual void CalculateProxySusceptibilityMatrix()
Calculates the susceptibility matrix of conserved charges proxies.
ThermalModelBase(ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters())
Construct a new ThermalModelBase object.
virtual void CalculateDensities()
Calculates the primordial and total (after decays) densities of all species.
@ GCE
Grand canonical ensemble.
virtual void CalculateTwoParticleFluctuationsDecays()
Computes particle number correlations and fluctuations for all final-state particles which are marked...
double DensityId(int i, double Pressure)
Calculate the ideal gas density of particle species i for the given pressure value.
double ExcludedVolume(int i) const
The excluded volume parameter of particle species i.
std::vector< double > m_densitiesidnoshift
void FillVirial(const std::vector< double > &ri=std::vector< double >(0))
Fills the vector of particle eigenvolume parameters.
ThermalModelEVDiagonal(ThermalParticleSystem *TPS, const ThermalModelParameters ¶ms=ThermalModelParameters())
Construct a new ThermalModelEVDiagonal object.
virtual double ParticleScalarDensity(int part)
virtual void CalculateFluctuations()
Computes the fluctuation observables.
virtual double ParticleKurtosis(int)
void SetRadius(double rad)
Set the same excluded volume radius parameter for all species.
virtual void CalculateTemperatureDerivatives()
Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron...
virtual double CalculateEnergyDensity()
void FillVirialEV(const std::vector< double > &vi=std::vector< double >(0))
Same as FillVirial() but uses the diagonal excluded-volume coefficients as input instead of radii.
virtual void WriteInteractionParameters(const std::string &filename)
Write the set of eigenvolume parameters for all particles to an external file.
void SetVirial(int i, int j, double b)
Sets the eigenvolume parameter of particle species i.
double CommonSuppressionFactor()
The density suppression factor , common for all species.
virtual std::vector< double > CalculateChargeFluctuations(const std::vector< double > &chgs, int order=4)
Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.
double PressureId(int i, double Pressure)
Calculate the ideal gas pressure of particle species i for the given pressure value.
std::vector< double > m_densitiesid
virtual double CalculatePressure()
virtual double MuShift(int i) const
The shift in the chemical potential of particle species i due to the excluded volume interactions.
virtual double CalculateEntropyDensity()
virtual void SolvePressure()
Solves the transcdental equation for the pressure.
std::vector< double > m_v
virtual void CalculatePrimordialDensities()
Calculates the primordial densities of all species.
virtual double ParticleSkewness(int)
virtual double CalculateEigenvolumeFraction()
double ScaledVarianceId(int i, double Pressure)
Calculate the ideal gas scaled variance of particle species i number fluctuations for the given press...
virtual void ReadInteractionParameters(const std::string &filename)
Read the set of eigenvolume parameters for all particles from an external file.
virtual double CalculateEnergyDensityDerivativeT()
virtual void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
double Pressure(double P)
Computes the l.h.s. of the transcdental equation for the pressure.
double VirialCoefficient(int i, int j) const
Return the eigenvolume parameter of particle species i.
virtual ~ThermalModelEVDiagonal(void)
Destroy the ThermalModelEVDiagonal object.
Class containing all information about a particle specie.
double Density(const ThermalModelParameters ¶ms, IdealGasFunctions::Quantity type=IdealGasFunctions::ParticleDensity, bool useWidth=0, double mu=0.) const
Computes a specified ideal gas thermodynamic function.
Class containing the particle list.
double vr(double r)
Computes the excluded volume parameter from a given radius parameter value.
std::vector< std::string > split(const std::string &s, char delim)
constexpr double GeVtoifm()
A constant to transform GeV into fm .
constexpr double GeVtoifm3()
A constant to transform GeV into fm .
constexpr double mnucleon()
Nucleon's mass. Value as in UrQMD.
The main namespace where all classes and functions of the Thermal-FIST library reside.
Structure containing all thermal parameters of the model.
Contains some extra mathematical functions used in the code.