38 for(
int i=0;i<6;++i)
m_chi[i].resize(3);
42 m_Virial.resize(m_densities.size(), vector<double>(m_densities.size(),0.));
47 m_TAG =
"ThermalModelVDW";
50 m_InteractionModel =
QvdW;
60 for(
size_t i = 0; i <
m_MuStar.size(); ++i)
67 for (
size_t i = 0; i <
m_MuStar.size(); ++i)
72 if (ri.size() != m_TPS->Particles().size()) {
73 std::cerr <<
"**WARNING** " << m_TAG <<
"::FillVirial(const vector<double> & ri): size of ri does not match number of hadrons in the list" << std::endl;
76 m_Virial.resize(m_TPS->Particles().size());
77 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
78 m_Virial[i].resize(m_TPS->Particles().size());
79 for (
int j = 0; j < m_TPS->ComponentsNumber(); ++j)
84 vector< vector<double> > fVirialTmp =
m_Virial;
85 for(
int i=0;i<m_TPS->ComponentsNumber();++i)
86 for(
int j=0;j<m_TPS->ComponentsNumber();++j) {
87 if (i==j)
m_Virial[i][j] = fVirialTmp[i][j];
88 else if ((fVirialTmp[i][i] + fVirialTmp[j][j]) > 0.0)
m_Virial[i][j] = 2. * fVirialTmp[i][j] * fVirialTmp[i][i] / (fVirialTmp[i][i] + fVirialTmp[j][j]);
91 m_VirialdT = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(),0.));
92 m_AttrdT = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(), 0.));
99 if (bij.size() != m_TPS->Particles().size()) {
100 std::cerr <<
"**WARNING** " << m_TAG <<
"::FillVirialEV(const vector<double> & bij): size of bij does not match number of hadrons in the list" << std::endl;
110 if (aij.size() != m_TPS->Particles().size()) {
111 std::cerr <<
"**WARNING** " << m_TAG <<
"::FillAttraction(const vector<double> & aij): size of aij does not match number of hadrons in the list" << std::endl;
121 m_Virial = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(), 0.));
122 m_Attr = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(), 0.));
124 ifstream fin(filename.c_str());
127 fin.getline(cc, 2000);
128 string tmp = string(cc);
130 if (elems.size() < 1)
132 istringstream iss(elems[0]);
133 long long pdgid1, pdgid2;
135 if (iss >> pdgid1 >> pdgid2 >> b) {
138 int ind1 = m_TPS->PdgToId(pdgid1);
139 int ind2 = m_TPS->PdgToId(pdgid2);
140 if (ind1 != -1 && ind2 != -1) {
153 ofstream fout(filename.c_str());
154 fout <<
"# List of the van dar Waals interaction parameters to be used in the QvdW-HRG model"
156 fout <<
"# Only particle pairs with a non-zero QvdW interaction are listed here"
163 fout <<
"#" <<
" " <<
"pdg_i"
165 <<
" " <<
"b_{ij}[fm^3]"
166 <<
" " <<
"a_{ij}[GeV*fm^3]"
168 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
169 for (
int j = 0; j < m_TPS->ComponentsNumber(); ++j) {
176 fout <<
" " << m_TPS->Particle(i).PdgId();
177 fout <<
" " << m_TPS->Particle(j).PdgId();
179 fout <<
" " <<
m_Attr[i][j];
189 m_Virial = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(), 0.));
190 m_Attr = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(), 0.));
191 m_VirialdT = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(), 0.));
192 m_AttrdT = vector< vector<double> >(m_TPS->Particles().size(), vector<double>(m_TPS->Particles().size(), 0.));
198 vector<double> ns(m_densities.size(), 0.);
199 for (
size_t i = 0; i < m_densities.size(); ++i)
207 int NN = m_densities.size();
210 MatrixXd densMatrix(NNdmu, NNdmu);
211 VectorXd solVector(NNdmu), xVector(NNdmu);
213 for (
int i = 0; i < NNdmu; ++i) {
214 for (
int j = 0; j < NNdmu; ++j) {
215 densMatrix(i, j) = 0.;
217 densMatrix(i, j) += 1.;
225 PartialPivLU<MatrixXd> decomp(densMatrix);
227 for (
int kp = 0; kp < NNdmu; ++kp) {
235 solVector = decomp.solve(xVector);
237 vector<double> ntildek(NNdmu, 0.);
238 for (
int i = 0; i < NNdmu; ++i)
239 ntildek[i] = solVector[i];
241 vector<double> np(m_densities.size(), 0.);
242 for (
int i = 0; i < NN; ++i) {
244 for (
int k = 0; k < NNdmu; ++k) {
247 np[i] = (1. - np[i]) * ns[i];
255 map< vector<double>,
int> MapVDW;
257 int NN = m_densities.size();
265 for (
int i = 0; i < NN; ++i) {
266 vector<double> VDWParam(0);
267 for (
int j = 0; j < NN; ++j) {
270 for (
int j = 0; j < NN; ++j) {
274 if (MapVDW.count(VDWParam) == 0) {
275 MapVDW[VDWParam] = tind;
292 int NN = m_densities.size();
295 vector<double> dmuscur(NNdmu, 0.);
296 for (
int i = 0; i < NNdmu; ++i)
299 BroydenEquationsVDW eqs(
this);
300 BroydenJacobianVDW jac(
this);
302 BroydenSolutionCriteriumVDW crit(
this);
304 dmuscur = broydn.
Solve(dmuscur, &crit);
312 vector<double> ret(NN);
313 for (
int i = 0; i < NN; ++i)
320 vector<double> csol(m_densities.size(), 0.);
325 double dmu = (muBmax - muBmin) / iters;
326 vector<double> curmust(m_densities.size(), 0.);
328 for(
int isol = 0; isol < iters; ++isol) {
329 double tmu = muBmin + (0.5 + isol) * dmu;
330 for(
size_t j = 0; j < curmust.size(); ++j) {
331 curmust[j] = m_Chem[j] + (tmu - m_Parameters.muB) * m_Chem[j] / m_Parameters.muB;
332 if (m_TPS->Particles()[j].Statistics()==-1 && curmust[j] > m_TPS->Particles()[j].Mass())
333 curmust[j] = 0.98 * m_TPS->Particles()[j].Mass();
339 for(
size_t i = 0; i < sol.size(); ++i)
340 if (sol[i] != sol[i]) fl =
false;
344 for(
int i = 0; i < m_TPS->ComponentsNumber(); ++i)
347 int NN = m_densities.size();
349 MatrixXd densMatrix(NN, NN);
350 VectorXd solVector(NN), xVector(NN);
352 for(
int i=0;i<NN;++i)
353 for(
int j=0;j<NN;++j) {
355 if (i==j) densMatrix(i,j) += 1.;
358 PartialPivLU<MatrixXd> decomp(densMatrix);
360 for(
int i=0;i<NN;++i)
362 solVector = decomp.solve(xVector);
363 for(
int i=0;i<NN;++i)
364 m_densities[i] = solVector[i];
367 for(
size_t i=0;i<m_densities.size();++i)
369 for(
size_t i=0;i<m_densities.size();++i)
370 for(
size_t j=0;j<m_densities.size();++j)
371 tP += -
m_Attr[i][j] * m_densities[i] * m_densities[j];
373 if (!solved || tP > Psol) {
388 vector<double> muStarInit =
m_MuStar;
390 for(
size_t i=0;i<muStarInit.size();++i) {
391 if (m_TPS->Particles()[i].Statistics()==-1 && muStarInit[i] > m_TPS->Particles()[i].Mass())
392 muStarInit[i] = 0.98 * m_TPS->Particles()[i].Mass();
402 CalculatePrimordialDensitiesNew();
406 void ThermalModelVDW::CalculatePrimordialDensitiesOld() {
407 m_FluctuationsCalculated =
false;
412 int NN = m_densities.size();
447 clock_t tbeg = clock();
451 vector<double> muStarInit =
m_MuStar;
453 for (
size_t i = 0; i<muStarInit.size(); ++i) {
454 if (m_TPS->Particles()[i].Statistics() == -1 && muStarInit[i] > m_TPS->Particles()[i].Mass())
455 muStarInit[i] = 0.98 * m_TPS->Particles()[i].Mass();
463 printf(
"Solution time = %lf ms\n", (clock() - tbeg) / (
double)(CLOCKS_PER_SEC) * 1.e3);
467 for(
int i=0;i<m_TPS->ComponentsNumber();++i)
470 MatrixXd densMatrix(NN, NN);
471 VectorXd solVector(NN), xVector(NN);
473 for(
int i=0;i<NN;++i)
474 for(
int j=0;j<NN;++j) {
476 if (i==j) densMatrix(i,j) += 1.;
479 PartialPivLU<MatrixXd> decomp(densMatrix);
482 solVector = decomp.solve(xVector);
483 for(
int i=0;i<NN;++i) m_densities[i] = solVector[i];
487 solVector = decomp.solve(xVector);
489 for(
int i=0;i<NN;++i)
m_scaldens[i] = solVector[i];
497 void ThermalModelVDW::CalculatePrimordialDensitiesNew() {
498 m_FluctuationsCalculated =
false;
503 int NN = m_densities.size();
505 clock_t tbeg = clock();
511 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i)
517 MatrixXd densMatrix(NNdmu, NNdmu);
518 VectorXd solVector(NNdmu), xVector(NNdmu);
520 for (
int i = 0; i < NNdmu; ++i) {
521 for (
int j = 0; j < NNdmu; ++j) {
522 densMatrix(i, j) = 0.;
524 densMatrix(i, j) += 1.;
532 PartialPivLU<MatrixXd> decomp(densMatrix);
534 for (
int kp = 0; kp < NNdmu; ++kp) {
541 solVector = decomp.solve(xVector);
543 vector<double> ntildek(NNdmu, 0.);
544 for (
int i = 0; i < NNdmu; ++i)
545 ntildek[i] = solVector[i];
548 for (
int i = 0; i < NN; ++i) {
550 for (
int k = 0; k < NNdmu; ++k) {
563 vector<double> ret(order + 1, 0.);
566 for(
size_t i=0;i<m_densities.size();++i)
567 ret[0] += chgs[i] * m_densities[i];
571 if (order<2)
return ret;
573 int NN = m_densities.size();
574 MatrixXd densMatrix(2*NN, 2*NN);
575 VectorXd solVector(2*NN), xVector(2*NN);
577 vector<double> chi2id(m_densities.size());
578 for(
int i=0;i<NN;++i)
579 chi2id[i] = m_TPS->Particles()[i].chi(2, m_Parameters, m_UseWidth,
m_MuStar[i]);
581 for(
int i=0;i<NN;++i)
582 for(
int j=0;j<NN;++j) {
584 if (i==j) densMatrix(i,j) += 1.;
587 for(
int i=0;i<NN;++i)
588 for(
int j=0;j<NN;++j)
589 densMatrix(i,NN+j) = 0.;
591 for(
int i=0;i<NN;++i) {
592 densMatrix(i,NN+i) = 0.;
593 for(
int k=0;k<NN;++k) {
594 densMatrix(i,NN+i) +=
m_Virial[k][i] * m_densities[k];
596 densMatrix(i,NN+i) = (densMatrix(i,NN+i) - 1.) * chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T;
599 for(
int i=0;i<NN;++i)
600 for(
int j=0;j<NN;++j) {
604 for(
int i=0;i<NN;++i)
605 for(
int j=0;j<NN;++j) {
607 if (i==j) densMatrix(NN+i,NN+j) += 1.;
611 PartialPivLU<MatrixXd> decomp(densMatrix);
614 vector<double> dni(NN, 0.), dmus(NN, 0.);
616 for(
int i=0;i<NN;++i) {
618 xVector[NN+i] = chgs[i];
621 solVector = decomp.solve(xVector);
623 for(
int i=0;i<NN;++i) {
624 dni[i] = solVector[i];
625 dmus[i] = solVector[NN+i];
628 for(
int i=0;i<NN;++i)
629 ret[1] += chgs[i] * dni[i];
633 if (order<3)
return ret;
636 vector<double> d2ni(NN, 0.), d2mus(NN, 0.);
638 vector<double> chi3id(m_densities.size());
639 for(
int i=0;i<NN;++i)
640 chi3id[i] = m_TPS->Particles()[i].chi(3, m_Parameters, m_UseWidth,
m_MuStar[i]);
642 for(
int i=0;i<NN;++i) {
646 for(
int j=0;j<NN;++j) tmp +=
m_Virial[j][i] * dni[j];
647 tmp = -2. * tmp * chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[i];
651 for(
int j=0;j<NN;++j) tmp +=
m_Virial[j][i] * m_densities[j];
652 tmp = -(tmp - 1.) * chi3id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * dmus[i] * dmus[i];
655 for(
int i=0;i<NN;++i) {
659 for(
int j=0;j<NN;++j) tmp += -
m_Virial[i][j] * dmus[j] * chi2id[j] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[j];
664 solVector = decomp.solve(xVector);
666 for(
int i=0;i<NN;++i) {
667 d2ni[i] = solVector[i];
668 d2mus[i] = solVector[NN+i];
671 for(
int i=0;i<NN;++i)
672 ret[2] += chgs[i] * d2ni[i];
677 if (order<4)
return ret;
680 vector<double> d3ni(NN, 0.), d3mus(NN, 0.);
682 vector<double> chi4id(m_densities.size());
683 for (
int i = 0; i < NN; ++i)
684 chi4id[i] = m_TPS->Particles()[i].chi(4, m_Parameters, m_UseWidth,
m_MuStar[i]);
686 vector<double> dnis(NN, 0.);
687 for(
int i=0;i<NN;++i) {
688 dnis[i] = chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[i];
691 vector<double> d2nis(NN, 0.);
692 for(
int i=0;i<NN;++i) {
693 d2nis[i] = chi3id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * dmus[i] * dmus[i] +
694 chi2id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * d2mus[i];
697 for(
int i=0;i<NN;++i) {
701 for(
int j=0;j<NN;++j) tmp +=
m_Virial[j][i] * dni[j];
702 tmp = -3. * tmp * d2nis[i];
706 for(
int j=0;j<NN;++j) tmp +=
m_Virial[j][i] * d2ni[j];
707 tmp = -3. * tmp * dnis[i];
711 for(
int j=0;j<NN;++j) tmps +=
m_Virial[j][i] * m_densities[j];
713 tmp = -(tmps - 1.) * chi3id[i] * pow(
xMath::GeVtoifm(), 3) * m_Parameters.T * d2mus[i] * 3. * dmus[i];
716 tmp = -(tmps - 1.) * chi4id[i] * pow(
xMath::GeVtoifm(), 3) * dmus[i] * dmus[i] * dmus[i];
719 for(
int i=0;i<NN;++i) {
723 for(
int j=0;j<NN;++j) tmp += -2. *
m_Virial[i][j] * d2mus[j] * dnis[j];
724 xVector[NN+i] += tmp;
727 for(
int j=0;j<NN;++j) tmp += -
m_Virial[i][j] * dmus[j] * d2nis[j];
728 xVector[NN+i] += tmp;
731 solVector = decomp.solve(xVector);
733 for(
int i=0;i<NN;++i) {
734 d3ni[i] = solVector[i];
735 d3mus[i] = solVector[NN+i];
738 for(
int i=0;i<NN;++i)
739 ret[3] += chgs[i] * d3ni[i];
748 if (order<1)
return m_chi;
750 vector<double> chgs(m_densities.size());
754 for(
size_t i=0;i<chgs.size();++i)
755 chgs[i] = m_TPS->Particles()[i].BaryonCharge();
757 for(
int i=0;i<order;++i)
m_chi[i][0] = chis[i];
760 for(
size_t i=0;i<chgs.size();++i)
761 chgs[i] = m_TPS->Particles()[i].ElectricCharge();
763 for(
int i=0;i<order;++i)
m_chi[i][1] = chis[i];
766 for(
size_t i=0;i<chgs.size();++i)
767 chgs[i] = m_TPS->Particles()[i].Strangeness();
769 for(
int i=0;i<order;++i)
m_chi[i][2] = chis[i];
772 for(
size_t i=0;i<chgs.size();++i)
773 chgs[i] = m_TPS->Particles()[i].ArbitraryCharge();
775 for(
int i=0;i<order;++i)
m_chiarb[i] = chis[i];
782 int NN = m_densities.size();
784 m_PrimCorrel.resize(NN);
785 for (
int i = 0; i < NN; ++i)
786 m_PrimCorrel[i].resize(NN);
787 m_dmusdmu = m_PrimCorrel;
788 m_TotalCorrel = m_PrimCorrel;
790 MatrixXd densMatrix(2 * NN, 2 * NN);
791 VectorXd solVector(2 * NN), xVector(2 * NN);
793 vector<double> chi2id(m_densities.size());
794 for (
int i = 0; i<NN; ++i)
796 chi2id[i] = m_TPS->Particles()[i].chiDimensionfull(2, m_Parameters, m_UseWidth,
m_MuStar[i]);
798 for (
int i = 0; i<NN; ++i)
799 for (
int j = 0; j<NN; ++j) {
801 if (i == j) densMatrix(i, j) += 1.;
804 for (
int i = 0; i<NN; ++i)
805 for (
int j = 0; j<NN; ++j)
806 densMatrix(i, NN + j) = 0.;
808 for (
int i = 0; i<NN; ++i) {
809 densMatrix(i, NN + i) = 0.;
810 for (
int k = 0; k<NN; ++k) {
811 densMatrix(i, NN + i) +=
m_Virial[k][i] * m_densities[k];
814 densMatrix(i, NN + i) = (densMatrix(i, NN + i) - 1.) * chi2id[i] * pow(
xMath::GeVtoifm(), 3);
817 for (
int i = 0; i<NN; ++i)
818 for (
int j = 0; j<NN; ++j) {
822 for (
int i = 0; i<NN; ++i)
823 for (
int j = 0; j<NN; ++j) {
825 if (i == j) densMatrix(NN + i, NN + j) += 1.;
828 PartialPivLU<MatrixXd> decomp(densMatrix);
830 for (
int k = 0; k < NN; ++k) {
831 vector<double> dni(NN, 0.), dmus(NN, 0.);
833 for (
int i = 0; i < NN; ++i) {
835 xVector[NN + i] =
static_cast<int>(i == k);
838 solVector = decomp.solve(xVector);
840 for (
int i = 0; i < NN; ++i) {
841 dni[i] = solVector[i];
842 dmus[i] = solVector[NN + i];
843 m_dmusdmu[i][k] = dmus[i];
846 for (
int j = 0; j < NN; ++j) {
847 m_PrimCorrel[j][k] = dni[j];
851 for (
int i = 0; i < NN; ++i) {
852 m_wprim[i] = m_PrimCorrel[i][i];
853 if (m_densities[i] > 0.) m_wprim[i] *= m_Parameters.T / m_densities[i];
854 else m_wprim[i] = 1.;
867 for (
size_t i = 0; i < m_wprim.size(); ++i) {
874 m_FluctuationsCalculated =
true;
876 if (IsTemperatureDerivativesCalculated()) {
877 m_TemperatureDerivativesCalculated =
false;
883 if (!IsTemperatureDerivativesCalculated())
889 for (
int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
891 double fi = 1., dvi = 0.;
893 for (
int k = 0; k < m_TPS->ComponentsNumber(); ++k) {
896 fi -=
m_Virial[k][i] * m_densities[k];
907 ret += dvi * m_dndT[i];
917 int N = m_TPS->ComponentsNumber();
918 m_dndT = vector<double>(N, 0.);
919 m_dmusdT = vector<double>(N, 0.);
920 m_PrimChi2sdT = vector<vector<double>>(N, vector<double>(N, 0.));
922 double T = m_Parameters.T;
924 vector<double> nids(N, 0.), dniddTs(N, 0.), chi2ids(N, 0.), dchi2idsdT(N, 0.), chi3ids(N, 0.);
925 for (
int i = 0; i < N; ++i) {
935 int NN = m_densities.size();
937 MatrixXd densMatrix(2 * NN, 2 * NN);
938 VectorXd solVector(2 * NN), xVector(2 * NN);
942 for (
int i = 0; i<NN; ++i)
943 for (
int j = 0; j<NN; ++j) {
945 if (i == j) densMatrix(i, j) += 1.;
948 for (
int i = 0; i<NN; ++i)
949 for (
int j = 0; j<NN; ++j)
950 densMatrix(i, NN + j) = 0.;
952 for (
int i = 0; i<NN; ++i) {
953 densMatrix(i, NN + i) = 0.;
954 for (
int k = 0; k<NN; ++k) {
955 densMatrix(i, NN + i) +=
m_Virial[k][i] * m_densities[k];
957 densMatrix(i, NN + i) = (densMatrix(i, NN + i) - 1.) * chi2ids[i];
960 for (
int i = 0; i<NN; ++i)
961 for (
int j = 0; j<NN; ++j) {
965 for (
int i = 0; i<NN; ++i)
966 for (
int j = 0; j<NN; ++j) {
968 if (i == j) densMatrix(NN + i, NN + j) += 1.;
971 PartialPivLU<MatrixXd> decomp(densMatrix);
973 for(
int i=0;i<NN;++i) {
975 for(
int k=0;k<NN;++k) {
976 fi -=
m_Virial[k][i] * m_densities[k];
978 xVector[i] = fi * dniddTs[i];
979 xVector[NN + i] = 0.;
980 for(
int j = 0; j < NN; ++j) {
985 solVector = decomp.solve(xVector);
987 for(
int i=0;i<NN;++i) {
988 m_dndT[i] = solVector[i];
989 m_dmusdT[i] = solVector[NN+i];
994 if (IsSusceptibilitiesCalculated()) {
1012 for (
int j = 0; j < NN; ++j) {
1017 for (
int i = 0; i < NN; ++i) {
1037 vector<double> dfik = vector<double>(NN, 0.);
1039 for (
int k = 0; k < NN; ++k) {
1041 fi -=
m_Virial[k][i] * m_densities[k];
1046 for (
int k = 0; k < NN; ++k) {
1047 a1 += m_PrimCorrel[k][j] * dfik[k] * (dniddTs[i] + chi2ids[i] * m_dmusdT[i]);
1048 a1 += m_dmusdmu[i][j] * dfik[k] * m_dndT[k] * chi2ids[i];
1050 a1 += m_dmusdmu[i][j] * fi * (dchi2idsdT[i] + chi3ids[i] * m_dmusdT[i]);
1055 for (
int k = 0; k < NN; ++k) {
1056 a2 += m_dmusdmu[k][j] * (-
m_Virial[i][k]) * (dniddTs[k] + chi2ids[k] * m_dmusdT[k]);
1058 xVector[i + NN] = a2;
1061 solVector = decomp.solve(xVector);
1063 for (
int i = 0; i < NN; ++i) {
1070 m_TemperatureDerivativesCalculated =
true;
1074 m_TemperatureDerivativesCalculated =
true;
1082 for(
size_t i=0;i<m_densities.size();++i)
1083 if (m_densities[i]>0.)
1085 for(
size_t i=0;i<m_densities.size();++i)
1086 for(
size_t j=0;j<m_densities.size();++j)
1087 ret += -
m_Attr[i][j] * m_densities[i] * m_densities[j];
1090 for (
size_t i = 0; i < m_densities.size(); ++i) {
1092 for (
size_t j = 0; j < m_densities.size(); ++j) {
1093 ret += -tPid * m_densities[j] * m_Parameters.T *
m_VirialdT[j][i];
1094 ret += m_Parameters.T *
m_AttrdT[i][j] * m_densities[i] * m_densities[j];
1105 for(
size_t i=0;i<m_densities.size();++i)
1106 if (m_densities[i]>0.)
1110 for (
size_t i = 0; i < m_densities.size(); ++i) {
1112 for (
size_t j = 0; j < m_densities.size(); ++j) {
1113 ret += -tPid * m_densities[j] *
m_VirialdT[j][i];
1114 ret +=
m_AttrdT[i][j] * m_densities[i] * m_densities[j];
1125 for(
size_t i=0;i<m_densities.size();++i)
1127 for(
size_t i=0;i<m_densities.size();++i)
1128 for(
size_t j=0;j<m_densities.size();++j)
1129 ret += -
m_Attr[i][j] * m_densities[i] * m_densities[j];
1142 if (
id >= 0. &&
id <
static_cast<int>(
m_Virial.size()))
1150 if (i<0 || i >=
static_cast<int>(
m_Virial.size()) || j < 0 || j >=
static_cast<int>(
m_Virial.size()))
1157 if (i<0 || i >=
static_cast<int>(
m_Attr.size()) || j < 0 || j >=
static_cast<int>(
m_Attr.size()))
1164 if (i<0 || i >=
static_cast<int>(
m_VirialdT.size()) || j < 0 || j >=
static_cast<int>(
m_VirialdT.size()))
1171 if (i<0 || i >=
static_cast<int>(
m_AttrdT.size()) || j < 0 || j >=
static_cast<int>(
m_AttrdT.size()))
1176 std::vector<double> ThermalModelVDW::BroydenEquationsVDW::Equations(
const std::vector<double>& x)
1178 int NN = m_THM->Densities().size();
1179 vector<double> Ps(NN, 0.);
1180 for (
int i = 0; i < NN; ++i) {
1181 Ps[i] = m_THM->TPS()->Particles()[i].Density(m_THM->
Parameters(),
1188 vector<double> ns(NN, 0.);
1189 for (
int i = 0; i < NN; ++i) {
1190 ns[i] = m_THM->TPS()->Particles()[i].Density(m_THM->
Parameters(),
1197 vector<double> np = m_THM->ComputeNp(x, ns);
1200 vector<double> ret(
m_N, 0.);
1201 for (
size_t i = 0; i < ret.size(); ++i) {
1203 for (
int j = 0; j < NN; ++j)
1204 ret[i] += m_THM->VirialCoefficient(m_THM->m_MapFromdMuStar[i], j) * Ps[j]
1205 - (m_THM->AttractionCoefficient(m_THM->m_MapFromdMuStar[i], j)
1206 + m_THM->AttractionCoefficient(j, m_THM->m_MapFromdMuStar[i])) * np[j];
1211 std::vector<double> ThermalModelVDW::BroydenJacobianVDW::Jacobian(
const std::vector<double>& x)
1213 int NN = m_THM->m_densities.size();
1214 int NNdmu = m_THM->m_MapFromdMuStar.size();
1216 bool attrfl =
false;
1217 for (
int i = 0; i < NN; ++i) {
1218 for (
int j = 0; j < NN; ++j) {
1219 if (m_THM->AttractionCoefficient(i, j) != 0.0) {
1227 MatrixXd densMatrix(NNdmu, NNdmu);
1228 VectorXd solVector(NNdmu), xVector(NNdmu);
1230 std::vector<double> ret(NNdmu*NNdmu, 0.);
1232 vector<double> Ps(NN, 0.);
1233 for (
int i = 0; i<NN; ++i)
1234 Ps[i] = m_THM->TPS()->Particles()[i].Density(m_THM->Parameters(),
1237 m_THM->ChemicalPotential(i) + x[m_THM->m_MapTodMuStar[i]]
1240 vector<double> ns(NN, 0.);
1241 for (
int i = 0; i<NN; ++i)
1242 ns[i] = m_THM->TPS()->Particles()[i].Density(m_THM->Parameters(),
1245 m_THM->ChemicalPotential(i) + x[m_THM->m_MapTodMuStar[i]]
1248 vector<double> chi2s(NN, 0.);
1249 for (
int i = 0; i<NN; ++i)
1250 chi2s[i] = m_THM->TPS()->Particles()[i].chiDimensionfull(2, m_THM->Parameters(),
1252 m_THM->ChemicalPotential(i) + x[m_THM->m_MapTodMuStar[i]]
1255 for (
int i = 0; i < NNdmu; ++i) {
1256 for (
int j = 0; j < NNdmu; ++j) {
1257 densMatrix(i, j) = 0.;
1259 densMatrix(i, j) += 1.;
1261 for (
size_t m = 0; m < m_THM->m_dMuStarIndices[i].size(); ++m) {
1262 densMatrix(i, j) += m_THM->m_Virial[m_THM->m_MapFromdMuStar[j]][m_THM->m_dMuStarIndices[i][m]] * ns[m_THM->m_dMuStarIndices[i][m]];
1267 PartialPivLU<MatrixXd> decomp(densMatrix);
1270 for (
int kp = 0; kp < NNdmu; ++kp) {
1272 for (
size_t m = 0; m < m_THM->m_dMuStarIndices[kp].size(); ++m) {
1273 xVector[kp] += ns[m_THM->m_dMuStarIndices[kp][m]];
1278 solVector = decomp.solve(xVector);
1280 vector<double> ntildek(NNdmu, 0.);
1281 for (
int i = 0; i < NNdmu; ++i)
1282 ntildek[i] = solVector[i];
1284 vector<double> np(NN, 0.);
1285 for (
int i = 0; i < NN; ++i) {
1287 for (
int k = 0; k < NNdmu; ++k) {
1288 np[i] += m_THM->m_Virial[m_THM->m_MapFromdMuStar[k]][i] * solVector[k];
1290 np[i] = (1. - np[i]) * ns[i];
1293 for (
int kp = 0; kp < NNdmu; ++kp) {
1296 for (
int l = 0; l < NNdmu; ++l) {
1298 for (
size_t m = 0; m < m_THM->m_dMuStarIndices[l].size(); ++m) {
1299 int ti = m_THM->m_dMuStarIndices[l][m];
1300 if (m_THM->m_MapTodMuStar[ti] != kp)
1305 tmps = np[ti] / ns[ti];
1310 solVector = decomp.solve(xVector);
1311 for (
int i = 0; i < NNdmu; ++i)
1312 if (solVector[i] > 1.) solVector[i] = 1.;
1315 vector<double> dnjdmukp(NN, 0.);
1317 for (
int j = 0; j < NN; ++j) {
1318 for (
int kk = 0; kk < NNdmu; ++kk) {
1319 dnjdmukp[j] += -m_THM->m_Virial[m_THM->m_MapFromdMuStar[kk]][j] * solVector[kk] * ns[j];
1322 if (m_THM->m_MapTodMuStar[j] == kp) {
1325 tmps = np[j] / ns[j];
1332 for (
int k = 0; k < NNdmu; ++k) {
1334 ret[k*NNdmu + kp] += 1.;
1335 for (
size_t m = 0; m < m_THM->m_dMuStarIndices[kp].size(); ++m) {
1336 int tj = m_THM->m_dMuStarIndices[kp][m];
1337 ret[k*NNdmu + kp] += m_THM->m_Virial[m_THM->m_MapFromdMuStar[k]][tj] * ns[tj];
1341 for (
int j = 0; j < NN; ++j) {
1342 ret[k*NNdmu + kp] += -(m_THM->m_Attr[m_THM->m_MapFromdMuStar[k]][j] + m_THM->m_Attr[j][m_THM->m_MapFromdMuStar[k]]) * dnjdmukp[j];
1353 bool ThermalModelVDW::BroydenSolutionCriteriumVDW::IsSolved(
const std::vector<double>& x,
const std::vector<double>& f,
const std::vector<double>& xdelta)
const
1355 if (xdelta.size() == x.size()) {
1356 double maxdiff = 0.;
1357 for (
size_t i = 0; i < xdelta.size(); ++i) {
1358 maxdiff = std::max(maxdiff, fabs(xdelta[i]));
1359 maxdiff = std::max(maxdiff, fabs(f[i]));
1361 return (maxdiff < m_MaximumError);
1373 for (
int i1 = 0; i1 < model->TPS()->Particles().size(); ++i1) {
1374 for (
int i2 = 0; i2 < model->TPS()->Particles().size(); ++i2) {
Contains some functions to deal with excluded volumes.
map< string, double > params
virtual bool IsSolved(const std::vector< double > &x, const std::vector< double > &f, const std::vector< double > &xdelta=std::vector< double >()) const
int m_N
The number of equations.
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
Abstract base class for an HRG model implementation.
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...
virtual void SetVirial(int, int, double)
Set the excluded volume coefficient .
virtual void FillChemicalPotentials()
Sets the chemical potentials of all particles.
@ QvdW
Quantum van der Waals model.
virtual void SetChemicalPotentials(const std::vector< double > &chem=std::vector< double >(0))
Sets the chemical potentials of all particles.
double ChemicalPotential(int i) const
Chemical potential of particle species i.
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.
bool UseWidth() const
Whether finite resonance widths are considered.
virtual void CalculateDensities()
Calculates the primordial and total (after decays) densities of all species.
virtual void SetAttraction(int, int, double)
Set the vdW mean field attraction coefficient .
const ThermalModelParameters & Parameters() const
@ GCE
Grand canonical ensemble.
virtual void CalculateTwoParticleFluctuationsDecays()
Computes particle number correlations and fluctuations for all final-state particles which are marked...
virtual void SetChemicalPotentials(const std::vector< double > &chem=std::vector< double >(0))
Sets the chemical potentials of all particles.
std::vector< double > m_chiarb
virtual double CalculateEnergyDensityDerivativeT()
virtual void WriteInteractionParameters(const std::string &filename)
Write the QvdW interaction parameters to a file.
double VirialCoefficientdT(int i, int j) const
The temperature derivative of the eigenvolume parameter .
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
bool m_VDWComponentMapCalculated
Whether the mapping to components with the same VDW parameters has been calculated.
virtual std::vector< double > CalculateChargeFluctuations(const std::vector< double > &chgs, int order=4)
Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.
virtual std::vector< double > SearchSingleSolution(const std::vector< double > &muStarInit)
Uses the Broyden method with a provided initial guess to determine the shifted chemical potentials by...
virtual double MuShift(int id) const
The shift in the chemical potential of particle species i due to the QvdW interactions.
std::vector< double > m_DensitiesId
Vector of ideal gas densities with shifted chemical potentials.
virtual void ReadInteractionParameters(const std::string &filename)
Reads the QvdW interaction parameters from a file.
void FillVirialEV(const std::vector< std::vector< double > > &bij=std::vector< std::vector< double > >(0))
Same as FillVirial() but uses the matrix of excluded-volume coefficients as input instead of radii.
double AttractionCoefficient(int i, int j) const
QvdW mean field attraction coefficient .
virtual double CalculateEntropyDensity()
std::vector< int > m_MapTodMuStar
void CalculateVDWComponentsMap()
Partitions particles species into sets that have identical VDW parameters.
virtual double CalculateEnergyDensity()
void FillChemicalPotentials()
Sets the chemical potentials of all particles.
virtual double CalculatePressure()
std::vector< std::vector< double > > m_AttrdT
void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
bool m_TemperatureDependentAB
std::vector< double > ComputeNp(const std::vector< double > &dmustar)
virtual void FillAttraction(const std::vector< std::vector< double > > &aij=std::vector< std::vector< double > >(0))
double AttractionCoefficientdT(int i, int j) const
The temperature derivative of the QvdW attraction parameter .
std::vector< std::vector< int > > m_dMuStarIndices
virtual ~ThermalModelVDW(void)
Destroy the ThermalModelVDW object.
virtual void CalculateTemperatureDerivatives()
Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron...
std::vector< int > m_MapFromdMuStar
std::vector< std::vector< double > > m_Virial
std::vector< double > m_MuStar
Vector of the shifted chemical potentials.
void CalculateFluctuations()
Computes the fluctuation observables.
virtual double ParticleScalarDensity(int part)
std::vector< double > m_scaldens
Vector of scalar densities. Not used.
virtual void CalculatePrimordialDensities()
Calculates the primordial densities of all species.
std::vector< std::vector< double > > m_chi
std::vector< std::vector< double > > m_VirialdT
bool m_SearchMultipleSolutions
Whether multiple solutions are considered.
ThermalModelVDW(ThermalParticleSystem *TPS_, const ThermalModelParameters ¶ms=ThermalModelParameters())
Construct a new ThermalModelVDW object.
bool m_LastBroydenSuccessFlag
Whether Broyden's method was successfull.
double VirialCoefficient(int i, int j) const
Excluded volume coefficient .
virtual void FillVirial(const std::vector< double > &ri=std::vector< double >(0))
Fills the excluded volume coefficients based on the provided radii parameters for all species.
std::vector< std::vector< double > > m_Attr
Matrix of the attractive QvdW coefficients .
std::vector< double > SearchMultipleSolutions(int iters=300)
Uses the Broyden method with different initial guesses to look for different possible solutions of th...
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.
std::vector< std::string > split(const std::string &s, char delim)
double brr(double r1, double r2)
Computes the symmetric 2nd virial coefficient of the classical hard spheres equation of state from t...
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.
std::vector< std::vector< double > > GetBaryonBaryonInteractionMatrix(const ThermalParticleSystem *TPS, double param)
Returns the matrix of attraction and repulsion parameters for baryon-baryon and antibaryon-antibaryon...
void SetVDWHRGInteractionParameters(ThermalModelBase *model, double a, double b)
Sets vdW interactions for baryon-baryon and antibaryon-antibaryon pairs as in https://arxiv....
Structure containing all thermal parameters of the model.
Contains some extra mathematical functions used in the code.