Thermal-FIST 1.5
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
ThermalModelEVDiagonal.cpp
Go to the documentation of this file.
1/*
2 * Thermal-FIST package
3 *
4 * Copyright (c) 2014-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 <fstream>
19#include <sstream>
20
21#include "HRGBase/xMath.h"
23
24#include <Eigen/Dense>
25
26using namespace Eigen;
27
28using namespace std;
29
30namespace thermalfist {
31
34 {
35 m_densitiesid.resize(m_TPS->Particles().size());
36 m_v.resize(m_TPS->Particles().size());
37 m_Volume = params.V;
38 m_TAG = "ThermalModelEVDiagonal";
39
40 m_Ensemble = GCE;
41 m_InteractionModel = DiagonalEV;
42 }
43
44
48
49
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) {
54 m_v[i] = CuteHRGHelper::vr(rad);
55 }
56 }
57
58 void ThermalModelEVDiagonal::FillVirial(const std::vector<double>& ri)
59 {
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;
62 return;
63 }
64 m_v.resize(m_TPS->Particles().size());
65 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
66 m_v[i] = CuteHRGHelper::vr(ri[i]);
67 }
68
69 void ThermalModelEVDiagonal::FillVirialEV(const std::vector<double>& vi)
70 {
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;
73 return;
74 }
75 m_v = vi;
76 }
77
78 void ThermalModelEVDiagonal::ReadInteractionParameters(const std::string & filename)
79 {
80 m_v = std::vector<double>(m_TPS->Particles().size(), 0.);
81
82 ifstream fin(filename.c_str());
83 char cc[2000];
84 while (!fin.eof()) {
85 fin.getline(cc, 2000);
86 string tmp = string(cc);
87 vector<string> elems = CuteHRGHelper::split(tmp, '#');
88 if (elems.size() < 1)
89 continue;
90 istringstream iss(elems[0]);
91 long long pdgid;
92 double b;
93 if (iss >> pdgid >> b) {
94 int ind = m_TPS->PdgToId(pdgid);
95 if (ind != -1)
96 m_v[ind] = b;
97 }
98 }
99 fin.close();
100 }
101
102 void ThermalModelEVDiagonal::WriteInteractionParameters(const std::string & filename)
103 {
104 ofstream fout(filename.c_str());
105 fout << "# List of eigenvolume parameters to be used in the Diagonal excluded-volume HRG model"
106 << std::endl;
107 fout << "# Only particles with a non-zero eigenvolume parameter are listed here"
108 << std::endl;
109 fout << "#" << std::setw(14) << "pdg"
110 << std::setw(15) << "v_i[fm^3]"
111 << std::endl;
112 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
113 if (m_v[i] != 0.) {
114 fout << std::setw(15) << m_TPS->Particle(i).PdgId();
115 fout << std::setw(15) << m_v[i];
116 fout << std::endl;
117 }
118 }
119 fout.close();
120 }
121
123 {
124 if (i < 0 || i >= static_cast<int>(m_v.size()))
125 return 0.;
126 return m_v[i];
127 }
128
129
132 m_densitiesid.resize(m_TPS->Particles().size());
133 }
134
135
137 double dMu = -m_v[i] * Pressure;
138
139 return m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i] + dMu);
140 }
141
143 double dMu = -m_v[i] * Pressure;
144
145 return m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i] + dMu);
146 }
147
149 double dMu = -m_v[i] * Pressure;
150
151 return m_TPS->Particles()[i].ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i] + dMu);
152 }
153
155 double ret = 0.;
156
157 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
158 double dMu = -m_v[i] * P;
159 ret += m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i] + dMu);
160 }
161
162 return ret;
163 }
164
166 BroydenEquationsDEV eqs(this);
167 BroydenJacobianDEV jac(this);
168 Broyden broydn(&eqs, &jac);
169 BroydenSolutionCriteriumDEV crit(this, 1.0E-8);
170
171 m_Pressure = Pressure(0.);
172 double mnc = pow(m_Parameters.T, 4.) * pow(xMath::GeVtoifm(), 3.);
173 if (m_Parameters.T == 0.0)
174 mnc = pow(xMath::mnucleon(), 4.) * pow(xMath::GeVtoifm(), 3.);
175 eqs.SetMnc(mnc);
176 jac.SetMnc(mnc);
177 double x0 = log(m_Pressure / mnc);
178 std::vector<double> x(1, x0);
179
180 x = broydn.Solve(x, &crit);
181
182
183
184 double PressureNew = mnc * exp(x[0]);
185
186 // UPDATE: Currently not used
187 //double current_precision = Broyden::TOL;
189 //while (abs(PressureNew) < current_precision && current_precision > 1.e-50 && abs(PressureNew /m_Pressure) < 1.e-5) {
190 // current_precision *= 1.e-10;
191 // x = broydn.Solve(x, &Broyden::BroydenSolutionCriterium(current_precision));
192 // PressureNew = mnc * exp(x[0]);
193 //}
194
195 m_Pressure = PressureNew;
196
197 if (broydn.Iterations() == broydn.MaxIterations())
198 m_LastCalculationSuccessFlag = false;
199 else m_LastCalculationSuccessFlag = true;
200
201 m_MaxDiff = broydn.MaxDifference();
202 }
203
205 m_FluctuationsCalculated = false;
206
208
209 m_wnSum = 0.;
210 m_Densityid = 0.;
211 m_TotalDensity = 0.;
212 m_Suppression = 0.;
213 double densityid = 0., suppression = 0.;
214
216
217#pragma omp parallel for reduction(+:densityid) reduction(+:suppression) if(useOpenMP)
218 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
219 double dMu = -m_v[i] * m_Pressure;
220 m_densitiesid[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i] + dMu);
221 densityid += m_densitiesid[i];
222 suppression += m_v[i] * m_densitiesid[i];
223
224 m_densitiesidnoshift[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
225 }
226
227 m_Densityid = densityid;
228 m_Suppression = suppression;
229
230 m_Suppression = 1. / (1. + m_Suppression);
231
232 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
233 m_densities[i] = m_densitiesid[i] * m_Suppression;
234 m_TotalDensity += m_densities[i];
235 m_wnSum += m_densities[i] * m_TPS->Particles()[i].ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure);
236 }
237
239
240 m_Calculated = true;
242 }
243
245 int NN = m_densities.size();
246 vector<double> tN(NN);
247 for (int i = 0; i < NN; ++i) tN[i] = DensityId(i, m_Pressure);
248
249 vector<double> chi2id(NN);
250 for (int i = 0; i < NN; ++i) {
251 chi2id[i] = m_TPS->Particles()[i].chiDimensionfull(2, m_Parameters, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure) * xMath::GeVtoifm3();
252 if (tN[i] > 0.0)
253 chi2id[i] *= m_densities[i] / tN[i];
254 }
255
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;
260
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];
267 double tmp = 0.;
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;
271
272 m_dmusdmu[i][j] = (i == j ? 1. : 0.) - m_v[i] * m_densities[j];
273 }
274 }
275
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.;
280 }
281
282 }
283
285 if (!IsTemperatureDerivativesCalculated())
287
288 // Compute dE/dT
289 double ret = 0.;
290
291 double eid = 0., dTbn = 0., bn = 0., dTeidterm = 0., dmueidterm = 0.;
292 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
293 const ThermalParticle &part = m_TPS->Particles()[i];
294 eid += part.Density(m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure);
295 dTbn += m_v[i] * m_dndT[i];
296 bn += m_v[i] * m_densities[i];
297 dTeidterm += part.Density(m_Parameters, IdealGasFunctions::dedT, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure);
298 dmueidterm += -m_v[i] * EntropyDensity() * part.Density(m_Parameters, IdealGasFunctions::dedmu, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure);
299 }
300
301 ret = -dTbn * eid + (1. - bn) * (dTeidterm + dmueidterm);
302
303 return ret;
304 }
305
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.));
311
312 double T = m_Parameters.T;
313
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) {
316 nids[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure);
317 dniddTs[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::dndT, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure);
318 chi2ids[i] = m_TPS->Particles()[i].chiDimensionfull(2, m_Parameters, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure) * xMath::GeVtoifm3();
319 dchi2idsdT[i] = (T * T * m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::dchi2dT, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure) * xMath::GeVtoifm3() + 2. * T * chi2ids[i] / T / T);
320 chi3ids[i] = m_TPS->Particles()[i].chiDimensionfull(3, m_Parameters, m_UseWidth, m_Chem[i] - m_v[i] * m_Pressure) * xMath::GeVtoifm3();
321 }
322
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];
333 }
334
335 double sum_dbndT = (1. - sum_bn) * (1. - sum_bn) * (sum_dTbns - sum_chi2idb2 * s);
336
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]);
339 }
340
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.;
345
346 double tval = 0.;
347 tval += -chi2ids[j] * m_densities[i] * m_v[j];
348 tval += sum_chi2idb2 * m_densities[j] * m_densities[i];
349 if (i == j)
350 tval += chi2ids[i];
351 tval += -chi2ids[i] * m_densities[j] * m_v[i];
352 m_PrimChi2sdT[i][j] += -sum_dbndT * tval;
353
354 tval = 0.;
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]);
361 if (i == j) {
362 tval += dchi2idsdT[i];
363 tval += chi3ids[i] * m_dmusdT[i];
364 }
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];
368
369 m_PrimChi2sdT[i][j] += (1. - sum_bn) * tval;
370
371 m_PrimChi2sdT[i][j] = m_PrimChi2sdT[i][j] / (xMath::GeVtoifm3() * T * T) - 2. * m_PrimCorrel[i][j] / T / T / T / xMath::GeVtoifm3();
372 }
373 }
374
375 m_TemperatureDerivativesCalculated = true;
377 }
378
379 m_TemperatureDerivativesCalculated = true;
380 }
381
382 // TODO include correlations
389
390 m_FluctuationsCalculated = true;
391
392 for (size_t i = 0; i < m_wprim.size(); ++i) {
393 m_skewprim[i] = ParticleSkewness(i);
394 m_kurtprim[i] = ParticleKurtosis(i);
395 }
396
397 for (size_t i = 0; i < m_wprim.size(); ++i) {
398 m_skewtot[i] = 1.;
399 m_kurttot[i] = 1.;
400 }
401
402 if (IsTemperatureDerivativesCalculated()) {
403 m_TemperatureDerivativesCalculated = false;
405 }
406 }
407
408
409 std::vector<double> ThermalModelEVDiagonal::CalculateChargeFluctuations(const std::vector<double>& chgs, int order)
410 {
411 vector<double> ret(order + 1, 0.);
412
413 // chi1
414 for (size_t i = 0; i < m_densities.size(); ++i)
415 ret[0] += chgs[i] * m_densities[i];
416
417 ret[0] /= pow(m_Parameters.T * xMath::GeVtoifm(), 3);
418
419 if (order < 2) return ret;
420
421
422
423 // Preparing matrix for system of linear equations
424 int NN = m_densities.size();
425
426 vector<double> MuStar(NN, 0.);
427 for (int i = 0; i < NN; ++i) {
428 MuStar[i] = m_Chem[i] + MuShift(i);
429 }
430
431
432 MatrixXd densMatrix(2 * NN, 2 * NN);
433 VectorXd solVector(2 * NN), xVector(2 * NN);
434
435 vector<double> DensitiesId(m_densities.size()), chi2id(m_densities.size());
436 for (int i = 0; i < NN; ++i) {
437 DensitiesId[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, MuStar[i]);
438 chi2id[i] = m_TPS->Particles()[i].chi(2, m_Parameters, m_UseWidth, MuStar[i]);
439 }
440
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.;
445 }
446
447 for (int i = 0; i < NN; ++i)
448 for (int j = 0; j < NN; ++j)
449 densMatrix(i, NN + j) = 0.;
450
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];
455 }
456 densMatrix(i, NN + i) = (densMatrix(i, NN + i) - 1.) * chi2id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * m_Parameters.T;
457 }
458
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.;
464 }
465
466
467 PartialPivLU<MatrixXd> decomp(densMatrix);
468
469 // chi2
470 vector<double> dni(NN, 0.), dmus(NN, 0.);
471
472 for (int i = 0; i < NN; ++i) {
473 xVector[i] = 0.;
474 xVector[NN + i] = chgs[i];
475 }
476
477 solVector = decomp.solve(xVector);
478
479 for (int i = 0; i < NN; ++i) {
480 dni[i] = solVector[i];
481 dmus[i] = solVector[NN + i];
482 }
483
484 for (int i = 0; i < NN; ++i)
485 ret[1] += chgs[i] * dni[i];
486
487 ret[1] /= pow(m_Parameters.T, 2) * pow(xMath::GeVtoifm(), 3);
488
489 if (order < 3) return ret;
490 // chi3
491 vector<double> d2ni(NN, 0.), d2mus(NN, 0.);
492
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]);
496
497 for (int i = 0; i < NN; ++i) {
498 xVector[i] = 0.;
499
500 double tmp = 0.;
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];
503 xVector[i] += tmp;
504
505 tmp = 0.;
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];
508 xVector[i] += tmp;
509 }
510 for (int i = 0; i < NN; ++i) {
511 xVector[NN + i] = 0.;
512
513 double tmp = 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];
515
516 xVector[NN + i] = tmp;
517 }
518
519 solVector = decomp.solve(xVector);
520
521 for (int i = 0; i < NN; ++i) {
522 d2ni[i] = solVector[i];
523 d2mus[i] = solVector[NN + i];
524 }
525
526 for (int i = 0; i < NN; ++i)
527 ret[2] += chgs[i] * d2ni[i];
528
529 ret[2] /= m_Parameters.T * pow(xMath::GeVtoifm(), 3);
530
531
532 if (order < 4) return ret;
533
534 // chi4
535 vector<double> d3ni(NN, 0.), d3mus(NN, 0.);
536
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]);
540
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];
544 }
545
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];
550 }
551
552 for (int i = 0; i < NN; ++i) {
553 xVector[i] = 0.;
554
555 double tmp = 0.;
556 for (int j = 0; j < NN; ++j) tmp += m_v[j] * dni[j];
557 tmp = -3. * tmp * d2nis[i];
558 xVector[i] += tmp;
559
560 tmp = 0.;
561 for (int j = 0; j < NN; ++j) tmp += m_v[j] * d2ni[j];
562 tmp = -3. * tmp * dnis[i];
563 xVector[i] += tmp;
564
565 double tmps = 0.;
566 for (int j = 0; j < NN; ++j) tmps += m_v[j] * m_densities[j];
567
568 tmp = -(tmps - 1.) * chi3id[i] * pow(xMath::GeVtoifm(), 3) * m_Parameters.T * d2mus[i] * 3. * dmus[i];
569 xVector[i] += tmp;
570
571 tmp = -(tmps - 1.) * chi4id[i] * pow(xMath::GeVtoifm(), 3) * dmus[i] * dmus[i] * dmus[i];
572 xVector[i] += tmp;
573 }
574 for (int i = 0; i < NN; ++i) {
575 xVector[NN + i] = 0.;
576
577 double tmp = 0.;
578 for (int j = 0; j < NN; ++j) tmp += -2. * m_v[i] * d2mus[j] * dnis[j];
579 xVector[NN + i] += tmp;
580
581 tmp = 0.;
582 for (int j = 0; j < NN; ++j) tmp += -m_v[i] * dmus[j] * d2nis[j];
583 xVector[NN + i] += tmp;
584 }
585
586 solVector = decomp.solve(xVector);
587
588 for (int i = 0; i < NN; ++i) {
589 d3ni[i] = solVector[i];
590 d3mus[i] = solVector[NN + i];
591 }
592
593 for (int i = 0; i < NN; ++i)
594 ret[3] += chgs[i] * d3ni[i];
595
596 ret[3] /= pow(xMath::GeVtoifm(), 3);
597
598 return ret;
599 }
600
601
603 if (!m_Calculated) CalculateDensities();
604 double ret = 0.;
605 double dMu = 0.;
606 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
607 dMu = -m_v[i] * m_Pressure;
608 ret += m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i] + dMu);
609 }
610 return ret * m_Suppression;
611 }
612
614 if (!m_Calculated) CalculateDensities();
615 double ret = 0.;
616 double dMu = 0.;
617 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
618 dMu = -m_v[i] * m_Pressure;
619 ret += m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::EntropyDensity, m_UseWidth, m_Chem[i] + dMu);
620 }
621 return ret * m_Suppression;
622 }
623
625 if (!m_Calculated) CalculateDensities();
626 return m_Pressure;
627 }
628
630 if (!m_Calculated) CalculateDensities();
631
632 double dMu = -m_v[part] * m_Pressure;
633 double ret = m_TPS->Particles()[part].Density(m_Parameters, IdealGasFunctions::ScalarDensity, m_UseWidth, m_Chem[part] + dMu);
634 return ret * m_Suppression;
635 }
636
638 {
639 if (!m_Calculated)
641 return m_Suppression;
642 }
643
645 {
646 if (id >= 0. && id < static_cast<int>(m_v.size()))
647 return -m_v[id] * m_Pressure;
648 else
649 return 0.0;
650 }
651
653 double tEV = 0.;
654 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
655 tEV += m_v[i] * m_densities[i] / 4.;
656 }
657 return tEV;
658 }
659
660 void ThermalModelEVDiagonal::SetRadius(int i, double rad)
661 {
662 if (i >= 0 && i < static_cast<int>(m_v.size()))
663 m_v[i] = CuteHRGHelper::vr(rad);
664 }
665
666 double ThermalModelEVDiagonal::VirialCoefficient(int i, int /*j*/) const
667 {
668 return m_v[i];
669 }
670
671 void ThermalModelEVDiagonal::SetVirial(int i, int j, double b)
672 {
673 if (i == j)
674 m_v[i] = b;
675 }
676
677 std::vector<double> ThermalModelEVDiagonal::BroydenEquationsDEV::Equations(const std::vector<double>& x)
678 {
679 std::vector<double> ret(1);
680 double pressure = m_mnc * exp(x[0]);
681
682 ret[0] = pressure - m_THM->Pressure(pressure);
683
684 return ret;
685 }
686
687 std::vector<double> ThermalModelEVDiagonal::BroydenJacobianDEV::Jacobian(const std::vector<double>& x)
688 {
689 double pressure = m_mnc * exp(x[0]);
690
691 double ret = 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);
694 ret += 1.;
695 ret *= pressure;
696
697 return std::vector<double>(1, ret);
698 }
699
700 std::vector<double> ThermalModelEVDiagonal::BroydenEquationsDEVOrig::Equations(const std::vector<double>& x)
701 {
702 std::vector<double> ret(1);
703 ret[0] = x[0] - m_THM->Pressure(x[0]);
704 return ret;
705 }
706
707 std::vector<double> ThermalModelEVDiagonal::BroydenJacobianDEVOrig::Jacobian(const std::vector<double>& x)
708 {
709 const double &pressure = x[0];
710
711 double ret = 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);
714 ret += 1.;
715
716 return std::vector<double>(1, ret);
717 }
718
719 bool ThermalModelEVDiagonal::BroydenSolutionCriteriumDEV::IsSolved(const std::vector<double>& x, const std::vector<double>& f, const std::vector<double>& /*xdelta*/) const
720 {
721 double maxdiff = 0.;
722 for (size_t i = 0; i < x.size(); ++i) {
723 maxdiff = std::max(maxdiff, fabs(f[i]) / m_THM->m_Pressure);
724 }
725 return (maxdiff < m_MaximumError);
726 }
727} // namespace thermalfist
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.
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
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 &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...
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.
void FillVirial(const std::vector< double > &ri=std::vector< double >(0))
Fills the vector of particle eigenvolume parameters.
ThermalModelEVDiagonal(ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
Construct a new ThermalModelEVDiagonal object.
virtual void CalculateFluctuations()
Computes the fluctuation observables.
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...
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.
virtual double MuShift(int i) const
The shift in the chemical potential of particle species i due to the excluded volume interactions.
virtual void SolvePressure()
Solves the transcdental equation for the pressure.
virtual void CalculatePrimordialDensities()
Calculates the primordial densities of all species.
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 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 &params, 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 .
Definition xMath.h:25
constexpr double GeVtoifm3()
A constant to transform GeV into fm .
Definition xMath.h:31
constexpr double mnucleon()
Nucleon's mass. Value as in UrQMD.
Definition xMath.h:34
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.