Thermal-FIST 1.6
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
ThermalModelCanonical.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#include <iostream>
11#include <cmath>
12#include <cstdlib>
13#include <algorithm>
14#include <cassert>
15
16#include "HRGBase/xMath.h"
18
19using namespace std;
20
21namespace thermalfist {
22
25 m_BCE(1),
26 m_QCE(1),
27 m_SCE(1),
28 m_CCE(1),
30 {
31
32 m_TAG = "ThermalModelCanonical";
33
34 m_Ensemble = CE;
35 m_InteractionModel = Ideal;
36
37 m_modelgce = NULL;
38
39 m_Banalyt = false;
42 }
43
44
45
47 {
48 CleanModelGCE();
49 }
50
54
56 {
57 m_BMAX = 0;
58 m_QMAX = 0;
59 m_SMAX = 0;
60 m_CMAX = 0;
61
62
63 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
64 ThermalParticle &part = m_TPS->Particle(i);
65
66 if (part.Statistics() != 0 && IsParticleCanonical(part)) {
68
69 m_BMAX = max(m_BMAX, abs(part.BaryonCharge() * part.ClusterExpansionOrder()));
70 m_QMAX = max(m_QMAX, abs(part.ElectricCharge() * part.ClusterExpansionOrder()));
71 m_SMAX = max(m_SMAX, abs(part.Strangeness() * part.ClusterExpansionOrder()));
72 m_CMAX = max(m_CMAX, abs(part.Charm() * part.ClusterExpansionOrder()));
73 }
74 else {
75 m_BMAX = max(m_BMAX, abs(part.BaryonCharge()));
76 m_QMAX = max(m_QMAX, abs(part.ElectricCharge()));
77 m_SMAX = max(m_SMAX, abs(part.Strangeness()));
78 m_CMAX = max(m_CMAX, abs(part.Charm()));
79 }
80 }
81
86
87 if (computeFluctuations) {
88 m_BMAX *= 2;
89 m_QMAX *= 2;
90 m_SMAX *= 2;
91 m_CMAX *= 2;
92 }
93
94 // Some charges may be treated grand-canonically
95 m_BMAX *= m_BCE;
96 m_QMAX *= m_QCE;
97 m_SMAX *= m_SCE;
98 m_CMAX *= m_CCE;
99
100 std::cout << "BMAX = " << m_BMAX << "\tQMAX = " << m_QMAX << "\tSMAX = " << m_SMAX << "\tCMAX = " << m_CMAX << std::endl;
101
102 m_QNMap.clear();
103 m_QNvec.resize(0);
104
105 m_Corr.resize(0);
106 m_PartialZ.resize(0);
107
108 int ind = 0;
109 for (int iB = -m_BMAX; iB <= m_BMAX; ++iB)
110 for (int iQ = -m_QMAX; iQ <= m_QMAX; ++iQ)
111 for (int iS = -m_SMAX; iS <= m_SMAX; ++iS)
112 for (int iC = -m_CMAX; iC <= m_CMAX; ++iC) {
113
114 QuantumNumbers qn(iB, iQ, iS, iC);
115 m_QNMap[qn] = ind;
116 m_QNvec.push_back(qn);
117
118 m_PartialZ.push_back(0.);
119 m_Corr.push_back(1.);
120
121 ind++;
122 }
123
124 m_PartialZCalculatedFlucts = computeFluctuations;
125 }
126
128 m_QuantumStats = stats;
129 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
130 m_TPS->Particle(i).UseStatistics(stats);
131
132 // Only cluster expansion method supported for particles with canonically conserved charges
133 if (stats && IsParticleCanonical(m_TPS->Particle(i)))
134 m_TPS->Particle(i).SetCalculationType(IdealGasFunctions::ClusterExpansion);
135 }
136 m_PartialZ.clear();
137 //CalculateQuantumNumbersRange();
138 }
139
141 {
142 m_ConstrainMuB &= !m_BCE;
143 m_ConstrainMuQ &= !m_QCE;
144 m_ConstrainMuS &= !m_SCE;
145 m_ConstrainMuC &= !m_CCE;
147 }
148
150 {
151 m_ConstrainMuB &= !m_BCE;
152 m_ConstrainMuQ &= !m_QCE;
153 m_ConstrainMuS &= !m_SCE;
154 m_ConstrainMuC &= !m_CCE;
156 }
157
158
160 assert(m_IGFExtraConfig.MagneticField.B == 0.); // No magnetic field supported currently
161
162 m_FluctuationsCalculated = false;
163
164 if (m_PartialZ.size() == 0)
166
167 if (m_BMAX_list == 1 && m_BCE && m_QCE && m_SCE && m_CCE && !UsePartialChemicalEquilibrium()) {
168 m_Banalyt = true;
169 m_Parameters.muB = 0.0;
170 m_Parameters.muQ = 0.0;
171 m_Parameters.muS = 0.0;
172 m_Parameters.muC = 0.0;
173 }
174 else {
175 m_Banalyt = false;
176 if (m_BCE)
177 m_Parameters.muB = 0.0;
178 if (m_QCE)
179 m_Parameters.muQ = 0.0;
180 if (m_SCE)
181 m_Parameters.muS = 0.0;
182 if (m_CCE)
183 m_Parameters.muC = 0.0;
184
185 //PrepareModelGCE(); // Plan B, may work better when quantum numbers are large
186 }
187
189
190 for (size_t i = 0; i < m_densities.size(); ++i) {
191 ThermalParticle &tpart = m_TPS->Particle(i);
192 m_densities[i] = 0.;
193
194 if (!IsParticleCanonical(tpart)) {
195 m_densities[i] = tpart.Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
196 }
197 else if (tpart.Statistics() == 0
199 {
200 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
201
202 if (ind < static_cast<int>(m_Corr.size()))
203 m_densities[i] = m_Corr[ind] * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
204 else {
205 cout << "ThermalModelCanonical::CalculatePrimordialDensities: Warning! No canonical partition function for this particle!" << endl;
206 cout << "B = " << m_BCE * tpart.BaryonCharge() << "\tQ = " << m_QCE * tpart.ElectricCharge() << "\tS = " << m_SCE * tpart.Strangeness() << "\tC = " << m_CCE * tpart.Charm() << endl;
207 }
208 }
209 else {
210 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
211 int ind = m_QNMap[QuantumNumbers(m_BCE*n*tpart.BaryonCharge(), m_QCE*n*tpart.ElectricCharge(), m_SCE*n*tpart.Strangeness(), m_CCE*n*tpart.Charm())];
212 if (ind < static_cast<int>(m_Corr.size()))
213 m_densities[i] += m_Corr[ind] * tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
214 }
215 }
216 }
217
218 m_Calculated = true;
220 }
221
223 {
225
226 char cc[1000];
227
228 double TOL = 1.e-4;
229
230 // Checking that the CE calculation is valid
231 if (m_BCE && m_BMAX != 0) {
232 double totB = CalculateBaryonDensity() * m_Parameters.SVc;
233 if (fabs(m_Parameters.B - totB) > TOL) {
234
235 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total baryon number.\
236\n\
237Expected: %d\n\
238Obtained: %lf\n\
239\n", m_Parameters.B, totB);
240
241 printf("%s", cc);
242
243 m_ValidityLog.append(cc);
244
245 m_LastCalculationSuccessFlag = false;
246 }
247 }
248
249
250 if (m_QCE && m_QMAX != 0) {
251 double totQ = CalculateChargeDensity() * m_Parameters.SVc;
252 if (fabs(m_Parameters.Q - totQ) > TOL) {
253 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total electric charge.\
254\n\
255Expected: %d\n\
256Obtained: %lf\n\
257\n", m_Parameters.Q, totQ);
258
259 printf("%s", cc);
260
261 m_ValidityLog.append(cc);
262
263 m_LastCalculationSuccessFlag = false;
264 }
265 }
266
267
268 if (m_SCE && m_SMAX != 0) {
269 double totS = CalculateStrangenessDensity() * m_Parameters.SVc;
270 if (fabs(m_Parameters.S - totS) > TOL) {
271 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total strangeness.\
272\n\
273Expected: %d\n\
274Obtained: %lf\n\
275\n", m_Parameters.S, totS);
276
277 printf("%s", cc);
278
279 m_ValidityLog.append(cc);
280
281 m_LastCalculationSuccessFlag = false;
282 }
283 }
284
285
286 if (m_CCE && m_CMAX != 0) {
287 double totC = CalculateCharmDensity() * m_Parameters.SVc;
288 if (fabs(m_Parameters.C - totC) > TOL) {
289 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total charm.\
290\n\
291Expected: %d\n\
292Obtained: %lf\n\
293\n", m_Parameters.C, totC);
294
295 printf("%s", cc);
296
297 m_ValidityLog.append(cc);
298
299 m_LastCalculationSuccessFlag = false;
300 }
301 }
302 }
303
305 {
306 if (Vc < 0.0)
307 Vc = m_Parameters.SVc;
308
311 else {
312 // Partial chemical equilibrium canonical ensemble currently works only if there is particle-antiparticle symmetry
313 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
314 int i2 = m_TPS->PdgToId(-m_TPS->Particle(i).PdgId());
315 if (i2 != -1) {
316 if (std::abs(m_Chem[i] - m_Chem[i2]) > 1.e-8) {
317 throw std::runtime_error("ThermalModelCanonical::CalculatePartitionFunctions: Partial chemical equilibrium canonical ensemble only supported if particle-antiparticle fugacities are symmetric!");
318 }
319 }
320 }
321 }
322
323 bool AllMuZero = true;
324 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
325 ThermalParticle &tpart = m_TPS->Particle(i);
326 if (IsParticleCanonical(tpart) && m_Chem[i] != 0.0)
327 {
328 AllMuZero = false;
329 break;
330 }
331 }
332
333 vector<double> Nsx(m_PartialZ.size(), 0.);
334 vector<double> Nsy(m_PartialZ.size(), 0.);
335
336 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
337 ThermalParticle &tpart = m_TPS->Particle(i);
338
339 if (!IsParticleCanonical(tpart)) {
340 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
341 if (ind != m_QNMap[QuantumNumbers(0, 0, 0, 0)]) {
342 throw std::invalid_argument("ThermalModelCanonical: neutral particle cannot have non-zero ce charges");
343 }
344 if (ind < static_cast<int>(Nsx.size()))
345 Nsx[ind] += tpart.Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
346 }
347 else if (tpart.Statistics() == 0
349 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
350
351 if (ind < static_cast<int>(Nsx.size())) {
353 double tdens = tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, 0.);
354 Nsx[ind] += tdens * cosh(m_Chem[i] / m_Parameters.T);
355 Nsy[ind] += tdens * sinh(m_Chem[i] / m_Parameters.T);
356 }
357 // Currently only works at mu = 0!!
358 else {
359 double tdens = tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
360 Nsx[ind] += tdens;
361 }
362 }
363 }
364 else {
365 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
366 int ind = m_QNMap[QuantumNumbers(m_BCE * n * tpart.BaryonCharge(), m_QCE * n * tpart.ElectricCharge(), m_SCE * n * tpart.Strangeness(), m_CCE * n * tpart.Charm())];
367 if (ind < static_cast<int>(Nsx.size())) {
369 double tdens = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, 0.) / static_cast<double>(n); // TODO: Check
370 Nsx[ind] += tdens * cosh(n * m_Chem[i] / m_Parameters.T);
371 Nsy[ind] += tdens * sinh(n * m_Chem[i] / m_Parameters.T);
372 }
373 // Currently only works at mu = 0!!
374 else {
375 double tdens = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]) / static_cast<double>(n); // TODO: Check
376 Nsx[ind] += tdens;
377 }
378 }
379 }
380 }
381 }
382
383 for (int i = 0; i < static_cast<int>(Nsx.size()); ++i) {
384 Nsx[i] *= Vc;
385 Nsy[i] *= Vc;
386 m_PartialZ[i] = 0.;
387 }
388
389 int nmax = max(3, (int)sqrt(m_Parameters.B*m_Parameters.B + m_Parameters.Q*m_Parameters.Q + m_Parameters.S*m_Parameters.S + m_Parameters.C*m_Parameters.C));
390 if (m_Parameters.B == 0 && m_Parameters.Q == 0 && m_Parameters.S == 0 && m_Parameters.C == 0)
391 nmax = 4;
392
393
394
395 int nmaxB = max(4, m_Parameters.B);
396 int nmaxQ = max(4, m_Parameters.Q);
397 int nmaxS = max(4, m_Parameters.S);
398 int nmaxC = max(4, m_Parameters.C);
399
400 // UPDATE: allow to increase the number of interations externally
402 nmaxB = nmaxQ = nmaxS = nmaxC = nmax;
403
404
405
406
407
408 m_MultExp = 0.;
409 m_MultExpBanalyt = 0.;
410 for (size_t i = 0; i < m_PartialZ.size(); ++i) {
411 if (!m_Banalyt || m_QNvec[i].B == 0)
412 m_MultExp += Nsx[i];
413 if (m_Banalyt && (m_QNvec[i].B == 1 || m_QNvec[i].B == -1))
414 m_MultExpBanalyt += Nsx[i];
415 }
416
417 double dphiB = xMath::Pi() / nmaxB;
418 int maxB = 2 * nmaxB;
419 if (m_BMAX == 0 || m_Banalyt)
420 maxB = 1;
421
422 for (int iB = 0; iB < maxB; ++iB) {
423
424 vector<double> xlegB, wlegB;
425
426 if (m_BMAX != 0 && !m_Banalyt) {
427 double aB = iB * dphiB;
428 if (iB >= nmaxB) aB = xMath::Pi() - (iB + 1) * dphiB;
429 double bB = aB + dphiB;
431 }
432 else {
433 xlegB.resize(1);
434 xlegB[0] = 0.;
435 wlegB.resize(1);
436 wlegB[0] = 1.;
437 }
438
439
440 double dphiS = xMath::Pi() / nmaxS;
441 int maxS = 2 * nmaxS;
442 if (m_SMAX == 0)
443 maxS = 1;
444
445 for (int iS = 0; iS < maxS; ++iS) {
446 vector<double> xlegS, wlegS;
447
448 if (m_SMAX != 0) {
449 double aS = iS * dphiS;
450 if (iS >= nmaxS) aS = xMath::Pi() - (iS + 1) * dphiS;
451 double bS = aS + dphiS;
453 }
454 else {
455 xlegS.resize(1);
456 xlegS[0] = 0.;
457 wlegS.resize(1);
458 wlegS[0] = 1.;
459 }
460
461 double dphiQ = xMath::Pi() / nmaxQ;
462 int maxQ = nmaxQ;
463 if (m_QMAX == 0)
464 maxQ = 1;
465
466 for (int iQ = 0; iQ < maxQ; ++iQ) {
467 vector<double> xlegQ, wlegQ;
468
469 if (m_QMAX != 0) {
470 double aQ = iQ * dphiQ;
471 double bQ = aQ + dphiQ;
473 }
474 else {
475 xlegQ.resize(1);
476 xlegQ[0] = 0.;
477 wlegQ.resize(1);
478 wlegQ[0] = 1.;
479 }
480
481 double dphiC = xMath::Pi() / nmaxC;
482 int maxC = 2 * nmaxC;
483 if (m_CMAX == 0)
484 maxC = 1;
485
486 for (int iC = 0; iC < maxC; ++iC) {
487 vector<double> xlegC, wlegC;
488 if (m_CMAX != 0) {
489 double aC = iC * dphiC;
490 if (iC >= nmaxC) aC = xMath::Pi() - (iC + 1) * dphiC;
491 double bC = aC + dphiC;
493 }
494 else {
495 xlegC.resize(1);
496 xlegC[0] = 0.;
497 wlegC.resize(1);
498 wlegC[0] = 1.;
499 }
500
501 for (size_t iBt = 0; iBt < xlegB.size(); ++iBt) {
502 for (size_t iSt = 0; iSt < xlegS.size(); ++iSt) {
503 for (size_t iQt = 0; iQt < xlegQ.size(); ++iQt) {
504 for (size_t iCt = 0; iCt < xlegC.size(); ++iCt) {
505 vector<double> cosph(m_PartialZ.size(), 0.), sinph(m_PartialZ.size(), 0.);
506 double wx = 0., wy = 0., mx = 0., my = 0.;
507 for (size_t i = 0; i < m_PartialZ.size(); ++i) {
508 int tB = m_QNvec[i].B;
509 int tQ = m_QNvec[i].Q;
510 int tS = m_QNvec[i].S;
511 int tC = m_QNvec[i].C;
512
513 if (m_Banalyt) {
514 cosph[i] = cos(tS*xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
515 sinph[i] = sin(tS*xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
516
517 if (m_QNvec[i].B == 1) {
518 wx += Nsx[i] * cosph[i];
519 wy += Nsx[i] * sinph[i];
520 }
521 else if (m_QNvec[i].B == 0) {
522 mx += Nsx[i] * (cosph[i] - 1.);
523 }
524 }
525 else {
526 cosph[i] = cos(tB*xlegB[iBt] + tS * xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
527 mx += Nsx[i] * (cosph[i] - 1.);
528
529 if (!AllMuZero) {
530 sinph[i] = sin(tB*xlegB[iBt] + tS * xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
531 my += Nsy[i] * sinph[i];
532 }
533 }
534 }
535
536 double wmod = 0.;
537 double warg = 0.;
538
539 if (m_Banalyt) {
540 wmod = sqrt(wx*wx + wy * wy);
541 warg = atan2(wy, wx);
542 }
543
544 for (size_t iN = 0; iN < m_PartialZ.size(); ++iN) {
545 int tBg = m_Parameters.B - m_QNvec[iN].B;
546 int tQg = m_Parameters.Q - m_QNvec[iN].Q;
547 int tSg = m_Parameters.S - m_QNvec[iN].S;
548 int tCg = m_Parameters.C - m_QNvec[iN].C;
549
550 if (m_Banalyt) {
551 //m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] * cos(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt] - tBg * warg) * exp(mx) * xMath::BesselI(tBg, 2. * wmod);
552 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] *
553 cos(tBg * xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt] - tBg * warg) *
554 exp(mx + 2. * wmod - m_MultExpBanalyt) *
555 xMath::BesselIexp(tBg, 2. * wmod);
556 }
557 else {
558 if (AllMuZero)
559 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] * cos(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * exp(mx);
560 else
561 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] * exp(mx)
562 * (cos(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * cos(my)
563 + sin(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * sin(my));
564 }
565 }
566 }
567 }
568 }
569 }
570
571 //int cind = iB * maxS * maxQ * maxC + iS * maxQ * maxC + iQ * maxC + iC;
572 //int tot = maxB * maxS * maxQ * maxC;
573 }
574 }
575 }
576 }
577
578 for (size_t iN = 0; iN < m_PartialZ.size(); ++iN) {
579 if (m_BMAX != 0 && m_BMAX != 1 && !m_Banalyt) // TODO: cross-check
580 m_PartialZ[iN] /= 2. * xMath::Pi();
581 if (m_QMAX != 0) {
582 m_PartialZ[iN] /= 2. * xMath::Pi();
583 m_PartialZ[iN] *= 2.;
584 }
585 if (m_SMAX != 0)
586 m_PartialZ[iN] /= 2. * xMath::Pi();
587 if (m_CMAX != 0)
588 m_PartialZ[iN] /= 2. * xMath::Pi();
589 }
590
591
592 m_Corr.resize(m_PartialZ.size());
593 for (size_t iN = 0; iN < m_PartialZ.size(); ++iN) {
594 m_Corr[iN] = m_PartialZ[iN] / m_PartialZ[m_QNMap[QuantumNumbers(0, 0, 0, 0)]];
595 }
596
598 }
599
601 {
602 ThermalParticle &tpart = m_TPS->Particle(part);
603 double ret1 = 0., ret2 = 0., ret3 = 0.;
604
605 if (!IsParticleCanonical(tpart)) {
606 return tpart.ScaledVariance(m_Parameters, m_UseWidth, m_Chem[part]);
607 }
608 else if (tpart.Statistics() == 0
610 {
611 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
612 int ind2 = m_QNMap[QuantumNumbers(m_BCE * 2 * tpart.BaryonCharge(), m_QCE * 2 * tpart.ElectricCharge(), m_SCE * 2 * tpart.Strangeness(), m_CCE * 2 * tpart.Charm())];
613
614 ret1 = 1.;
615 if (ind < static_cast<int>(m_Corr.size()) && ind2 < static_cast<int>(m_Corr.size()))
616 ret2 = m_Corr[ind2] / m_Corr[ind] * m_Parameters.SVc * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
617
618 if (ind < static_cast<int>(m_Corr.size()))
619 ret3 = -m_Corr[ind] * m_Parameters.SVc * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
620 }
621 else {
622 double ret1num = 0., ret1zn = 0.;
623 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
624 int ind = m_QNMap[QuantumNumbers(m_BCE*n*tpart.BaryonCharge(), m_QCE*n*tpart.ElectricCharge(), m_SCE*n*tpart.Strangeness(), m_CCE*n*tpart.Charm())];
625
626 double densityClusterN = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
627
628 if (ind < static_cast<int>(m_Corr.size())) {
629 ret1num += m_Corr[ind] * n * densityClusterN;
630 ret1zn += m_Corr[ind] * densityClusterN;
631 }
632
633 for (int n2 = 1; n2 <= tpart.ClusterExpansionOrder(); ++n2) {
634 if (m_QNMap.count(QuantumNumbers(m_BCE * (n + n2)*tpart.BaryonCharge(), m_QCE * (n + n2)*tpart.ElectricCharge(), m_SCE * (n + n2)*tpart.Strangeness(), m_CCE * (n + n2)*tpart.Charm())) != 0) {
635 int ind2 = m_QNMap[QuantumNumbers(m_BCE * (n + n2)*tpart.BaryonCharge(), m_QCE * (n + n2)*tpart.ElectricCharge(), m_SCE * (n + n2)*tpart.Strangeness(), m_CCE * (n + n2)*tpart.Charm())];
636 if (ind < static_cast<int>(m_Corr.size()) && ind2 < static_cast<int>(m_Corr.size()))
637 ret2 += densityClusterN * m_Corr[ind2] * m_Parameters.SVc * tpart.DensityCluster(n2, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
638 }
639 }
640 }
641
642 if (ret1zn == 0.0)
643 return 1.;
644
645 ret1 = ret1num / ret1zn;
646 ret2 = ret2 / ret1zn;
647 ret3 = -ret1zn * m_Parameters.SVc;
648 }
649 return ret1 + ret2 + ret3;
650 }
651
653 int NN = m_densities.size();
654
655 vector<double> yld(NN, 0);
656 vector<double> ret1num(NN, 0);
657 vector< vector<double> > ret2num(NN, vector<double>(NN, 0.));
658
659 for (int i = 0; i < NN; ++i)
660 yld[i] = m_densities[i] * m_Parameters.SVc;
661
662 for (int i = 0; i < NN; ++i) {
663 ThermalParticle &tpart = m_TPS->Particle(i);
664 if (!IsParticleCanonical(tpart)) {
665 ret1num[i] = tpart.ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i]) * yld[i];
666 }
667 else if (tpart.Statistics() == 0
669 {
670 ret1num[i] = yld[i];
671 }
672 else {
673 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
674 int ind = m_QNMap[QuantumNumbers(m_BCE * n * tpart.BaryonCharge(), m_QCE * n * tpart.ElectricCharge(), m_SCE * n * tpart.Strangeness(), m_CCE * n * tpart.Charm())];
675
676 double densityClusterN = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
677
678 if (ind < static_cast<int>(m_Corr.size()))
679 ret1num[i] += m_Corr[ind] * n * densityClusterN * m_Parameters.SVc;
680 }
681 }
682 }
683
684 for (int i = 0; i < NN; ++i) {
685 for (int j = 0; j < NN; ++j) {
686 ThermalParticle &tpart1 = m_TPS->Particle(i);
687 ThermalParticle &tpart2 = m_TPS->Particle(j);
688
689 int n1max = tpart1.ClusterExpansionOrder();
690 int n2max = tpart2.ClusterExpansionOrder();
691
692 if (tpart1.Statistics() == 0 || tpart1.CalculationType() != IdealGasFunctions::ClusterExpansion)
693 n1max = 1;
694 if (tpart2.Statistics() == 0 || tpart2.CalculationType() != IdealGasFunctions::ClusterExpansion)
695 n2max = 1;
696
697 if (!IsParticleCanonical(tpart1) || !IsParticleCanonical(tpart2)) {
698 ret2num[i][j] = yld[i] * yld[j];
699 }
700 else {
701 for (int n1 = 1; n1 <= n1max; ++n1) {
702 for (int n2 = 1; n2 <= n2max; ++n2) {
703 int ind = m_QNMap[QuantumNumbers(
704 m_BCE*(n1*tpart1.BaryonCharge() + n2 * tpart2.BaryonCharge()),
705 m_QCE*(n1*tpart1.ElectricCharge() + n2 * tpart2.ElectricCharge()),
706 m_SCE*(n1*tpart1.Strangeness() + n2 * tpart2.Strangeness()),
707 m_CCE*(n1*tpart1.Charm() + n2 * tpart2.Charm()))];
708
709 double densityClusterN1 = tpart1.DensityCluster(n1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
710 double densityClusterN2 = tpart2.DensityCluster(n2, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[j]);
711
712 if (ind < static_cast<int>(m_Corr.size()))
713 ret2num[i][j] += m_Corr[ind] * densityClusterN1 * densityClusterN2 * m_Parameters.SVc * m_Parameters.SVc;
714 }
715 }
716 }
717 }
718 }
719
720
721 m_PrimCorrel.resize(NN);
722 for (int i = 0; i < NN; ++i)
723 m_PrimCorrel[i].resize(NN);
724 m_TotalCorrel = m_PrimCorrel;
725
726 for (int i = 0; i < NN; ++i) {
727 for (int j = 0; j < NN; ++j) {
728 m_PrimCorrel[i][j] = 0.;
729 if (i == j)
730 m_PrimCorrel[i][j] += ret1num[i] / yld[i];
731 m_PrimCorrel[i][j] += ret2num[i][j] / yld[i];
732 m_PrimCorrel[i][j] += -yld[j];
733 m_PrimCorrel[i][j] *= yld[i] / m_Parameters.SVc / m_Parameters.T;
734 if (yld[i] == 0.0)
735 m_PrimCorrel[i][j] = 0.;
736 }
737 }
738
739 for (int i = 0; i < NN; ++i) {
740 m_wprim[i] = m_PrimCorrel[i][i];
741 if (m_densities[i] > 0.) m_wprim[i] *= m_Parameters.T / m_densities[i];
742 else m_wprim[i] = 1.;
743 if (m_wprim[i] != m_wprim[i]) m_wprim[i] = 1.;
744 }
745
746 m_TwoParticleCorrelationsCalculated = true;
747 }
748
753 }
754
760
761 m_FluctuationsCalculated = true;
762
763 for (size_t i = 0; i < m_wprim.size(); ++i) {
764 m_skewprim[i] = 1.;
765 m_kurtprim[i] = 1.;
766 m_skewtot[i] = 1.;
767 m_kurttot[i] = 1.;
768 }
769 }
770
772 if (!m_Calculated) CalculateDensities();
773 double ret = 0.;
774
775 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
776 ThermalParticle &tpart = m_TPS->Particle(i);
777 {
778 if (!IsParticleCanonical(tpart)) {
779 ret += tpart.Density(m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i]);
780 }
781 else if (tpart.Statistics() == 0
783 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
784
785 if (ind < static_cast<int>(m_Corr.size()))
786 ret += m_Corr[ind] * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i]);
787 }
788 else {
789 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
790 int ind = m_QNMap[QuantumNumbers(m_BCE * n * tpart.BaryonCharge(), m_QCE * n * tpart.ElectricCharge(), m_SCE * n * tpart.Strangeness(), m_CCE * n * tpart.Charm())];
791 if (ind < static_cast<int>(m_Corr.size()))
792 ret += m_Corr[ind] * tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i]);
793 }
794 }
795 }
796 }
797
798 return ret;
799 }
800
802 if (!m_Calculated) CalculateDensities();
803 double ret = 0.;
804
805 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
806 ThermalParticle &tpart = m_TPS->Particle(i);
807 {
808 if (!IsParticleCanonical(tpart)) {
809 ret += tpart.Density(m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i]);
810 }
811 else if (tpart.Statistics() == 0
813 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
814
815 if (ind < static_cast<int>(m_Corr.size()))
816 ret += m_Corr[ind] * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i]);
817 }
818 else {
819 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
820 int ind = m_QNMap[QuantumNumbers(m_BCE * n * tpart.BaryonCharge(), m_QCE * n * tpart.ElectricCharge(), m_SCE * n * tpart.Strangeness(), m_CCE * n * tpart.Charm())];
821
822 if (ind < static_cast<int>(m_Corr.size()))
823 ret += m_Corr[ind] * tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i]);
824 }
825 }
826 }
827 }
828
829 return ret;
830 }
831
833 {
834 double ret = (CalculateEnergyDensity() / m_Parameters.T) + (m_MultExp + m_MultExpBanalyt + log(m_PartialZ[m_QNMap[QuantumNumbers(0, 0, 0, 0)]])) / m_Parameters.SVc;
835
836 if (m_BCE)
837 ret += -m_Parameters.muB / m_Parameters.T * m_Parameters.B / m_Parameters.SVc;
838 else
839 ret += -m_Parameters.muB * CalculateBaryonDensity() / m_Parameters.T;
840
841 if (m_QCE)
842 ret += -m_Parameters.muQ / m_Parameters.T * m_Parameters.Q / m_Parameters.SVc;
843 else
844 ret += -m_Parameters.muQ * CalculateChargeDensity() / m_Parameters.T;
845
846 if (m_SCE)
847 ret += -m_Parameters.muS / m_Parameters.T * m_Parameters.S / m_Parameters.SVc;
848 else
849 ret += -m_Parameters.muS * CalculateStrangenessDensity() / m_Parameters.T;
850
851 if (m_CCE)
852 ret += -m_Parameters.muC / m_Parameters.T * m_Parameters.C / m_Parameters.SVc;
853 else
854 ret += -m_Parameters.muC * CalculateCharmDensity() / m_Parameters.T;
855
856 return ret;
857 }
858
860 {
861 return m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
862 }
863
865 {
866 return !(
867 (part.BaryonCharge() == 0 || m_BCE == 0)
868 && (part.ElectricCharge() == 0 || m_QCE == 0)
869 && (part.Strangeness() == 0 || m_SCE == 0)
870 && (part.Charm() == 0 || m_CCE == 0)
871 );
872 }
873
875 {
876 if (charge == ConservedCharge::BaryonCharge)
877 return (m_BCE != 0);
878 else if (charge == ConservedCharge::ElectricCharge)
879 return (m_QCE != 0);
880 else if (charge == ConservedCharge::StrangenessCharge)
881 return (m_SCE != 0);
882 else if (charge == ConservedCharge::CharmCharge)
883 return (m_CCE != 0);
884 return 0;
885 }
886
887 void ThermalModelCanonical::PrepareModelGCE()
888 {
889 CleanModelGCE();
890
891 m_modelgce = new ThermalModelIdeal(m_TPS, m_Parameters);
892 m_modelgce->SetUseWidth(m_UseWidth);
894
895 if (m_BCE)
896 m_Parameters.muB = 0.0;
897
898 if (!m_BCE && m_SCE)
899 m_Parameters.muS = m_Parameters.muB / 3.;
900
901 if (!m_BCE && m_QCE)
902 m_Parameters.muQ = -m_Parameters.muB / 30.;
903
904 m_Parameters.muC = 0.;
905
906 m_modelgce->SolveChemicalPotentials(m_Parameters.B, m_Parameters.Q, m_Parameters.S, m_Parameters.C,
907 m_Parameters.muB, m_Parameters.muQ, m_Parameters.muS, m_Parameters.muC,
908 static_cast<bool>(m_BCE),
909 static_cast<bool>(m_QCE),
910 static_cast<bool>(m_SCE),
911 static_cast<bool>(m_CCE));
912
913 m_Parameters.muB = m_modelgce->Parameters().muB;
914 m_Parameters.muQ = m_modelgce->Parameters().muQ;
915 m_Parameters.muS = m_modelgce->Parameters().muS;
916 m_Parameters.muC = m_modelgce->Parameters().muC;
917
918
919 // Possible alternative below
920 //double tdens = 0.;
921 //for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
922 // ThermalParticle &part = m_TPS->Particle(i);
923 // if (part.BaryonCharge() == 1)
924 // tdens += part.Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, 0.);
925 //}
926 //m_Parameters.muB = m_Parameters.T * asinh(m_Parameters.B / m_Parameters.SVc / 2. / tdens);
927 //m_Parameters.muS = m_Parameters.muB / 3.;
928 //m_Parameters.muQ = -m_Parameters.muB / 30.;
929 }
930
931 void ThermalModelCanonical::CleanModelGCE()
932 {
933 if (m_modelgce != NULL) {
934 delete m_modelgce;
935 m_modelgce = NULL;
936 }
937 }
938
939} // namespace thermalfist
map< string, double > params
virtual void CalculateSusceptibilityMatrix()
Calculates the conserved charges susceptibility matrix.
virtual void FixParametersNoReset()
Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
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...
bool UsePartialChemicalEquilibrium()
Whether partial chemical equilibrium with additional chemical potentials is used.
virtual bool SolveChemicalPotentials(double totB=0., double totQ=0., double totS=0., double totC=0., double muBinit=0., double muQinit=0., double muSinit=0., double muCinit=0., bool ConstrMuB=true, bool ConstrMuQ=true, bool ConstrMuS=true, bool ConstrMuC=true)
The procedure which calculates the chemical potentials which reproduce the specified total baryon,...
virtual void FillChemicalPotentials()
Sets the chemical potentials of all particles.
virtual void FixParameters()
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
void SetUseWidth(bool useWidth)
Sets whether finite resonance widths are used. Deprecated.
virtual void SetChemicalPotentials(const std::vector< double > &chem=std::vector< double >(0))
Sets the chemical potentials of all particles.
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.
const ThermalModelParameters & Parameters() const
virtual void CalculateTwoParticleFluctuationsDecays()
Computes particle number correlations and fluctuations for all final-state particles which are marked...
std::map< QuantumNumbers, int > m_QNMap
Maps QuantumNumbers combinations to a 1-dimensional index.
double m_MultExpBanalyt
Exponential multiplier for analytical baryon fugacity calculations.
int m_IntegrationIterationsMultiplier
A multiplier to increase the number of iterations during the numerical integration used to calculate ...
std::vector< double > m_PartialZ
The computed canonical partition function.
int m_QCE
Flag indicating if electric charge is conserved canonically.
void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
virtual void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
int m_CCE
Flag indicating if charm is conserved canonically.
virtual void CalculateFluctuations()
Computes the fluctuation observables.
double m_MultExp
Exponential multiplier for canonical partition function calculations.
virtual double ParticleScaledVariance(int part)
std::vector< QuantumNumbers > m_QNvec
A set of QuantumNumbers combinations for which it is necessary to compute the canonical partition fun...
std::vector< double > m_Corr
A vector of chemical factors.
virtual double GetGCEDensity(int i) const
Density of particle species i in the grand-canonical ensemble.
virtual bool IsParticleCanonical(const ThermalParticle &part)
Determines whether the specified ThermalParticle is treat canonically or grand-canonically in the pre...
virtual void CalculatePartitionFunctions(double Vc=-1.)
Calculates all necessary canonical partition functions.
virtual ~ThermalModelCanonical(void)
Destroy the ThermalModelCanonical object.
virtual void CalculateQuantumNumbersRange(bool computeFluctuations=false)
Calculates the range of quantum numbers values for which it is necessary to compute the canonical par...
virtual void FixParametersNoReset()
Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
int m_BCE
Flag indicating if baryon charge is conserved canonically.
ThermalModelIdeal * m_modelgce
Pointer to a ThermalModelIdeal object used for GCE calculations.
virtual void ValidateCalculation()
Checks whether issues have occured during the calculation of particle densities in the CalculateDensi...
virtual bool IsConservedChargeCanonical(ConservedCharge::Name charge) const
bool m_Banalyt
Flag indicating whether the analytical calculation of baryon fugacity is used.
int m_SCE
Flag indicating if strangeness is conserved canonically.
virtual void FixParameters()
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
virtual void CalculatePrimordialDensities()
Calculates the primordial densities of all species.
ThermalModelCanonical(ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
Construct a new ThermalModelCanonical object.
Class implementing the Ideal HRG model.
Class containing all information about a particle specie.
int BaryonCharge() const
Particle's baryon number.
int Statistics() const
Particle's statistics.
int Strangeness() const
Particle's strangeness.
void SetCalculationType(IdealGasFunctions::QStatsCalculationType type)
Sets the CalculationType() method to evaluate quantum statistics.
double Density(const ThermalModelParameters &params, IdealGasFunctions::Quantity type=IdealGasFunctions::ParticleDensity, bool useWidth=0, double mu=0.) const
Computes a specified ideal gas thermodynamic function.
IdealGasFunctions::QStatsCalculationType CalculationType() const
Method to evaluate quantum statistics.
int ElectricCharge() const
Particle's electric charge.
int Charm() const
Particle's charm.
int ClusterExpansionOrder() const
Number of terms in the cluster expansion method.
double DensityCluster(int n, const ThermalModelParameters &params, IdealGasFunctions::Quantity type=IdealGasFunctions::ParticleDensity, bool useWidth=0, double mu=0.) const
double ScaledVariance(const ThermalModelParameters &params, bool useWidth=0, double mu=0.) const
Computes the scaled variance of particle number fluctuations in the ideal gas. Computes the scaled va...
Class containing the particle list.
void GetCoefsIntegrateLegendre10(double a, double b, std::vector< double > *x, std::vector< double > *w)
constexpr double Pi()
Pi constant.
Definition xMath.h:23
double BesselIexp(int n, double x)
Modified Bessel function I_n(x), divided by exponential factor.
Definition xMath.cpp:791
The main namespace where all classes and functions of the Thermal-FIST library reside.
Definition CosmicEoS.h:9
Name
Set of all conserved charges considered.
@ ElectricCharge
Electric charge.
Struct containing a set of quantum numbers: Baryon number, electric charge, strangeness,...
Structure containing all thermal parameters of the model.
Contains some extra mathematical functions used in the code.