Thermal-FIST 1.5
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;
40 }
41
42
43
45 {
46 CleanModelGCE();
47 }
48
52
54 {
55 m_BMAX = 0;
56 m_QMAX = 0;
57 m_SMAX = 0;
58 m_CMAX = 0;
59
60
61 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
62 ThermalParticle &part = m_TPS->Particle(i);
63
64 if (part.Statistics() != 0 && IsParticleCanonical(part)) {
66
67 m_BMAX = max(m_BMAX, abs(part.BaryonCharge() * part.ClusterExpansionOrder()));
68 m_QMAX = max(m_QMAX, abs(part.ElectricCharge() * part.ClusterExpansionOrder()));
69 m_SMAX = max(m_SMAX, abs(part.Strangeness() * part.ClusterExpansionOrder()));
70 m_CMAX = max(m_CMAX, abs(part.Charm() * part.ClusterExpansionOrder()));
71 }
72 else {
73 m_BMAX = max(m_BMAX, abs(part.BaryonCharge()));
74 m_QMAX = max(m_QMAX, abs(part.ElectricCharge()));
75 m_SMAX = max(m_SMAX, abs(part.Strangeness()));
76 m_CMAX = max(m_CMAX, abs(part.Charm()));
77 }
78 }
79
84
85 if (computeFluctuations) {
86 m_BMAX *= 2;
87 m_QMAX *= 2;
88 m_SMAX *= 2;
89 m_CMAX *= 2;
90 }
91
92 // Some charges may be treated grand-canonically
93 m_BMAX *= m_BCE;
94 m_QMAX *= m_QCE;
95 m_SMAX *= m_SCE;
96 m_CMAX *= m_CCE;
97
98 std::cout << "BMAX = " << m_BMAX << "\tQMAX = " << m_QMAX << "\tSMAX = " << m_SMAX << "\tCMAX = " << m_CMAX << std::endl;
99
100 m_QNMap.clear();
101 m_QNvec.resize(0);
102
103 m_Corr.resize(0);
104 m_PartialZ.resize(0);
105
106 int ind = 0;
107 for (int iB = -m_BMAX; iB <= m_BMAX; ++iB)
108 for (int iQ = -m_QMAX; iQ <= m_QMAX; ++iQ)
109 for (int iS = -m_SMAX; iS <= m_SMAX; ++iS)
110 for (int iC = -m_CMAX; iC <= m_CMAX; ++iC) {
111
112 QuantumNumbers qn(iB, iQ, iS, iC);
113 m_QNMap[qn] = ind;
114 m_QNvec.push_back(qn);
115
116 m_PartialZ.push_back(0.);
117 m_Corr.push_back(1.);
118
119 ind++;
120 }
121
122 }
123
125 m_QuantumStats = stats;
126 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
127 m_TPS->Particle(i).UseStatistics(stats);
128
129 // Only cluster expansion method supported for particles with canonically conserved charges
130 if (stats && IsParticleCanonical(m_TPS->Particle(i)))
131 m_TPS->Particle(i).SetCalculationType(IdealGasFunctions::ClusterExpansion);
132 }
133 m_PartialZ.clear();
134 //CalculateQuantumNumbersRange();
135 }
136
138 {
139 m_ConstrainMuB &= !m_BCE;
140 m_ConstrainMuQ &= !m_QCE;
141 m_ConstrainMuS &= !m_SCE;
142 m_ConstrainMuC &= !m_CCE;
144 }
145
147 {
148 m_ConstrainMuB &= !m_BCE;
149 m_ConstrainMuQ &= !m_QCE;
150 m_ConstrainMuS &= !m_SCE;
151 m_ConstrainMuC &= !m_CCE;
153 }
154
155
157 assert(m_IGFExtraConfig.MagneticField.B == 0.); // No magnetic field supported currently
158
159 m_FluctuationsCalculated = false;
160
161 if (m_PartialZ.size() == 0)
163
164 if (m_BMAX_list == 1 && m_BCE && m_QCE && m_SCE && m_CCE && !UsePartialChemicalEquilibrium()) {
165 m_Banalyt = true;
166 m_Parameters.muB = 0.0;
167 m_Parameters.muQ = 0.0;
168 m_Parameters.muS = 0.0;
169 m_Parameters.muC = 0.0;
170 }
171 else {
172 m_Banalyt = false;
173 if (m_BCE)
174 m_Parameters.muB = 0.0;
175 if (m_QCE)
176 m_Parameters.muQ = 0.0;
177 if (m_SCE)
178 m_Parameters.muS = 0.0;
179 if (m_CCE)
180 m_Parameters.muC = 0.0;
181
182 //PrepareModelGCE(); // Plan B, may work better when quantum numbers are large
183 }
184
186
187 for (size_t i = 0; i < m_densities.size(); ++i) {
188 ThermalParticle &tpart = m_TPS->Particle(i);
189 m_densities[i] = 0.;
190
191 if (!IsParticleCanonical(tpart)) {
192 m_densities[i] = tpart.Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
193 }
194 else if (tpart.Statistics() == 0
196 {
197 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
198
199 if (ind < static_cast<int>(m_Corr.size()))
200 m_densities[i] = m_Corr[ind] * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
201 }
202 else {
203 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
204 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())];
205 if (ind < static_cast<int>(m_Corr.size()))
206 m_densities[i] += m_Corr[ind] * tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
207 }
208 }
209 }
210
211 m_Calculated = true;
213 }
214
216 {
218
219 char cc[1000];
220
221 double TOL = 1.e-4;
222
223 // Checking that the CE calculation is valid
224 if (m_BCE && m_BMAX != 0) {
225 double totB = CalculateBaryonDensity() * m_Parameters.SVc;
226 if (fabs(m_Parameters.B - totB) > TOL) {
227
228 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total baryon number.\
229\n\
230Expected: %d\n\
231Obtained: %lf\n\
232\n", m_Parameters.B, totB);
233
234 printf("%s", cc);
235
236 m_ValidityLog.append(cc);
237
238 m_LastCalculationSuccessFlag = false;
239 }
240 }
241
242
243 if (m_QCE && m_QMAX != 0) {
244 double totQ = CalculateChargeDensity() * m_Parameters.SVc;
245 if (fabs(m_Parameters.Q - totQ) > TOL) {
246 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total electric charge.\
247\n\
248Expected: %d\n\
249Obtained: %lf\n\
250\n", m_Parameters.Q, totQ);
251
252 printf("%s", cc);
253
254 m_ValidityLog.append(cc);
255
256 m_LastCalculationSuccessFlag = false;
257 }
258 }
259
260
261 if (m_SCE && m_SMAX != 0) {
262 double totS = CalculateStrangenessDensity() * m_Parameters.SVc;
263 if (fabs(m_Parameters.S - totS) > TOL) {
264 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total strangeness.\
265\n\
266Expected: %d\n\
267Obtained: %lf\n\
268\n", m_Parameters.S, totS);
269
270 printf("%s", cc);
271
272 m_ValidityLog.append(cc);
273
274 m_LastCalculationSuccessFlag = false;
275 }
276 }
277
278
279 if (m_CCE && m_CMAX != 0) {
280 double totC = CalculateCharmDensity() * m_Parameters.SVc;
281 if (fabs(m_Parameters.C - totC) > TOL) {
282 sprintf(cc, "**WARNING** ThermalModelCanonical: Inaccurate calculation of total charm.\
283\n\
284Expected: %d\n\
285Obtained: %lf\n\
286\n", m_Parameters.C, totC);
287
288 printf("%s", cc);
289
290 m_ValidityLog.append(cc);
291
292 m_LastCalculationSuccessFlag = false;
293 }
294 }
295 }
296
298 {
299 if (Vc < 0.0)
300 Vc = m_Parameters.SVc;
301
304 else {
305 // Partial chemical equilibrium canonical ensemble currently works only if there is particle-antiparticle symmetry
306 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
307 int i2 = m_TPS->PdgToId(-m_TPS->Particle(i).PdgId());
308 if (i2 != -1) {
309 if (std::abs(m_Chem[i] - m_Chem[i2]) > 1.e-8) {
310 throw std::runtime_error("ThermalModelCanonical::CalculatePartitionFunctions: Partial chemical equilibrium canonical ensemble only supported if particle-antiparticle fugacities are symmetric!");
311 }
312 }
313 }
314 }
315
316 bool AllMuZero = true;
317 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
318 ThermalParticle &tpart = m_TPS->Particle(i);
319 if (IsParticleCanonical(tpart) && m_Chem[i] != 0.0)
320 {
321 AllMuZero = false;
322 break;
323 }
324 }
325
326 vector<double> Nsx(m_PartialZ.size(), 0.);
327 vector<double> Nsy(m_PartialZ.size(), 0.);
328
329 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
330 ThermalParticle &tpart = m_TPS->Particle(i);
331
332 if (!IsParticleCanonical(tpart)) {
333 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
334 if (ind != m_QNMap[QuantumNumbers(0, 0, 0, 0)]) {
335 throw std::invalid_argument("ThermalModelCanonical: neutral particle cannot have non-zero ce charges");
336 }
337 if (ind < static_cast<int>(Nsx.size()))
338 Nsx[ind] += tpart.Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
339 }
340 else if (tpart.Statistics() == 0
342 int ind = m_QNMap[QuantumNumbers(m_BCE * tpart.BaryonCharge(), m_QCE * tpart.ElectricCharge(), m_SCE * tpart.Strangeness(), m_CCE * tpart.Charm())];
343
344 if (ind < static_cast<int>(Nsx.size())) {
346 double tdens = tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, 0.);
347 Nsx[ind] += tdens * cosh(m_Chem[i] / m_Parameters.T);
348 Nsy[ind] += tdens * sinh(m_Chem[i] / m_Parameters.T);
349 }
350 // Currently only works at mu = 0!!
351 else {
352 double tdens = tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
353 Nsx[ind] += tdens;
354 }
355 }
356 }
357 else {
358 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
359 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())];
360 if (ind < static_cast<int>(Nsx.size())) {
362 double tdens = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, 0.) / static_cast<double>(n); // TODO: Check
363 Nsx[ind] += tdens * cosh(n * m_Chem[i] / m_Parameters.T);
364 Nsy[ind] += tdens * sinh(n * m_Chem[i] / m_Parameters.T);
365 }
366 // Currently only works at mu = 0!!
367 else {
368 double tdens = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]) / static_cast<double>(n); // TODO: Check
369 Nsx[ind] += tdens;
370 }
371 }
372 }
373 }
374 }
375
376 for (int i = 0; i < static_cast<int>(Nsx.size()); ++i) {
377 Nsx[i] *= Vc;
378 Nsy[i] *= Vc;
379 m_PartialZ[i] = 0.;
380 }
381
382 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));
383 if (m_Parameters.B == 0 && m_Parameters.Q == 0 && m_Parameters.S == 0 && m_Parameters.C == 0)
384 nmax = 4;
385
386
387
388 int nmaxB = max(4, m_Parameters.B);
389 int nmaxQ = max(4, m_Parameters.Q);
390 int nmaxS = max(4, m_Parameters.S);
391 int nmaxC = max(4, m_Parameters.C);
392
393 // UPDATE: allow to increase the number of interations externally
395 nmaxB = nmaxQ = nmaxS = nmaxC = nmax;
396
397
398
399
400
401 m_MultExp = 0.;
402 m_MultExpBanalyt = 0.;
403 for (size_t i = 0; i < m_PartialZ.size(); ++i) {
404 if (!m_Banalyt || m_QNvec[i].B == 0)
405 m_MultExp += Nsx[i];
406 if (m_Banalyt && (m_QNvec[i].B == 1 || m_QNvec[i].B == -1))
407 m_MultExpBanalyt += Nsx[i];
408 }
409
410 double dphiB = xMath::Pi() / nmaxB;
411 int maxB = 2 * nmaxB;
412 if (m_BMAX == 0 || m_Banalyt)
413 maxB = 1;
414
415 for (int iB = 0; iB < maxB; ++iB) {
416
417 vector<double> xlegB, wlegB;
418
419 if (m_BMAX != 0 && !m_Banalyt) {
420 double aB = iB * dphiB;
421 if (iB >= nmaxB) aB = xMath::Pi() - (iB + 1) * dphiB;
422 double bB = aB + dphiB;
424 }
425 else {
426 xlegB.resize(1);
427 xlegB[0] = 0.;
428 wlegB.resize(1);
429 wlegB[0] = 1.;
430 }
431
432
433 double dphiS = xMath::Pi() / nmaxS;
434 int maxS = 2 * nmaxS;
435 if (m_SMAX == 0)
436 maxS = 1;
437
438 for (int iS = 0; iS < maxS; ++iS) {
439 vector<double> xlegS, wlegS;
440
441 if (m_SMAX != 0) {
442 double aS = iS * dphiS;
443 if (iS >= nmaxS) aS = xMath::Pi() - (iS + 1) * dphiS;
444 double bS = aS + dphiS;
446 }
447 else {
448 xlegS.resize(1);
449 xlegS[0] = 0.;
450 wlegS.resize(1);
451 wlegS[0] = 1.;
452 }
453
454 double dphiQ = xMath::Pi() / nmaxQ;
455 int maxQ = nmaxQ;
456 if (m_QMAX == 0)
457 maxQ = 1;
458
459 for (int iQ = 0; iQ < maxQ; ++iQ) {
460 vector<double> xlegQ, wlegQ;
461
462 if (m_QMAX != 0) {
463 double aQ = iQ * dphiQ;
464 double bQ = aQ + dphiQ;
466 }
467 else {
468 xlegQ.resize(1);
469 xlegQ[0] = 0.;
470 wlegQ.resize(1);
471 wlegQ[0] = 1.;
472 }
473
474 double dphiC = xMath::Pi() / nmaxC;
475 int maxC = 2 * nmaxC;
476 if (m_CMAX == 0)
477 maxC = 1;
478
479 for (int iC = 0; iC < maxC; ++iC) {
480 vector<double> xlegC, wlegC;
481 if (m_CMAX != 0) {
482 double aC = iC * dphiC;
483 if (iC >= nmaxC) aC = xMath::Pi() - (iC + 1) * dphiC;
484 double bC = aC + dphiC;
486 }
487 else {
488 xlegC.resize(1);
489 xlegC[0] = 0.;
490 wlegC.resize(1);
491 wlegC[0] = 1.;
492 }
493
494 for (size_t iBt = 0; iBt < xlegB.size(); ++iBt) {
495 for (size_t iSt = 0; iSt < xlegS.size(); ++iSt) {
496 for (size_t iQt = 0; iQt < xlegQ.size(); ++iQt) {
497 for (size_t iCt = 0; iCt < xlegC.size(); ++iCt) {
498 vector<double> cosph(m_PartialZ.size(), 0.), sinph(m_PartialZ.size(), 0.);
499 double wx = 0., wy = 0., mx = 0., my = 0.;
500 for (size_t i = 0; i < m_PartialZ.size(); ++i) {
501 int tB = m_QNvec[i].B;
502 int tQ = m_QNvec[i].Q;
503 int tS = m_QNvec[i].S;
504 int tC = m_QNvec[i].C;
505
506 if (m_Banalyt) {
507 cosph[i] = cos(tS*xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
508 sinph[i] = sin(tS*xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
509
510 if (m_QNvec[i].B == 1) {
511 wx += Nsx[i] * cosph[i];
512 wy += Nsx[i] * sinph[i];
513 }
514 else if (m_QNvec[i].B == 0) {
515 mx += Nsx[i] * (cosph[i] - 1.);
516 }
517 }
518 else {
519 cosph[i] = cos(tB*xlegB[iBt] + tS * xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
520 mx += Nsx[i] * (cosph[i] - 1.);
521
522 if (!AllMuZero) {
523 sinph[i] = sin(tB*xlegB[iBt] + tS * xlegS[iSt] + tQ * xlegQ[iQt] + tC * xlegC[iCt]);
524 my += Nsy[i] * sinph[i];
525 }
526 }
527 }
528
529 double wmod = 0.;
530 double warg = 0.;
531
532 if (m_Banalyt) {
533 wmod = sqrt(wx*wx + wy * wy);
534 warg = atan2(wy, wx);
535 }
536
537 for (size_t iN = 0; iN < m_PartialZ.size(); ++iN) {
538 int tBg = m_Parameters.B - m_QNvec[iN].B;
539 int tQg = m_Parameters.Q - m_QNvec[iN].Q;
540 int tSg = m_Parameters.S - m_QNvec[iN].S;
541 int tCg = m_Parameters.C - m_QNvec[iN].C;
542
543 if (m_Banalyt) {
544 //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);
545 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] *
546 cos(tBg * xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt] - tBg * warg) *
547 exp(mx + 2. * wmod - m_MultExpBanalyt) *
548 xMath::BesselIexp(tBg, 2. * wmod);
549 }
550 else {
551 if (AllMuZero)
552 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);
553 else
554 m_PartialZ[iN] += wlegB[iBt] * wlegS[iSt] * wlegQ[iQt] * wlegC[iCt] * exp(mx)
555 * (cos(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * cos(my)
556 + sin(tBg*xlegB[iBt] + tSg * xlegS[iSt] + tQg * xlegQ[iQt] + tCg * xlegC[iCt]) * sin(my));
557 }
558 }
559 }
560 }
561 }
562 }
563
564 //int cind = iB * maxS * maxQ * maxC + iS * maxQ * maxC + iQ * maxC + iC;
565 //int tot = maxB * maxS * maxQ * maxC;
566 }
567 }
568 }
569 }
570
571 for (size_t iN = 0; iN < m_PartialZ.size(); ++iN) {
572 if (m_BMAX != 0 && m_BMAX != 1 && !m_Banalyt) // TODO: cross-check
573 m_PartialZ[iN] /= 2. * xMath::Pi();
574 if (m_QMAX != 0) {
575 m_PartialZ[iN] /= 2. * xMath::Pi();
576 m_PartialZ[iN] *= 2.;
577 }
578 if (m_SMAX != 0)
579 m_PartialZ[iN] /= 2. * xMath::Pi();
580 if (m_CMAX != 0)
581 m_PartialZ[iN] /= 2. * xMath::Pi();
582 }
583
584
585 m_Corr.resize(m_PartialZ.size());
586 for (size_t iN = 0; iN < m_PartialZ.size(); ++iN) {
587 m_Corr[iN] = m_PartialZ[iN] / m_PartialZ[m_QNMap[QuantumNumbers(0, 0, 0, 0)]];
588 }
589 }
590
592 {
593 ThermalParticle &tpart = m_TPS->Particle(part);
594 double ret1 = 0., ret2 = 0., ret3 = 0.;
595
596 if (!IsParticleCanonical(tpart)) {
597 return tpart.ScaledVariance(m_Parameters, m_UseWidth, m_Chem[part]);
598 }
599 else if (tpart.Statistics() == 0
601 {
602 int ind = m_QNMap[QuantumNumbers(m_BCE*tpart.BaryonCharge(), m_QCE*tpart.ElectricCharge(), m_SCE*tpart.Strangeness(), m_CCE*tpart.Charm())];
603 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())];
604
605 ret1 = 1.;
606 if (ind < static_cast<int>(m_Corr.size()) && ind2 < static_cast<int>(m_Corr.size()))
607 ret2 = m_Corr[ind2] / m_Corr[ind] * m_Parameters.SVc * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
608
609 if (ind < static_cast<int>(m_Corr.size()))
610 ret3 = -m_Corr[ind] * m_Parameters.SVc * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
611 }
612 else {
613 double ret1num = 0., ret1zn = 0.;
614 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
615 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())];
616
617 double densityClusterN = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
618
619 if (ind < static_cast<int>(m_Corr.size())) {
620 ret1num += m_Corr[ind] * n * densityClusterN;
621 ret1zn += m_Corr[ind] * densityClusterN;
622 }
623
624 for (int n2 = 1; n2 <= tpart.ClusterExpansionOrder(); ++n2) {
625 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) {
626 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())];
627 if (ind < static_cast<int>(m_Corr.size()) && ind2 < static_cast<int>(m_Corr.size()))
628 ret2 += densityClusterN * m_Corr[ind2] * m_Parameters.SVc * tpart.DensityCluster(n2, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[part]);
629 }
630 }
631 }
632
633 if (ret1zn == 0.0)
634 return 1.;
635
636 ret1 = ret1num / ret1zn;
637 ret2 = ret2 / ret1zn;
638 ret3 = -ret1zn * m_Parameters.SVc;
639 }
640 return ret1 + ret2 + ret3;
641 }
642
644 int NN = m_densities.size();
645
646 vector<double> yld(NN, 0);
647 vector<double> ret1num(NN, 0);
648 vector< vector<double> > ret2num(NN, vector<double>(NN, 0.));
649
650 for (int i = 0; i < NN; ++i)
651 yld[i] = m_densities[i] * m_Parameters.SVc;
652
653 for (int i = 0; i < NN; ++i) {
654 ThermalParticle &tpart = m_TPS->Particle(i);
655 if (!IsParticleCanonical(tpart)) {
656 ret1num[i] = tpart.ScaledVariance(m_Parameters, m_UseWidth, m_Chem[i]) * yld[i];
657 }
658 else if (tpart.Statistics() == 0
660 {
661 ret1num[i] = yld[i];
662 }
663 else {
664 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
665 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())];
666
667 double densityClusterN = tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
668
669 if (ind < static_cast<int>(m_Corr.size()))
670 ret1num[i] += m_Corr[ind] * n * densityClusterN * m_Parameters.SVc;
671 }
672 }
673 }
674
675 for (int i = 0; i < NN; ++i) {
676 for (int j = 0; j < NN; ++j) {
677 ThermalParticle &tpart1 = m_TPS->Particle(i);
678 ThermalParticle &tpart2 = m_TPS->Particle(j);
679
680 int n1max = tpart1.ClusterExpansionOrder();
681 int n2max = tpart2.ClusterExpansionOrder();
682
683 if (tpart1.Statistics() == 0 || tpart1.CalculationType() != IdealGasFunctions::ClusterExpansion)
684 n1max = 1;
685 if (tpart2.Statistics() == 0 || tpart2.CalculationType() != IdealGasFunctions::ClusterExpansion)
686 n2max = 1;
687
688 if (!IsParticleCanonical(tpart1) || !IsParticleCanonical(tpart2)) {
689 ret2num[i][j] = yld[i] * yld[j];
690 }
691 else {
692 for (int n1 = 1; n1 <= n1max; ++n1) {
693 for (int n2 = 1; n2 <= n2max; ++n2) {
694 int ind = m_QNMap[QuantumNumbers(
695 m_BCE*(n1*tpart1.BaryonCharge() + n2 * tpart2.BaryonCharge()),
696 m_QCE*(n1*tpart1.ElectricCharge() + n2 * tpart2.ElectricCharge()),
697 m_SCE*(n1*tpart1.Strangeness() + n2 * tpart2.Strangeness()),
698 m_CCE*(n1*tpart1.Charm() + n2 * tpart2.Charm()))];
699
700 double densityClusterN1 = tpart1.DensityCluster(n1, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
701 double densityClusterN2 = tpart2.DensityCluster(n2, m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[j]);
702
703 if (ind < static_cast<int>(m_Corr.size()))
704 ret2num[i][j] += m_Corr[ind] * densityClusterN1 * densityClusterN2 * m_Parameters.SVc * m_Parameters.SVc;
705 }
706 }
707 }
708 }
709 }
710
711
712 m_PrimCorrel.resize(NN);
713 for (int i = 0; i < NN; ++i)
714 m_PrimCorrel[i].resize(NN);
715 m_TotalCorrel = m_PrimCorrel;
716
717 for (int i = 0; i < NN; ++i) {
718 for (int j = 0; j < NN; ++j) {
719 m_PrimCorrel[i][j] = 0.;
720 if (i == j)
721 m_PrimCorrel[i][j] += ret1num[i] / yld[i];
722 m_PrimCorrel[i][j] += ret2num[i][j] / yld[i];
723 m_PrimCorrel[i][j] += -yld[j];
724 m_PrimCorrel[i][j] *= yld[i] / m_Parameters.SVc / m_Parameters.T;
725 if (yld[i] == 0.0)
726 m_PrimCorrel[i][j] = 0.;
727 }
728 }
729
730 for (int i = 0; i < NN; ++i) {
731 m_wprim[i] = m_PrimCorrel[i][i];
732 if (m_densities[i] > 0.) m_wprim[i] *= m_Parameters.T / m_densities[i];
733 else m_wprim[i] = 1.;
734 }
735
736 }
737
744
745 m_FluctuationsCalculated = true;
746
747 for (size_t i = 0; i < m_wprim.size(); ++i) {
748 m_skewprim[i] = 1.;
749 m_kurtprim[i] = 1.;
750 m_skewtot[i] = 1.;
751 m_kurttot[i] = 1.;
752 }
753 }
754
756 if (!m_Calculated) CalculateDensities();
757 double ret = 0.;
758
759 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
760 ThermalParticle &tpart = m_TPS->Particle(i);
761 {
762 if (!IsParticleCanonical(tpart)) {
763 ret += tpart.Density(m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i]);
764 }
765 else if (tpart.Statistics() == 0
767 int ind = m_QNMap[QuantumNumbers(tpart.BaryonCharge(), tpart.ElectricCharge(), tpart.Strangeness(), tpart.Charm())];
768
769 if (ind < static_cast<int>(m_Corr.size()))
770 ret += m_Corr[ind] * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i]);
771 }
772 else {
773 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
774 int ind = m_QNMap[QuantumNumbers(n*tpart.BaryonCharge(), n*tpart.ElectricCharge(), n*tpart.Strangeness(), n*tpart.Charm())];
775 if (ind < static_cast<int>(m_Corr.size()))
776 ret += m_Corr[ind] * tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::EnergyDensity, m_UseWidth, m_Chem[i]);
777 }
778 }
779 }
780 }
781
782 return ret;
783 }
784
786 if (!m_Calculated) CalculateDensities();
787 double ret = 0.;
788
789 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
790 ThermalParticle &tpart = m_TPS->Particle(i);
791 {
792 if (!IsParticleCanonical(tpart)) {
793 ret += tpart.Density(m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i]);
794 }
795 else if (tpart.Statistics() == 0
797 int ind = m_QNMap[QuantumNumbers(tpart.BaryonCharge(), tpart.ElectricCharge(), tpart.Strangeness(), tpart.Charm())];
798
799 if (ind < static_cast<int>(m_Corr.size()))
800 ret += m_Corr[ind] * tpart.DensityCluster(1, m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i]);
801 }
802 else {
803 for (int n = 1; n <= tpart.ClusterExpansionOrder(); ++n) {
804 int ind = m_QNMap[QuantumNumbers(n*tpart.BaryonCharge(), n*tpart.ElectricCharge(), n*tpart.Strangeness(), n*tpart.Charm())];
805
806 if (ind < static_cast<int>(m_Corr.size()))
807 ret += m_Corr[ind] * tpart.DensityCluster(n, m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, m_Chem[i]);
808 }
809 }
810 }
811 }
812
813 return ret;
814 }
815
817 {
818 double ret = (CalculateEnergyDensity() / m_Parameters.T) + (m_MultExp + m_MultExpBanalyt + log(m_PartialZ[m_QNMap[QuantumNumbers(0, 0, 0, 0)]])) / m_Parameters.SVc;
819
820 if (m_BCE)
821 ret += -m_Parameters.muB / m_Parameters.T * m_Parameters.B / m_Parameters.SVc;
822 else
823 ret += -m_Parameters.muB * CalculateBaryonDensity() / m_Parameters.T;
824
825 if (m_QCE)
826 ret += -m_Parameters.muQ / m_Parameters.T * m_Parameters.Q / m_Parameters.SVc;
827 else
828 ret += -m_Parameters.muQ * CalculateChargeDensity() / m_Parameters.T;
829
830 if (m_SCE)
831 ret += -m_Parameters.muS / m_Parameters.T * m_Parameters.S / m_Parameters.SVc;
832 else
833 ret += -m_Parameters.muS * CalculateStrangenessDensity() / m_Parameters.T;
834
835 if (m_CCE)
836 ret += -m_Parameters.muC / m_Parameters.T * m_Parameters.C / m_Parameters.SVc;
837 else
838 ret += -m_Parameters.muC * CalculateCharmDensity() / m_Parameters.T;
839
840 return ret;
841 }
842
844 {
845 return m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
846 }
847
849 {
850 return !(
851 (part.BaryonCharge() == 0 || m_BCE == 0)
852 && (part.ElectricCharge() == 0 || m_QCE == 0)
853 && (part.Strangeness() == 0 || m_SCE == 0)
854 && (part.Charm() == 0 || m_CCE == 0)
855 );
856 }
857
859 {
860 if (charge == ConservedCharge::BaryonCharge)
861 return (m_BCE != 0);
862 else if (charge == ConservedCharge::ElectricCharge)
863 return (m_QCE != 0);
864 else if (charge == ConservedCharge::StrangenessCharge)
865 return (m_SCE != 0);
866 else if (charge == ConservedCharge::CharmCharge)
867 return (m_CCE != 0);
868 return 0;
869 }
870
871 void ThermalModelCanonical::PrepareModelGCE()
872 {
873 CleanModelGCE();
874
875 m_modelgce = new ThermalModelIdeal(m_TPS, m_Parameters);
876 m_modelgce->SetUseWidth(m_UseWidth);
878
879 if (m_BCE)
880 m_Parameters.muB = 0.0;
881
882 if (!m_BCE && m_SCE)
883 m_Parameters.muS = m_Parameters.muB / 3.;
884
885 if (!m_BCE && m_QCE)
886 m_Parameters.muQ = -m_Parameters.muB / 30.;
887
888 m_Parameters.muC = 0.;
889
890 m_modelgce->SolveChemicalPotentials(m_Parameters.B, m_Parameters.Q, m_Parameters.S, m_Parameters.C,
891 m_Parameters.muB, m_Parameters.muQ, m_Parameters.muS, m_Parameters.muC,
892 static_cast<bool>(m_BCE),
893 static_cast<bool>(m_QCE),
894 static_cast<bool>(m_SCE),
895 static_cast<bool>(m_CCE));
896
897 m_Parameters.muB = m_modelgce->Parameters().muB;
898 m_Parameters.muQ = m_modelgce->Parameters().muQ;
899 m_Parameters.muS = m_modelgce->Parameters().muS;
900 m_Parameters.muC = m_modelgce->Parameters().muC;
901
902
903 // Possible alternative below
904 //double tdens = 0.;
905 //for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
906 // ThermalParticle &part = m_TPS->Particle(i);
907 // if (part.BaryonCharge() == 1)
908 // tdens += part.Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, 0.);
909 //}
910 //m_Parameters.muB = m_Parameters.T * asinh(m_Parameters.B / m_Parameters.SVc / 2. / tdens);
911 //m_Parameters.muS = m_Parameters.muB / 3.;
912 //m_Parameters.muQ = -m_Parameters.muB / 30.;
913 }
914
915 void ThermalModelCanonical::CleanModelGCE()
916 {
917 if (m_modelgce != NULL) {
918 delete m_modelgce;
919 m_modelgce = NULL;
920 }
921 }
922
923} // 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.