Thermal-FIST 1.6
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
ThermalModelEVCrosstermsLegacy.cpp
Go to the documentation of this file.
1/*
2 * Thermal-FIST package
3 *
4 * Copyright (c) 2015-2019 Volodymyr Vovchenko
5 *
6 * GNU General Public License (GPLv3 or later)
7 */
9
10#ifdef USE_OPENMP
11#include <omp.h>
12#endif
13
14#include <vector>
15#include <cmath>
16#include <iostream>
17#include <iomanip>
18#include <algorithm>
19#include <fstream>
20#include <sstream>
21
22#include "HRGBase/xMath.h"
24
25#include <Eigen/Dense>
26
27using namespace Eigen;
28
29using namespace std;
30
31namespace thermalfist {
32
35 {
36 m_densitiesid.resize(m_TPS->Particles().size());
37 m_Volume = params.V;
38 m_Ps.resize(m_TPS->Particles().size());
39 FillVirial(std::vector<double>(m_TPS->Particles().size(), 0.));
40 m_TAG = "ThermalModelEVCrosstermsLegacy";
41
42 m_Ensemble = GCE;
43 m_InteractionModel = CrosstermsEV;
44 }
45
49
50 void ThermalModelEVCrosstermsLegacy::FillVirial(const std::vector<double> & ri) {
51 if (ri.size() != m_TPS->Particles().size()) {
52 std::cerr << "**WARNING** " << m_TAG << "::FillVirial(const std::vector<double> & ri): size " << static_cast<int>(ri.size()) << " of ri does not match number of hadrons " << static_cast<int>(m_TPS->Particles().size()) << " in the list" << std::endl;
53 return;
54 }
55 m_Virial.resize(m_TPS->Particles().size());
56 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
57 m_Virial[i].resize(m_TPS->Particles().size());
58 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j)
59 m_Virial[i][j] = CuteHRGHelper::brr(ri[i], ri[j]);
60 }
61
62 // Correction for non-diagonal terms R1=R2 and R2=0
63 std::vector< std::vector<double> > fVirialTmp = m_Virial;
64 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
65 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) {
66 if (i == j) m_Virial[i][j] = fVirialTmp[i][j];
67 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]);
68 }
69 }
70
72 {
73 m_Virial = std::vector< std::vector<double> >(m_TPS->Particles().size(), std::vector<double>(m_TPS->Particles().size(), 0.));
74
75 ifstream fin(filename.c_str());
76 char cc[2000];
77 while (!fin.eof()) {
78 fin.getline(cc, 2000);
79 string tmp = string(cc);
80 vector<string> elems = CuteHRGHelper::split(tmp, '#');
81 if (elems.size() < 1)
82 continue;
83 istringstream iss(elems[0]);
84 int pdgid1, pdgid2;
85 double b;
86 if (iss >> pdgid1 >> pdgid2 >> b) {
87 int ind1 = m_TPS->PdgToId(pdgid1);
88 int ind2 = m_TPS->PdgToId(pdgid2);
89 if (ind1 != -1 && ind2 != -1)
90 m_Virial[ind1][ind2] = b;
91 }
92 }
93 fin.close();
94 }
95
97 {
98 ofstream fout(filename.c_str());
99 fout << "# List of crossterms parameters to be used in the Crossterms excluded-volume HRG model"
100 << std::endl;
101 fout << "# Only particle pairs with a non-zero eigenvolume parameter are listed here"
102 << std::endl;
103 /*fout << "#" << std::setw(14) << "pdg_i"
104 << std::setw(15) << "pdg_j"
105 << std::setw(15) << "b_{ij}[fm^3]"
106 << std::endl;*/
107 fout << "#" << " " << "pdg_i"
108 << " " << "pdg_j"
109 << " " << "b_{ij}[fm^3]"
110 << std::endl;
111 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
112 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) {
113 if (m_Virial[i][j] != 0.) {
114 //fout << std::setw(15) << m_TPS->Particle(i).PdgId();
115 //fout << std::setw(15) << m_TPS->Particle(j).PdgId();
116 //fout << std::setw(15) << m_Virial[i][j];
117 fout << " " << m_TPS->Particle(i).PdgId();
118 fout << " " << m_TPS->Particle(j).PdgId();
119 fout << " " << m_Virial[i][j];
120 fout << std::endl;
121 }
122 }
123 }
124 fout.close();
125 }
126
128 FillVirial(vector<double>(m_TPS->Particles().size(), rad));
129 }
130
132 if (i < 0 || i >= static_cast<int>(m_Virial.size()) || j < 0 || j > static_cast<int>(m_Virial.size()))
133 return 0.;
134 return m_Virial[i][j];
135 }
136
137 void ThermalModelEVCrosstermsLegacy::SetVirial(int i, int j, double b) {
138 if (i >= 0 && i < static_cast<int>(m_Virial.size()) && j >= 0 && j < static_cast<int>(m_Virial.size()))
139 m_Virial[i][j] = b;
140 else std::cerr << "**WARNING** Index overflow in ThermalModelEVCrosstermsLegacy::SetVirial" << std::endl;
141 }
142
145 m_densitiesid.resize(m_TPS->Particles().size());
146 FillVirial();
147 }
148
149 double ThermalModelEVCrosstermsLegacy::DensityId(int i, const std::vector<double>& pstars)
150 {
151 double dMu = 0.;
152 if (pstars.size() == m_TPS->Particles().size())
153 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j)
154 dMu += -m_Virial[i][j] * pstars[j];
155 else
156 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j)
157 dMu += -m_Virial[i][j] * m_Ps[j];
158
159 return m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i] + dMu);
160 }
161
162 double ThermalModelEVCrosstermsLegacy::Pressure(int i, const std::vector<double>& pstars)
163 {
164 double dMu = 0.;
165 if (pstars.size() == m_TPS->Particles().size())
166 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j)
167 dMu += -m_Virial[i][j] * pstars[j];
168 else
169 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j)
170 dMu += -m_Virial[i][j] * m_Ps[j];
171
172 return m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i] + dMu);
173 }
174
176 double dMu = 0.;
177 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) dMu += -m_Virial[i][j] * m_Ps[j];
178
179 return m_TPS->Particles()[i].ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i] + dMu);
180 }
181
183 double dMu = 0.;
184 dMu += -m_Virial[i][i] * P;
185
186 return m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i] + dMu);
187 }
188
189
191 double ret = 0.;
192 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
193 ret += PartialPressureDiagonal(i, P);
194 return ret;
195 }
196
198 BroydenEquationsCRSDEV eqs(this);
199 BroydenJacobian jac(&eqs);
200 jac.SetDx(1.0E-8);
201 Broyden broydn(&eqs, &jac);
203
204 m_Pressure = 0.;
205 double x0 = m_Pressure;
206 std::vector<double> x(1, x0);
207
208 x = broydn.Solve(x, &crit);
209
210 m_Pressure = x[0];
211 for (size_t i = 0; i < m_Ps.size(); ++i)
213 }
214
215
217 if (resetPartials) {
218 m_Ps.resize(m_TPS->Particles().size());
219 for (size_t i = 0; i < m_Ps.size(); ++i) m_Ps[i] = 0.;
221 }
222
223 BroydenEquationsCRS eqs(this);
224 BroydenJacobianCRS jac(this);
225 Broyden broydn(&eqs, &jac);
226 BroydenSolutionCriteriumCRS crit(this);
227
228 m_Ps = broydn.Solve(m_Ps, &crit);
229 m_Pressure = 0.;
230 for (size_t i = 0; i < m_Ps.size(); ++i)
231 m_Pressure += m_Ps[i];
232
233 if (broydn.Iterations() == broydn.MaxIterations())
234 m_LastCalculationSuccessFlag = false;
235 else m_LastCalculationSuccessFlag = true;
236
237 m_MaxDiff = broydn.MaxDifference();
238 }
239
241 m_FluctuationsCalculated = false;
242
243 map< vector<double>, int> m_MapEVcomponent;
244
245 {
246 int NN = m_densities.size();
247 m_MapToEVComponent.resize(NN);
248 m_MapFromEVComponent.clear();
249 m_MapEVcomponent.clear();
250 m_EVComponentIndices.clear();
251
252 int tind = 0;
253 for (int i = 0; i < NN; ++i) {
254 vector<double> EVParam(0);
255 for (int j = 0; j < NN; ++j) {
256 EVParam.push_back(m_Virial[i][j]);
257 EVParam.push_back(m_Virial[j][i]);
258 }
259
260 if (m_MapEVcomponent.count(EVParam) == 0) {
261 m_MapEVcomponent[EVParam] = tind;
262 m_MapToEVComponent[i] = tind;
263 m_MapFromEVComponent.push_back(i);
264 m_EVComponentIndices.push_back(vector<int>(1, i));
265 tind++;
266 }
267 else {
268 m_MapToEVComponent[i] = m_MapEVcomponent[EVParam];
269 m_EVComponentIndices[m_MapEVcomponent[EVParam]].push_back(i);
270 }
271 }
272 }
273
274 // Pressure
276 vector<double> tN(m_densities.size());
277 for (size_t i = 0; i < m_Ps.size(); ++i)
278 tN[i] = DensityId(i);
279
280 // Densities
281
282 int NN = m_densities.size();
283
284 MatrixXd densMatrix(NN, NN);
285 VectorXd solVector(NN), xVector(NN);
286
287 for (int i = 0; i < NN; ++i)
288 for (int j = 0; j < NN; ++j) {
289 densMatrix(i, j) = m_Virial[i][j] * tN[i];
290 if (i == j) densMatrix(i, j) += 1.;
291 }
292
293 PartialPivLU<MatrixXd> decomp(densMatrix);
294
295 for (int i = 0; i < NN; ++i) xVector[i] = 0.;
296 for (int i = 0; i < NN; ++i) {
297 xVector[i] = tN[i];
298 solVector = decomp.solve(xVector);
299 if (1) {
300 m_densities[i] = 0.;
301 for (int j = 0; j < NN; ++j)
302 m_densities[i] += solVector[j];
303 }
304 else {
305 std::cerr << "**WARNING** Could not recover m_densities from partial pressures!\n";
306 }
307 xVector[i] = 0.;
308 }
309
310 std::vector<double> fEntropyP(m_densities.size());
311 for (int i = 0; i < NN; ++i) {
312 double dMu = 0.;
313 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) dMu += -m_Virial[i][j] * m_Ps[j];
314 xVector[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::EntropyDensity, m_UseWidth, m_Chem[i] + dMu);
315 }
316
317 solVector = decomp.solve(xVector);
318
319 if (1) {
321 for (int i = 0; i < NN; ++i)
322 m_TotalEntropyDensity += solVector[i];
323 }
324 else {
325 std::cerr << "**WARNING** Could not recover m_densities from partial pressures!\n";
326 }
327
328 m_Calculated = true;
330 }
331
333 // Pressure
334 SolvePressure(false);
335 vector<double> tN(m_densities.size());
336 for (size_t i = 0; i < m_Ps.size(); ++i)
337 tN[i] = DensityId(i);
338
339 // Densities
340
341 int NN = m_densities.size();
342
343 MatrixXd densMatrix(NN, NN);
344 VectorXd solVector(NN), xVector(NN);
345
346 for (int i = 0; i < NN; ++i)
347 for (int j = 0; j < NN; ++j) {
348 densMatrix(i, j) = m_Virial[i][j] * tN[i];
349 if (i == j) densMatrix(i, j) += 1.;
350 }
351
352 PartialPivLU<MatrixXd> decomp(densMatrix);
353
354 for (int i = 0; i < NN; ++i) xVector[i] = 0.;
355 for (int i = 0; i < NN; ++i) {
356 xVector[i] = tN[i];//m_Ps[i] / m_Parameters.T;
357 //solVector = lu.Solve(xVector, ok);
358 solVector = decomp.solve(xVector);
359 //if (ok) {
360 if (1) {
361 //if (decomp.info()==Eigen::Success) {
362 m_densities[i] = 0.;
363 for (int j = 0; j < NN; ++j)
364 m_densities[i] += solVector[j];
365 }
366 else {
367 std::cerr << "**WARNING** Could not recover m_densities from partial pressures!\n";
368 return;
369 }
370 xVector[i] = 0.;
371 }
372
373 std::vector<double> fEntropyP(m_densities.size());
374 for (int i = 0; i < NN; ++i) {
375 double dMu = 0.;
376 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) dMu += -m_Virial[i][j] * m_Ps[j];
377 xVector[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::EntropyDensity, m_UseWidth, m_Chem[i] + dMu);
378 }
379
380 solVector = decomp.solve(xVector);
381
382 //if (ok) {
383 if (1) {
384 //if (decomp.info()==Eigen::Success) {
386 for (int i = 0; i < NN; ++i)
387 m_TotalEntropyDensity += solVector[i];
388 }
389 else {
390 std::cerr << "**WARNING** Could not recover m_densities from partial pressures!\n";
391 return;
392 }
393
394 // Decays
395
397
398 m_Calculated = true;
399 }
400
402 m_Ps.resize(m_TPS->Particles().size());
403 for (size_t i = 0; i < m_Ps.size(); ++i)
404 m_Ps[i] = 0.;
406 vector<double> Pstmp = m_Ps;
407 int iter = 0;
408 double maxdiff = 0.;
409 for (iter = 0; iter < 1000; ++iter)
410 {
411 maxdiff = 0.;
412 for (size_t i = 0; i < m_Ps.size(); ++i) {
413 Pstmp[i] = Pressure(i);
414 maxdiff = max(maxdiff, fabs((Pstmp[i] - m_Ps[i]) / Pstmp[i]));
415 }
416 m_Ps = Pstmp;
417 //cout << iter << "\t" << maxdiff << "\n";
418 if (maxdiff < 1.e-10) break;
419 }
420 if (iter == 1000) std::cerr << iter << "\t" << maxdiff << "\n";
421 m_Pressure = 0.;
422 for (size_t i = 0; i < m_Ps.size(); ++i)
423 m_Pressure += m_Ps[i];
424 }
425
427 // Pressure
429
430 int NN = m_densities.size();
431 vector<double> tN(NN);
432 for (int i = 0; i < NN; ++i) tN[i] = DensityId(i);
433
434 MatrixXd densMatrix(NN, NN);
435 for (int i = 0; i < NN; ++i)
436 for (int j = 0; j < NN; ++j) {
437 densMatrix(i, j) = m_Virial[j][i] * tN[i];
438 if (i == j) densMatrix(i, j) += 1.;
439 }
440 //densMatrix.SetMatrixArray(matr.GetArray());
441
442 VectorXd solVector(NN), xVector(NN);
443 for (int i = 0; i < NN; ++i) xVector[i] = tN[i];
444
445 PartialPivLU<MatrixXd> decomp(densMatrix);
446
447 solVector = decomp.solve(xVector);
448
449 if (1) {
450 //if (decomp.info()==Eigen::Success) {
451 for (int i = 0; i < NN; ++i)
452 m_densities[i] = solVector[i];
453 }
454 else {
455 std::cerr << "**WARNING** Could not recover m_densities from partial pressures!\n";
456 return;
457 }
458
459 // Entropy
460 for (int i = 0; i < NN; ++i)
461 for (int j = 0; j < NN; ++j) {
462 densMatrix(i, j) = m_Virial[i][j] * tN[i];
463 if (i == j) densMatrix(i, j) += 1.;
464 }
465
466 std::vector<double> fEntropyP(m_densities.size());
467 for (int i = 0; i < NN; ++i) {
468 double dMu = 0.;
469 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) dMu += -m_Virial[i][j] * m_Ps[j];
470 xVector[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::EntropyDensity, m_UseWidth, m_Chem[i] + dMu);
471 }
472
473 decomp = PartialPivLU<MatrixXd>(densMatrix);
474 solVector = decomp.solve(xVector);
475
476 if (1) {
477 //if (decomp.info()==Eigen::Success) {
479 for (int i = 0; i < NN; ++i)
480 m_TotalEntropyDensity += solVector[i];
481 }
482 else {
483 std::cerr << "Could not recover entropy m_densities from partial pressures!\n";
484 }
485
486 m_Calculated = true;
487 }
488
490 int NN = m_densities.size();
491 vector<double> tN(NN);// , tW(NN);
492 for (int i = 0; i < NN; ++i) tN[i] = DensityId(i);
493 //for (int i = 0; i < NN; ++i) tW[i] = ScaledVarianceId(i);
494 vector<double> chi2id(NN);
495 for (int i = 0; i < NN; ++i) {
496 double dMu = 0.;
497 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) dMu += -m_Virial[i][j] * m_Ps[j];
498 chi2id[i] = m_densities[i] / tN[i] * m_TPS->Particles()[i].chiDimensionfull(2, m_Parameters, m_UseWidth, m_Chem[i] + dMu) * xMath::GeVtoifm3();
499 }
500
501 MatrixXd densMatrix(NN, NN);
502 VectorXd solVector(NN), xVector(NN), xVector2(NN);
503
504 for (int i = 0; i < NN; ++i)
505 for (int j = 0; j < NN; ++j) {
506 densMatrix(i, j) = m_Virial[i][j] * tN[i];
507 if (i == j) densMatrix(i, j) += 1.;
508 }
509
510 PartialPivLU<MatrixXd> decomp(densMatrix);
511
512 vector< vector<double> > ders, coefs;
513
514 ders.resize(NN);
515 coefs.resize(NN);
516
517 for (int i = 0; i < NN; ++i) {
518 ders[i].resize(NN);
519 coefs[i].resize(NN);
520 }
521
522 for (int i = 0; i < NN; ++i) xVector[i] = 0.;
523 for (int i = 0; i < NN; ++i) {
524 xVector[i] = tN[i];
525 solVector = decomp.solve(xVector);
526 if (1) {
527 //if (decomp.info()==Eigen::Success) {
528 for (int j = 0; j < NN; ++j) {
529 ders[j][i] = solVector[j];
530 }
531
532 for (int l = 0; l < NN; ++l) {
533 coefs[l][i] = 0.;
534 for (int k = 0; k < NN; ++k) {
535 coefs[l][i] += -m_Virial[l][k] * ders[k][i];
536 }
537 if (l == i) coefs[l][i] += 1.;
538 }
539 }
540 else {
541 std::cerr << "**WARNING** Could not recover fluctuations!\n";
542 }
543 xVector[i] = 0.;
544 }
545
546
547 m_PrimCorrel.resize(NN);
548 for (int i = 0; i < NN; ++i) m_PrimCorrel[i].resize(NN);
549 m_TotalCorrel = m_PrimCorrel;
550
551 for (int i = 0; i < NN; ++i)
552 for (int j = i; j < NN; ++j) {
553 for (int l = 0; l < NN; ++l)
554 //xVector[l] = tN[l] / m_Parameters.T * tW[l] * coefs[l][i] * coefs[l][j];
555 xVector[l] = chi2id[l] * coefs[l][i] * coefs[l][j];
556 solVector = decomp.solve(xVector);
557 if (1) {
558 //if (decomp.info()==Eigen::Success) {
559 m_PrimCorrel[i][j] = 0.;
560 for (int k = 0; k < NN; ++k) {
561 m_PrimCorrel[i][j] += solVector[k];
562 }
563 m_PrimCorrel[j][i] = m_PrimCorrel[i][j];
564 }
565 else {
566 std::cerr << "**WARNING** Could not recover fluctuations!\n";
567 }
568 }
569
570 //cout << "Primaries solved!\n";
571
572 for (int i = 0; i < NN; ++i) {
573 m_wprim[i] = m_PrimCorrel[i][i];
574 if (m_densities[i] > 0.) m_wprim[i] *= m_Parameters.T / m_densities[i];
575 else m_wprim[i] = 1.;
576 if (m_wprim[i] != m_wprim[i]) m_wprim[i] = 1.;
577 }
578
579 m_TwoParticleCorrelationsCalculated = true;
580 }
581
582 // TODO include correlations
589
590 m_FluctuationsCalculated = true;
591
592 for (size_t i = 0; i < m_wprim.size(); ++i) {
593 m_skewprim[i] = 1.;
594 m_kurtprim[i] = 1.;
595 m_skewtot[i] = 1.;
596 m_kurttot[i] = 1.;
597 }
598 }
599
600 std::vector<double> ThermalModelEVCrosstermsLegacy::CalculateChargeFluctuations(const std::vector<double>& chgs, int order, bool dimensionfull)
601 {
602 vector<double> ret(order + 1, 0.);
603
604 // chi1
605 for (size_t i = 0; i < m_densities.size(); ++i)
606 ret[0] += chgs[i] * m_densities[i];
607
608 if (!dimensionfull)
609 ret[0] /= pow(m_Parameters.T * xMath::GeVtoifm(), 3);
610 else
611 ret[0] /= pow(xMath::GeVtoifm(), 3);
612
613 if (order < 2) return ret;
614 // Preparing matrix for system of linear equations
615 int NN = m_densities.size();
616
617 vector<double> MuStar(NN, 0.);
618 for (int i = 0; i < NN; ++i) {
619 MuStar[i] = m_Chem[i] + MuShift(i);
620 }
621
622 MatrixXd densMatrix(2 * NN, 2 * NN);
623 VectorXd solVector(2 * NN), xVector(2 * NN);
624
625 vector<double> DensitiesId(m_densities.size()), chi2id(m_densities.size());
626 for (int i = 0; i < NN; ++i) {
627 DensitiesId[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, MuStar[i]);
628 chi2id[i] = m_TPS->Particles()[i].chi(2, m_Parameters, m_UseWidth, MuStar[i]);
629 }
630
631 for (int i = 0; i < NN; ++i)
632 for (int j = 0; j < NN; ++j) {
633 densMatrix(i, j) = m_Virial[j][i] * DensitiesId[i];
634 if (i == j) densMatrix(i, j) += 1.;
635 }
636
637 for (int i = 0; i < NN; ++i)
638 for (int j = 0; j < NN; ++j)
639 densMatrix(i, NN + j) = 0.;
640
641 for (int i = 0; i < NN; ++i) {
642 densMatrix(i, NN + i) = 0.;
643 for (int k = 0; k < NN; ++k) {
644 densMatrix(i, NN + i) += m_Virial[k][i] * m_densities[k];
645 }
646 densMatrix(i, NN + i) = (densMatrix(i, NN + i) - 1.) * chi2id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T;
647 }
648
649 for (int i = 0; i < NN; ++i)
650 for (int j = 0; j < NN; ++j) {
651 densMatrix(NN + i, NN + j) = m_Virial[i][j] * DensitiesId[j];
652 if (i == j) densMatrix(NN + i, NN + j) += 1.;
653 }
654
655
656 PartialPivLU<MatrixXd> decomp(densMatrix);
657
658 // chi2
659 vector<double> dni(NN, 0.), dmus(NN, 0.);
660
661 for (int i = 0; i < NN; ++i) {
662 xVector[i] = 0.;
663 xVector[NN + i] = chgs[i];
664 }
665
666 solVector = decomp.solve(xVector);
667
668 for (int i = 0; i < NN; ++i) {
669 dni[i] = solVector[i];
670 dmus[i] = solVector[NN + i];
671 }
672
673 for (int i = 0; i < NN; ++i)
674 ret[1] += chgs[i] * dni[i];
675
676 if (!dimensionfull)
677 ret[1] /= pow(m_Parameters.T, 2) * pow(xMath::GeVtoifm(), 3);
678 else
679 ret[1] /= pow(xMath::GeVtoifm(), 3);
680
681 if (order < 3) return ret;
682 // chi3
683 vector<double> d2ni(NN, 0.), d2mus(NN, 0.);
684
685 vector<double> chi3id(m_densities.size());
686 for (int i = 0; i < NN; ++i)
687 chi3id[i] = m_TPS->Particles()[i].chi(3, m_Parameters, m_UseWidth, MuStar[i]);
688
689 for (int i = 0; i < NN; ++i) {
690 xVector[i] = 0.;
691
692 double tmp = 0.;
693 for (int j = 0; j < NN; ++j) tmp += m_Virial[j][i] * dni[j];
694 tmp = -2. * tmp * chi2id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[i];
695 xVector[i] += tmp;
696
697 tmp = 0.;
698 for (int j = 0; j < NN; ++j) tmp += m_Virial[j][i] * m_densities[j];
699 tmp = -(tmp - 1.) * chi3id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * dmus[i] * dmus[i];
700 xVector[i] += tmp;
701 }
702 for (int i = 0; i < NN; ++i) {
703 xVector[NN + i] = 0.;
704
705 double tmp = 0.;
706 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];
707
708 xVector[NN + i] = tmp;
709 }
710
711 solVector = decomp.solve(xVector);
712
713 for (int i = 0; i < NN; ++i) {
714 d2ni[i] = solVector[i];
715 d2mus[i] = solVector[NN + i];
716 }
717
718 for (int i = 0; i < NN; ++i)
719 ret[2] += chgs[i] * d2ni[i];
720
721 if (!dimensionfull)
722 ret[2] /= m_Parameters.T * pow(xMath::GeVtoifm(), 3);
723 else
724 ret[2] /= pow(xMath::GeVtoifm(), 3);
725
726
727 if (order < 4) return ret;
728
729 // chi4
730 vector<double> d3ni(NN, 0.), d3mus(NN, 0.);
731
732 vector<double> chi4id(m_densities.size());
733 for (int i = 0; i < NN; ++i)
734 chi4id[i] = m_TPS->Particles()[i].chi(4, m_Parameters, m_UseWidth, MuStar[i]);
735
736 vector<double> dnis(NN, 0.);
737 for (int i = 0; i < NN; ++i) {
738 dnis[i] = chi2id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * dmus[i];
739 }
740
741 vector<double> d2nis(NN, 0.);
742 for (int i = 0; i < NN; ++i) {
743 d2nis[i] = chi3id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * dmus[i] * dmus[i] +
744 chi2id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T * d2mus[i];
745 }
746
747 for (int i = 0; i < NN; ++i) {
748 xVector[i] = 0.;
749
750 double tmp = 0.;
751 for (int j = 0; j < NN; ++j) tmp += m_Virial[j][i] * dni[j];
752 tmp = -3. * tmp * d2nis[i];
753 xVector[i] += tmp;
754
755 tmp = 0.;
756 for (int j = 0; j < NN; ++j) tmp += m_Virial[j][i] * d2ni[j];
757 tmp = -3. * tmp * dnis[i];
758 xVector[i] += tmp;
759
760 double tmps = 0.;
761 for (int j = 0; j < NN; ++j) tmps += m_Virial[j][i] * m_densities[j];
762
763 tmp = -(tmps - 1.) * chi3id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * d2mus[i] * 3. * dmus[i];
764 xVector[i] += tmp;
765
766 tmp = -(tmps - 1.) * chi4id[i] * pow(xMath::GeVtoifm(), 3) * dmus[i] * dmus[i] * dmus[i];
767 xVector[i] += tmp;
768 }
769 for (int i = 0; i < NN; ++i) {
770 xVector[NN + i] = 0.;
771
772 double tmp = 0.;
773 for (int j = 0; j < NN; ++j) tmp += -2. * m_Virial[i][j] * d2mus[j] * dnis[j];
774 xVector[NN + i] += tmp;
775
776 tmp = 0.;
777 for (int j = 0; j < NN; ++j) tmp += -m_Virial[i][j] * dmus[j] * d2nis[j];
778 xVector[NN + i] += tmp;
779 }
780
781 solVector = decomp.solve(xVector);
782
783 for (int i = 0; i < NN; ++i) {
784 d3ni[i] = solVector[i];
785 d3mus[i] = solVector[NN + i];
786 }
787
788 for (int i = 0; i < NN; ++i)
789 ret[3] += chgs[i] * d3ni[i];
790
791 ret[3] /= pow(xMath::GeVtoifm(), 3);
792
793 return ret;
794 }
795
797 double ret = 0.;
798 ret += m_Parameters.T * CalculateEntropyDensity();
799 ret += -CalculatePressure();
800 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
801 ret += m_Chem[i] * m_densities[i];
802 return ret;
803 }
804
809
811 if (!m_Calculated) CalculateDensities();
812 return m_Pressure;
813 }
814
815
817 {
818 if (id >= 0. && id < static_cast<int>(m_Virial.size())) {
819 double dMu = 0.;
820 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j)
821 dMu += -m_Virial[id][j] * m_Ps[j];
822 return dMu;
823 }
824 else
825 return 0.0;
826 }
827
828 std::vector<double> ThermalModelEVCrosstermsLegacy::BroydenEquationsCRS::Equations(const std::vector<double>& x)
829 {
830 std::vector<double> ret(m_N);
831 for (size_t i = 0; i < x.size(); ++i)
832 ret[i] = x[i] - m_THM->Pressure(i, x);
833 return ret;
834 }
835
836 std::vector<double> ThermalModelEVCrosstermsLegacy::BroydenJacobianCRS::Jacobian(const std::vector<double>& x)
837 {
838 int N = x.size();
839
840 vector<double> tN(N);
841 for (int i = 0; i < N; ++i) {
842 tN[i] = m_THM->DensityId(i, x);
843 }
844
845 vector<double> ret(N*N, 0.);
846 for (int i = 0; i < N; ++i) {
847 for (int j = 0; j < N; ++j) {
848 if (i == j)
849 ret[i*N + j] += 1.;
850 ret[i*N + j] += m_THM->VirialCoefficient(i, j) * tN[i];
851 }
852 }
853
854 return ret;
855 }
856
857 bool ThermalModelEVCrosstermsLegacy::BroydenSolutionCriteriumCRS::IsSolved(const std::vector<double>& x, const std::vector<double>& f, const std::vector<double>& /*xdelta*/) const
858 {
859 double maxdiff = 0.;
860 for (size_t i = 0; i < x.size(); ++i) {
861 maxdiff = std::max(maxdiff, fabs(f[i]) / x[i]);
862 }
863 return (maxdiff < m_MaximumError);
864 }
865
866 std::vector<double> ThermalModelEVCrosstermsLegacy::BroydenEquationsCRSDEV::Equations(const std::vector<double>& x)
867 {
868 std::vector<double> ret(1);
869 ret[0] = x[0] - m_THM->PressureDiagonalTotal(x[0]);
870 return ret;
871 }
872} // namespace thermalfist
Contains some functions to deal with excluded volumes.
map< string, double > params
Sub-class where it is determined whether the required accuracy is achieved in the Broyden's method.
Definition Broyden.h:144
int m_N
The number of equations.
Definition Broyden.h:66
Class implementing the Broyden method to solve a system of non-linear equations.
Definition Broyden.h:132
double MaxDifference() const
Definition Broyden.h:239
virtual std::vector< double > Solve(const std::vector< double > &x0, BroydenSolutionCriterium *solcrit=NULL, int max_iterations=MAX_ITERS)
Definition Broyden.cpp:83
int Iterations() const
Definition Broyden.h:229
int MaxIterations() const
Definition Broyden.h:234
Class which implements calculation of the Jacobian needed for the Broyden's method.
Definition Broyden.h:77
void SetDx(double dx)
Set the finite variable difference value used for calculating the Jacobian numerically.
Definition Broyden.h:114
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...
@ CrosstermsEV
Crossterms 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 &params=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...
virtual double Pressure(int i, const std::vector< double > &pstars=std::vector< double >())
Calculate the ideal gas pressure of particle species i for the given values of partial pressures.
virtual void WriteInteractionParameters(const std::string &filename)
Write the QvdW interaction parameters to a file.
double VirialCoefficient(int i, int j) const
Excluded volume coefficient .
virtual std::vector< double > CalculateChargeFluctuations(const std::vector< double > &chgs, int order=4, bool dimensionfull=false)
Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.
void SetVirial(int i, int j, double b)
Set the excluded volume coefficient .
void CalculateFluctuations()
Computes the fluctuation observables.
double PressureDiagonalTotal(double P)
The total pressure for the given input pressure in the diagonal model.
virtual double MuShift(int i) const
The shift in the chemical potential of particle species i due to the excluded volume interactions.
virtual ~ThermalModelEVCrosstermsLegacy(void)
Destroy the ThermalModelEVCrossterms object.
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.
virtual double DensityId(int i, const std::vector< double > &pstars=std::vector< double >())
Calculate the ideal gas density of particle species i for the given values of partial pressures.
void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
virtual void SolvePressure(bool resetPartials=true)
Solves the system of transcdental equations for the pressure using the Broyden's method.
virtual void ReadInteractionParameters(const std::string &filename)
Reads the QvdW interaction parameters from a file.
virtual void SolveDiagonal()
Solves the transcendental equation of the corresponding diagonal EV model.
void SetRadius(double rad)
Set the same excluded volume radius parameter for all species.
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
ThermalModelEVCrosstermsLegacy(ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
Construct a new ThermalModelEVCrossterms object.
double PartialPressureDiagonal(int i, double P)
The "partial pressure" of hadron species i for the given total pressure in the diagonal model.
virtual void CalculatePrimordialDensities()
Calculates the primordial densities of all species.
double ScaledVarianceId(int ind)
Calculate the ideal gas scaled variance of particle species i number fluctuations for the given value...
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 .
Definition xMath.h:25
constexpr double GeVtoifm3()
A constant to transform GeV into fm .
Definition xMath.h:31
The main namespace where all classes and functions of the Thermal-FIST library reside.
Definition CosmicEoS.h:9
Structure containing all thermal parameters of the model.
Contains some extra mathematical functions used in the code.