Thermal-FIST 1.5
Package for hadron resonance gas model applications
Loading...
Searching...
No Matches
ThermalModelBase.cpp
Go to the documentation of this file.
1/*
2 * Thermal-FIST package
3 *
4 * Copyright (c) 2014-2025 Volodymyr Vovchenko
5 *
6 * GNU General Public License (GPLv3 or later)
7 */
9
10#include <cstdio>
11#include <algorithm>
12
13#include <Eigen/Dense>
14
15#include "HRGBase/Utility.h"
17
18using namespace Eigen;
19
20using namespace std;
21
22namespace thermalfist {
23
25 m_TPS(TPS_),
26 m_Parameters(params),
27 m_UseWidth(false),
28 m_PCE(false),
29 m_Calculated(false),
30 m_FeeddownCalculated(false),
31 m_FluctuationsCalculated(false),
32 m_TemperatureDerivativesCalculated(false),
33 m_GCECalculated(false),
34 m_NormBratio(false),
35 m_QuantumStats(true),
36 m_MaxDiff(0.),
37 m_useOpenMP(0),
38 m_IGFExtraConfig()
39 {
42 }
43
44 m_QBgoal = 0.4;
45 m_SBgoal = 50.;
46 m_Chem.resize(m_TPS->Particles().size());
47 m_Volume = params.V;
48 m_densities.resize(m_TPS->Particles().size());
49 m_densitiestotal.resize(m_TPS->Particles().size());
50 m_densitiesbyfeeddown = std::vector< std::vector<double> >(ParticleDecayType::NumberOfDecayTypes, m_densitiestotal);
51
52 m_wprim.resize(m_TPS->Particles().size());
53 m_wtot.resize(m_TPS->Particles().size());
54 m_skewprim.resize(m_TPS->Particles().size());
55 m_skewtot.resize(m_TPS->Particles().size());
56 m_kurtprim.resize(m_TPS->Particles().size());
57 m_kurttot.resize(m_TPS->Particles().size());
58
59 m_ConstrainMuB = false;
60 m_ConstrainMuC = m_ConstrainMuQ = m_ConstrainMuS = true;
61
62 m_Susc.resize(4);
63 for (int i = 0; i < 4; ++i) m_Susc[i].resize(4);
64 m_dSuscdT = m_Susc;
65
66 m_NormBratio = false;
67
68 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
69 const ThermalParticle &tpart = m_TPS->Particles()[i];
70 for (size_t j = 0; j < tpart.Decays().size(); ++j) {
71 if (tpart.DecaysOriginal().size() == tpart.Decays().size() && tpart.Decays()[j].mBratio != tpart.DecaysOriginal()[j].mBratio)
72 m_NormBratio = true;
73 }
74 }
75
76 m_PrimCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(TPS()->ComponentsNumber(), 0.));
77 m_TotalCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(TPS()->ComponentsNumber(), 0.));
78 m_PrimChargesCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(4, 0.));
79 m_FinalChargesCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(4, 0.));
80
81 m_Ensemble = GCE;
82 m_InteractionModel = Ideal;
83
84 //SetStatistics(m_QuantumStats);
85 //SetCalculationType(IdealGasFunctions::Quadratures);
86 SetUseWidth(TPS()->ResonanceWidthIntegrationType());
87
88 ResetCalculatedFlags();
89
90 m_ValidityLog = "";
91 }
92
93
94 void ThermalModelBase::FillVirial(const std::vector<double>& /*ri*/)
95 {
96 }
97
99 {
100 if (!useWidth && m_TPS->ResonanceWidthIntegrationType() != ThermalParticle::ZeroWidth) {
101 m_TPS->SetResonanceWidthIntegrationType(ThermalParticle::ZeroWidth);
102 //m_TPS->ProcessDecays();
103 }
104 if (useWidth && m_TPS->ResonanceWidthIntegrationType() == ThermalParticle::ZeroWidth) {
105 m_TPS->SetResonanceWidthIntegrationType(ThermalParticle::BWTwoGamma);
106 }
107 m_UseWidth = useWidth;
108 }
109
111 {
112 m_UseWidth = (type != ThermalParticle::ZeroWidth);
113 m_TPS->SetResonanceWidthIntegrationType(type);
114 }
115
116
117 void ThermalModelBase::SetNormBratio(bool normBratio) {
118 if (normBratio != m_NormBratio) {
119 m_NormBratio = normBratio;
120 if (m_NormBratio) {
121 m_TPS->NormalizeBranchingRatios();
122 }
123 else {
124 m_TPS->RestoreBranchingRatios();
125 }
126 }
127 }
128
129
130 void ThermalModelBase::ResetChemicalPotentials() {
131 m_Parameters.muS = m_Parameters.muB / 5.;
132 m_Parameters.muQ = -m_Parameters.muB / 50.;
133 m_Parameters.muC = 0.;
134 }
135
136
138 m_Parameters = params;
139 m_Volume = m_Parameters.V;
140 ResetCalculatedFlags();
141 }
142
144 {
145 m_Parameters.T = T;
146 ResetCalculatedFlags();
147 }
148
150 {
151 m_Parameters.muB = muB;
153 ResetCalculatedFlags();
154 }
155
157 {
158 m_Parameters.muQ = muQ;
160 ResetCalculatedFlags();
161 }
162
164 {
165 m_Parameters.muS = muS;
167 ResetCalculatedFlags();
168 }
169
171 {
172 m_Parameters.muC = muC;
174 ResetCalculatedFlags();
175 }
176
177 void ThermalModelBase::SetGammaS(double gammaS)
178 {
179 m_Parameters.gammaS = gammaS;
180 ResetCalculatedFlags();
181 }
182
183 void ThermalModelBase::SetGammaC(double gammaC)
184 {
185 m_Parameters.gammaC = gammaC;
186 ResetCalculatedFlags();
187 }
188
190 {
191 m_Parameters.B = B;
192 ResetCalculatedFlags();
193 }
194
196 {
197 m_Parameters.Q = Q;
198 ResetCalculatedFlags();
199 }
200
202 {
203 m_Parameters.S = S;
204 ResetCalculatedFlags();
205 }
206
208 {
209 m_Parameters.C = C;
210 ResetCalculatedFlags();
211 }
212
214 {
215 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
216 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
217 const ThermalParticle &part1 = TPS()->Particles()[i];
218 const ThermalParticle &part2 = TPS()->Particles()[j];
219 if (part1.BaryonCharge() == 0 && part2.BaryonCharge() == 0)
220 SetVirial(i, j, 0.);
221 }
222 }
223 }
224
226 {
227 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
228 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
229 const ThermalParticle &part1 = TPS()->Particles()[i];
230 const ThermalParticle &part2 = TPS()->Particles()[j];
231 if (part1.BaryonCharge() == 0 && part2.BaryonCharge() == 0)
232 SetAttraction(i, j, 0.);
233 }
234 }
235 }
236
238 {
239 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
240 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
241 const ThermalParticle &part1 = TPS()->Particles()[i];
242 const ThermalParticle &part2 = TPS()->Particles()[j];
243 if ((part1.BaryonCharge() == 0 && part2.BaryonCharge() != 0)
244 || (part1.BaryonCharge() != 0 && part2.BaryonCharge() == 0))
245 SetVirial(i, j, 0.);
246 }
247 }
248 }
249
251 {
252 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
253 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
254 const ThermalParticle &part1 = TPS()->Particles()[i];
255 const ThermalParticle &part2 = TPS()->Particles()[j];
256 if ((part1.BaryonCharge() == 0 && part2.BaryonCharge() != 0)
257 || (part1.BaryonCharge() != 0 && part2.BaryonCharge() == 0))
258 SetAttraction(i, j, 0.);
259 }
260 }
261 }
262
264 {
265 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
266 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
267 const ThermalParticle &part1 = TPS()->Particles()[i];
268 const ThermalParticle &part2 = TPS()->Particles()[j];
269 if ((part1.BaryonCharge() > 0 && part2.BaryonCharge() > 0)
270 || (part1.BaryonCharge() < 0 && part2.BaryonCharge() < 0))
271 SetVirial(i, j, 0.);
272 }
273 }
274 }
275
277 {
278 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
279 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
280 const ThermalParticle &part1 = TPS()->Particles()[i];
281 const ThermalParticle &part2 = TPS()->Particles()[j];
282 if ((part1.BaryonCharge() > 0 && part2.BaryonCharge() > 0)
283 || (part1.BaryonCharge() < 0 && part2.BaryonCharge() < 0))
284 SetAttraction(i, j, 0.);
285 }
286 }
287 }
288
290 {
291 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
292 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
293 const ThermalParticle &part1 = TPS()->Particles()[i];
294 const ThermalParticle &part2 = TPS()->Particles()[j];
295 if ((part1.BaryonCharge() > 0 && part2.BaryonCharge() < 0)
296 || (part1.BaryonCharge() < 0 && part2.BaryonCharge() > 0))
297 SetVirial(i, j, 0.);
298 }
299 }
300 }
301
303 {
304 for (int i = 0; i < TPS()->ComponentsNumber(); ++i) {
305 for (int j = 0; j < TPS()->ComponentsNumber(); ++j) {
306 const ThermalParticle &part1 = TPS()->Particles()[i];
307 const ThermalParticle &part2 = TPS()->Particles()[j];
308 if ((part1.BaryonCharge() > 0 && part2.BaryonCharge() < 0)
309 || (part1.BaryonCharge() < 0 && part2.BaryonCharge() > 0))
310 SetAttraction(i, j, 0.);
311 }
312 }
313 }
314
315 void ThermalModelBase::SetGammaq(double gammaq)
316 {
317 m_Parameters.gammaq = gammaq;
318 ResetCalculatedFlags();
319 }
320
321
323 m_TPS = TPS_;
324 m_Chem.resize(m_TPS->Particles().size());
325 m_densities.resize(m_TPS->Particles().size());
326 m_densitiestotal.resize(m_TPS->Particles().size());
327 m_densitiesbyfeeddown = std::vector< std::vector<double> >(ParticleDecayType::NumberOfDecayTypes, m_densitiestotal);
328 m_wprim.resize(m_TPS->Particles().size());
329 m_wtot.resize(m_TPS->Particles().size());
330 m_skewprim.resize(m_TPS->Particles().size());
331 m_skewtot.resize(m_TPS->Particles().size());
332 m_kurtprim.resize(m_TPS->Particles().size());
333 m_kurttot.resize(m_TPS->Particles().size());
334 m_PrimCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(TPS()->ComponentsNumber(), 0.));
335 m_TotalCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(TPS()->ComponentsNumber(), 0.));
336 m_PrimChargesCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(4, 0.));
337 m_FinalChargesCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(4, 0.));
338 ResetCalculatedFlags();
339 }
340
342 m_QuantumStats = stats;
343 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
344 m_TPS->Particle(i).UseStatistics(stats);
345 }
346
348 {
349 if (!m_UseWidth) {
350 std::cerr << "**WARNING** ThermalModelBase::SetResonanceWidthIntegrationType: Using resonance widths is switched off!" << std::endl;
351 m_TPS->SetResonanceWidthIntegrationType(ThermalParticle::BWTwoGamma);
352 }
353 else
354 m_TPS->SetResonanceWidthIntegrationType(type);
355 }
356
358 m_Chem.resize(m_TPS->Particles().size());
359 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
360 m_Chem[i] = m_TPS->Particles()[i].BaryonCharge() * m_Parameters.muB + m_TPS->Particles()[i].Strangeness() * m_Parameters.muS + m_TPS->Particles()[i].ElectricCharge() * m_Parameters.muQ + m_TPS->Particles()[i].Charm() * m_Parameters.muC;
361 }
362
363 void ThermalModelBase::SetChemicalPotentials(const std::vector<double>& chem)
364 {
365 if (chem.size() != m_TPS->Particles().size()) {
366 std::cerr << "**WARNING** " << m_TAG << "::SetChemicalPotentials(const std::vector<double> & chem): size of chem does not match number of hadrons in the list" << std::endl;
367 return;
368 }
369 m_Chem = chem;
370 }
371
373 {
374 if (i < 0 || i >= static_cast<int>(m_Chem.size())) {
375 throw std::out_of_range("ThermalModelBase::ChemicalPotential: i is out of bounds!");
376 }
377 return m_Chem[i];
378 }
379
381 {
382 if (i < 0 || i >= static_cast<int>(m_Chem.size())) {
383 throw std::out_of_range("ThermalModelBase::SetChemicalPotential: i is out of bounds!");
384 }
385 m_Chem[i] = chem;
386 }
387
389 {
390 if (i < 0 || i >= static_cast<int>(m_Chem.size())) {
391 throw std::out_of_range("ThermalModelBase::FullIdealChemicalPotential: i is out of bounds!");
392 }
393
394 double ret = ChemicalPotential(i);
395
396 ret += MuShift(i);
397
398 const ThermalParticle& part = m_TPS->Particles()[i];
399
400 if (!(m_Parameters.gammaq == 1.)) ret += log(m_Parameters.gammaq) * part.AbsoluteQuark() * m_Parameters.T;
401 if (!(m_Parameters.gammaS == 1. || part.AbsoluteStrangeness() == 0.)) ret += log(m_Parameters.gammaS) * part.AbsoluteStrangeness() * m_Parameters.T;
402 if (!(m_Parameters.gammaC == 1. || part.AbsoluteCharm() == 0.)) ret += log(m_Parameters.gammaC) * part.AbsoluteCharm() * m_Parameters.T;
403
404 return ret;
405 }
406
408 if (m_UseWidth && m_TPS->ResonanceWidthIntegrationType() == ThermalParticle::eBW) {
409 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
410 m_TPS->Particle(i).CalculateThermalBranchingRatios(m_Parameters, m_UseWidth, m_Chem[i] + MuShift(i));
411 }
412 m_TPS->ProcessDecays();
413 }
414
415 // Primordial
416 m_densitiesbyfeeddown[static_cast<int>(Feeddown::Primordial)] = m_densities;
417
418 // According to stability flags
419 int feed_index = static_cast<int>(Feeddown::StabilityFlag);
420 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
421 m_densitiestotal[i] = m_densities[i];
422 //const std::vector< std::pair<double, int> >& decayContributions = m_TPS->Particles()[i].DecayContributionsByFeeddown()[feed_index];
423 const ThermalParticleSystem::DecayContributionsToParticle& decayContributions = m_TPS->DecayContributionsByFeeddown()[feed_index][i];
424 for (size_t j = 0; j < decayContributions.size(); ++j)
425 if (i != decayContributions[j].second)
426 m_densitiestotal[i] += decayContributions[j].first * m_densities[decayContributions[j].second];
427 }
428
429 m_densitiesbyfeeddown[feed_index] = m_densitiestotal;
430
431 // Weak, EM, strong
432 for (feed_index = static_cast<int>(Feeddown::Weak); feed_index <= static_cast<int>(Feeddown::Strong); ++feed_index) {
433 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
434 m_densitiesbyfeeddown[feed_index][i] = m_densities[i];
435 //const std::vector< std::pair<double, int> >& decayContributions = m_TPS->Particles()[i].DecayContributionsByFeeddown()[feed_index];
436 const ThermalParticleSystem::DecayContributionsToParticle& decayContributions = m_TPS->DecayContributionsByFeeddown()[feed_index][i];
437 for (size_t j = 0; j < decayContributions.size(); ++j)
438 if (i != decayContributions[j].second)
439 m_densitiesbyfeeddown[feed_index][i] += decayContributions[j].first * m_densities[decayContributions[j].second];
440 }
441 }
442
443 m_FeeddownCalculated = true;
444 }
445
446
448 {
449 if (resetInitialValues)
451 else
453 }
454
456 if (fabs(m_Parameters.muB) < 1e-6 && !m_ConstrainMuB) {
457 if (m_ConstrainMuS)
458 m_Parameters.muS = 0.;
459 if (m_ConstrainMuQ)
460 m_Parameters.muQ = 0.;
461 if (m_ConstrainMuC)
462 m_Parameters.muC = 0.;
465 return;
466 }
467 if (m_ConstrainMuB) {
468 m_Parameters.muB = xMath::mnucleon() / 2.;
469 }
470 double suppr = 10;
471 if (m_Parameters.muB > 0.150) suppr = 8.;
472 if (m_Parameters.muB > 0.300) suppr = 7.;
473 if (m_Parameters.muB > 0.450) suppr = 6.;
474 if (m_Parameters.muB > 0.600) suppr = 6.;
475 if (m_Parameters.muB > 0.750) suppr = 5.;
476 if (m_Parameters.muB > 0.900) suppr = 4.;
477 if (m_Parameters.muB > 1.000) suppr = 3.;
478 if (m_ConstrainMuS)
479 m_Parameters.muS = m_Parameters.muB / suppr;
480 if (m_ConstrainMuQ)
481 m_Parameters.muQ = -m_Parameters.muB / suppr / 10.;
482 if (m_ConstrainMuC)
483 m_Parameters.muC = -m_Parameters.muS;
484
486 }
487
489 if (fabs(m_Parameters.muB) < 1e-6 && !m_ConstrainMuB) {
490 if (m_ConstrainMuS)
491 m_Parameters.muS = 0.;
492 if (m_ConstrainMuQ)
493 m_Parameters.muQ = 0.;
494 if (m_ConstrainMuC)
495 m_Parameters.muC = 0.;
496 //m_Parameters.muS = m_Parameters.muQ = m_Parameters.muC = 0.;
499 return;
500 }
501
502 m_ConstrainMuB &= m_TPS->hasBaryons();
503 m_ConstrainMuQ &= (m_TPS->hasCharged() && m_TPS->hasBaryons());
504 m_ConstrainMuS &= m_TPS->hasStrange();
505 m_ConstrainMuC &= m_TPS->hasCharmed();
506
507 vector<double> x22(4);
508 x22[0] = m_Parameters.muB;
509 x22[1] = m_Parameters.muQ;
510 x22[2] = m_Parameters.muS;
511 x22[3] = m_Parameters.muC;
512 vector<double> x2(4), xinit(4);
513 xinit[0] = x2[0] = m_Parameters.muB;
514 xinit[1] = x2[1] = m_Parameters.muQ;
515 xinit[2] = x2[2] = m_Parameters.muS;
516 xinit[3] = x2[3] = m_Parameters.muC;
517 int iter = 0, iterMAX = 2;
518 while (iter < iterMAX) {
519 BroydenEquationsChem eqs(this);
520 BroydenJacobianChem jaco(this);
521 BroydenChem broydn(this, &eqs, &jaco);
523 broydn.Solve(x22, &crit);
524 //std::cerr << "Broyden iters: " << broydn.Iterations() << std::endl;
525 break;
526 }
527 }
528
529 bool ThermalModelBase::FixChemicalPotentialsThroughDensities(double rhoB, double rhoQ, double rhoS, double rhoC,
530 double muBinit, double muQinit, double muSinit, double muCinit,
531 bool ConstrMuB, bool ConstrMuQ, bool ConstrMuS, bool ConstrMuC) {
532 double V = Parameters().V;
534 V * rhoB, V * rhoQ, V * rhoS, V * rhoC,
535 muBinit, muQinit, muSinit, muCinit,
536 ConstrMuB, ConstrMuQ, ConstrMuS, ConstrMuC);
537 }
538
539 bool ThermalModelBase::SolveChemicalPotentials(double totB, double totQ, double totS, double totC,
540 double muBinit, double muQinit, double muSinit, double muCinit,
541 bool ConstrMuB, bool ConstrMuQ, bool ConstrMuS, bool ConstrMuC) {
543 std::cerr << "**WARNING** PCE enabled, cannot assume chemical equilibrium to do optimization..." << std::endl;
544 return false;
545 }
546
547 // TODO: Stability analysis for small densities/numbers
548 // if (abs(totB) < 1.e-6)
549 // totB = 0.;
550 // if (abs(totQ) < 1.e-6)
551 // totQ = 0.;
552 // if (abs(totS) < 1.e-6)
553 // totS = 0.;
554 // if (abs(totC) < 1.e-6)
555 // totC = 0.;
556
557 if (ConstrMuB)
558 m_Parameters.muB = muBinit;
559 if (ConstrMuS)
560 m_Parameters.muS = muSinit;
561 if (ConstrMuQ)
562 m_Parameters.muQ = muQinit;
563 if (ConstrMuC)
564 m_Parameters.muC = muCinit;
565 bool allzero = true;
566 allzero &= (totB == 0.0 && ConstrMuB) || (muBinit == 0 && !ConstrMuB);
567 allzero &= (totQ == 0.0 && ConstrMuQ) || (muQinit == 0 && !ConstrMuQ);
568 allzero &= (totS == 0.0 && ConstrMuS) || (muSinit == 0 && !ConstrMuS);
569 allzero &= (totC == 0.0 && ConstrMuC) || (muCinit == 0 && !ConstrMuC);
570
571
572 if (allzero) {
573 m_Parameters.muB = 0.;
574 m_Parameters.muS = 0.;
575 m_Parameters.muQ = 0.;
576 m_Parameters.muC = 0.;
579 return true;
580 }
581 vector<int> vConstr(4, 1);
582 vector<int> vType(4, 0);
583
584
585 vConstr[0] = m_TPS->hasBaryons() && ConstrMuB;
586 vConstr[1] = m_TPS->hasCharged() && ConstrMuQ;
587 vConstr[2] = m_TPS->hasStrange() && ConstrMuS;
588 vConstr[3] = m_TPS->hasCharmed() && ConstrMuC;
589
590 vType[0] = (int)(totB == 0.0);
591 vType[1] = (int)(totQ == 0.0);
592 vType[2] = (int)(totS == 0.0);
593 vType[3] = (int)(totC == 0.0);
594
595 vector<double> vTotals(4);
596 vTotals[0] = totB;
597 vTotals[1] = totQ;
598 vTotals[2] = totS;
599 vTotals[3] = totC;
600
601 vector<double> xin(4, 0.);
602 xin[0] = muBinit;
603 xin[1] = muQinit;
604 xin[2] = muSinit;
605 xin[3] = muCinit;
606
607 vector<double> xinactual;
608 for (int i = 0; i < 4; ++i) {
609 if (vConstr[i]) {
610 xinactual.push_back(xin[i]);
611 }
612 }
613
614 BroydenEquationsChemTotals eqs(vConstr, vType, vTotals, this);
615 BroydenJacobianChemTotals jaco(vConstr, vType, vTotals, this);
616 Broyden broydn(&eqs, &jaco);
618 broydn.Solve(xinactual, &crit);
619
620
621 //std::cerr << BaryonDensity() * Volume() << std::endl;
622
623 return (broydn.Iterations() < broydn.MaxIterations());
624 }
625
632
634 {
635 m_ValidityLog = "";
636
637 char cc[1000];
638
639 m_LastCalculationSuccessFlag = true;
640 for (size_t i = 0; i < m_densities.size(); ++i) {
641 if (m_densities[i] != m_densities[i]) {
642 m_LastCalculationSuccessFlag = false;
643
644 std::cerr << "**WARNING** Density for particle " << m_TPS->Particle(i).PdgId() << " (" << m_TPS->Particle(i).Name() << ") is NaN!\n\n";
645
646 m_ValidityLog.append(cc);
647 }
648 //m_LastCalculationSuccessFlag &= (m_densities[i] == m_densities[i]);
649 }
650 }
651
652 std::vector<double> ThermalModelBase::CalculateChargeFluctuations(const std::vector<double>& /*chgs*/, int /*order*/)
653 {
654 std::cerr << "**WARNING** " << m_TAG << "::CalculateChargeFluctuations(const std::vector<double>& chgs, int order) not implemented!" << std::endl;
655 return std::vector<double>();
656 }
657
658 std::vector<double> ThermalModelBase::CalculateGeneralizedSusceptibilities(const std::vector<std::vector<double>> &/*chgs*/)
659 {
660 throw std::runtime_error("ThermalModelBase::CalculateGeneralizedSusceptibilities: not implemented!");
661 }
662
663 double ThermalModelBase::CalculateHadronDensity() {
664 if (!m_Calculated) CalculateDensities();
665 double ret = 0.;
666 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
667 ret += m_densities[i];
668
669 return ret;
670 }
671
672 double ThermalModelBase::CalculateBaryonDensity() {
673 if (!m_Calculated) CalculateDensities();
674 double ret = 0.;
675 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
676 ret += m_TPS->Particles()[i].BaryonCharge() * m_densities[i];
677
678 return ret;
679 }
680
681 double ThermalModelBase::CalculateChargeDensity() {
682 if (!m_Calculated) CalculateDensities();
683 double ret = 0.;
684 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
685 ret += m_TPS->Particles()[i].ElectricCharge() * m_densities[i];
686
687 return ret;
688 }
689
690 double ThermalModelBase::CalculateStrangenessDensity() {
691 if (!m_Calculated) CalculateDensities();
692 double ret = 0.;
693 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
694 ret += m_TPS->Particles()[i].Strangeness() * m_densities[i];
695
696 return ret;
697 }
698
699 double ThermalModelBase::CalculateCharmDensity() {
700 if (!m_Calculated) CalculateDensities();
701 double ret = 0.;
702 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
703 ret += m_TPS->Particles()[i].Charm() * m_densities[i];
704 return ret;
705 }
706
707 double ThermalModelBase::CalculateAbsoluteBaryonDensity() {
708 if (!m_Calculated) CalculateDensities();
709 double ret = 0.;
710 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
711 ret += fabs((double)m_TPS->Particles()[i].BaryonCharge()) * m_densities[i];
712 return ret;
713 }
714
715 double ThermalModelBase::CalculateAbsoluteChargeDensity() {
716 if (!m_Calculated) CalculateDensities();
717 double ret = 0.;
718 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
719 ret += fabs((double)m_TPS->Particles()[i].ElectricCharge()) * m_densities[i];
720 return ret;
721 }
722
723 double ThermalModelBase::CalculateAbsoluteStrangenessDensity() {
724 if (!m_Calculated) CalculateDensities();
725 double ret = 0.;
726 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
727 ret += m_TPS->Particles()[i].AbsoluteStrangeness() * m_densities[i];
728 return ret;
729 }
730
731 double ThermalModelBase::CalculateAbsoluteStrangenessDensityModulo() {
732 if (!m_Calculated) CalculateDensities();
733 double ret = 0.;
734 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
735 ret += fabs((double)m_TPS->Particles()[i].Strangeness()) * m_densities[i];
736 return ret;
737 }
738
739 double ThermalModelBase::CalculateAbsoluteCharmDensity() {
740 if (!m_Calculated) CalculateDensities();
741 double ret = 0.;
742 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
743 ret += m_TPS->Particles()[i].AbsoluteCharm() * m_densities[i];
744 return ret;
745 }
746
747 double ThermalModelBase::CalculateAbsoluteCharmDensityModulo() {
748 if (!m_Calculated) CalculateDensities();
749 double ret = 0.;
750 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
751 ret += fabs((double)m_TPS->Particles()[i].Charm()) * m_densities[i];
752 return ret;
753 }
754
755 double ThermalModelBase::CalculateArbitraryChargeDensity() {
756 if (!m_Calculated) CalculateDensities();
757 double ret = 0.;
758 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i)
759 ret += m_TPS->Particles()[i].ArbitraryCharge() * m_densities[i];
760 return ret;
761 }
762
763
764 double ThermalModelBase::GetDensity(long long PDGID, const std::vector<double> *dens)
765 {
766 if (m_TPS->PdgToId(PDGID) != -1)
767 return dens->operator[](m_TPS->PdgToId(PDGID));
768
769 // 1 - Npart
770 if (PDGID == 1) return CalculateBaryonDensity();
771
772 // K0S or K0L
773 if (PDGID == 310 || PDGID == 130)
774 if (m_TPS->PdgToId(311) != -1 && m_TPS->PdgToId(-311) != -1)
775 return (dens->operator[](m_TPS->PdgToId(311)) + dens->operator[](m_TPS->PdgToId(-311))) / 2.;
776
777 // Id Pdg code has a trailing zero, try to construct a particle + anti-particle yield
778 if (PDGID % 10 == 0) {
779 long long tpdgid = PDGID / 10;
780 if (m_TPS->PdgToId(tpdgid) != -1 && m_TPS->PdgToId(-tpdgid) != -1)
781 return dens->operator[](m_TPS->PdgToId(tpdgid)) + dens->operator[](m_TPS->PdgToId(-tpdgid));
782 }
783
784 // 22122112 - nucleons
785 if (PDGID == 22122112 && m_TPS->PdgToId(2212) != -1 && m_TPS->PdgToId(2112) != -1)
786 return dens->operator[](m_TPS->PdgToId(2212)) + dens->operator[](m_TPS->PdgToId(2112));
787
788 std::cerr << "**WARNING** " << m_TAG << ": Density with PDG ID " << PDGID << " not found!" << std::endl;
789
790 return 0.;
791 }
792
793 double ThermalModelBase::GetDensity(long long PDGID, Feeddown::Type feeddown)
794 {
795 std::vector<double> *dens = NULL;
796 if (feeddown == Feeddown::Primordial)
797 dens = &m_densities;
798 else if (feeddown == Feeddown::StabilityFlag)
799 dens = &m_densitiestotal;
800 else if (static_cast<size_t>(feeddown) < m_densitiesbyfeeddown.size())
801 dens = &m_densitiesbyfeeddown[static_cast<int>(feeddown)];
802 else {
803 std::cerr << "**WARNING** " << m_TAG << ": GetDensity: Unknown feeddown: " << static_cast<int>(feeddown) << std::endl;
804 return 0.;
805 }
806
807 if (!m_Calculated)
809
810 if (feeddown != Feeddown::Primordial && !m_FeeddownCalculated)
812
813 double ret = GetDensity(PDGID, dens);
814
815 // Weak decay contributions from K0S if this particle is not in the list
816 if (feeddown == Feeddown::Weak && m_TPS->PdgToId(310) == -1) {
817 // pi0
818 if (PDGID == 111) {
819 ret += 2. * 0.308 * GetDensity(310, dens);
820 }
821 // pi+,-
822 if (PDGID == 211 || PDGID == -211) {
823 ret += 0.692 * GetDensity(310, dens);
824 }
825 }
826
827 return ret;
828 }
829
830
831 std::vector<double> ThermalModelBase::GetIdealGasDensities() const {
832 std::vector<double> ret = m_densities;
833 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
834 ret[i] = m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::ParticleDensity, m_UseWidth, m_Chem[i]);
835 }
836 return ret;
837 }
838
839 void ThermalModelBase::ResetCalculatedFlags()
840 {
841 m_Calculated = false;
842 m_FeeddownCalculated = false;
843 m_SusceptibilitiesCalculated = false;
844 m_FluctuationsCalculated = false;
845 m_TemperatureDerivativesCalculated = false;
846 m_GCECalculated = false;
847 }
848
849 double ThermalModelBase::ConservedChargeDensity(ConservedCharge::Name chg)
850 {
852 return BaryonDensity();
854 return ElectricChargeDensity();
856 return StrangenessDensity();
858 return CharmDensity();
859 return 0.0;
860 }
861
862 double ThermalModelBase::ConservedChargeDensitydT(ConservedCharge::Name chg)
863 {
864 if (!IsTemperatureDerivativesCalculated())
866
867 double ret = 0.0;
868 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
869 ret += m_TPS->Particles()[i].ConservedCharge(chg) * m_dndT[i];
870 }
871
872 return ret;
873 }
874
875 double ThermalModelBase::ChargedMultiplicity(int type)
876 {
877 if (!m_Calculated) CalculateDensities();
878 double ret = 0.0;
879 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
880 int tQ = m_TPS->Particles()[i].ElectricCharge();
881 bool fl = false;
882 if (type == 0 && tQ != 0)
883 fl = true;
884 if (type == 1 && tQ > 0)
885 fl = true;
886 if (type == -1 && tQ < 0)
887 fl = true;
888 if (fl)
889 ret += m_densities[i];
890 }
891 return ret * Volume();
892 }
893
894 double ThermalModelBase::ChargedScaledVariance(int type)
895 {
896 if (!m_FluctuationsCalculated) {
897 std::cerr << "**WARNING** " << m_TAG << ": ChargedScaledVariance(int): Fluctuations were not calculated\n";
898 return 1.;
899 }
900 double ret = 0.0;
901 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
902 int tQ = m_TPS->Particles()[i].ElectricCharge();
903 bool fl = false;
904 if (type == 0 && tQ != 0)
905 fl = true;
906 if (type == 1 && tQ > 0)
907 fl = true;
908 if (type == -1 && tQ < 0)
909 fl = true;
910 if (fl) {
911 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) {
912 int tQ2 = m_TPS->Particles()[j].ElectricCharge();
913 bool fl2 = false;
914 if (type == 0 && tQ2 != 0)
915 fl2 = true;
916 if (type == 1 && tQ2 > 0)
917 fl2 = true;
918 if (type == -1 && tQ2 < 0)
919 fl2 = true;
920
921 if (fl2) {
922 ret += m_PrimCorrel[i][j];
923 }
924 }
925 }
926 }
927 return ret * m_Parameters.T * Volume() / ChargedMultiplicity(type);
928 }
929
930 double ThermalModelBase::ChargedMultiplicityFinal(int type)
931 {
932 if (!m_Calculated) CalculateDensities();
933
934 int op = type;
935 if (type == -1)
936 op = 2;
937
938 double ret = 0.0;
939 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
940 ret += m_densities[i] * m_TPS->Particles()[i].Nch()[op];
941 }
942 return ret * Volume();
943 }
944
945 double ThermalModelBase::ChargedScaledVarianceFinal(int type)
946 {
947 if (!m_FluctuationsCalculated) {
948 std::cerr << "**WARNING** " << m_TAG << ": ChargedScaledVarianceFinal(int): Fluctuations were not calculated" << std::endl;
949 return 1.;
950 }
951 int op = type;
952 if (type == -1)
953 op = 2;
954 double ret = 0.0;
955 for (int i = 0; i < m_TPS->ComponentsNumber(); ++i) {
956 ret += m_densities[i] * Volume() * m_TPS->Particles()[i].DeltaNch()[op];
957 for (int j = 0; j < m_TPS->ComponentsNumber(); ++j) {
958 ret += m_PrimCorrel[i][j] * m_Parameters.T * Volume() * m_TPS->Particles()[i].Nch()[op] * m_TPS->Particles()[j].Nch()[op];
959 }
960 }
961 return ret / ChargedMultiplicityFinal(type);
962 }
963
965 throw std::runtime_error("ThermalModelBase::CalculateTwoParticleCorrelations: Calculation of two-particle correlations and fluctuations is not implemented");
966 }
967
968
970 {
971 // Decay contributions here are done according to Eq. (47) in nucl-th/0606036
972
973 int NN = m_densities.size();
974
975 // Fluctuations for all
976 for (int i = 0; i < NN; ++i)
977 //for(int j=0;j<NN;++j)
978 {
979 m_TotalCorrel[i][i] = m_PrimCorrel[i][i];
980 //for (int r = 0; r < m_TPS->Particles()[i].DecayContributions().size(); ++r) {
981 const ThermalParticleSystem::DecayContributionsToParticle& decayContributions = m_TPS->DecayContributionsByFeeddown()[Feeddown::StabilityFlag][i];
982 for (size_t r = 0; r < decayContributions.size(); ++r) {
983 int rr = decayContributions[r].second;
984
985 m_TotalCorrel[i][i] += 2. * m_PrimCorrel[i][rr] * decayContributions[r].first;
986
987 for (size_t r2 = 0; r2 < decayContributions.size(); ++r2) {
988 int rr2 = decayContributions[r2].second;
989 m_TotalCorrel[i][i] += m_PrimCorrel[rr][rr2] * decayContributions[r].first * decayContributions[r2].first;
990 }
991
992 // Probabilistic decays
993 m_TotalCorrel[i][i] += m_densities[rr] / m_Parameters.T * m_TPS->DecayCumulants()[i][r].first[1];
994 //m_TotalCorrel[i][i] += m_densities[rr] / m_Parameters.T * m_TPS->Particles()[i].DecayCumulants()[r].first[1];
995 //if (rr != i) { // && !m_TPS->Particles()[r].IsStable()) {
996 // double nij = 0., ni = 0., nj = 0., dnij = 0.;
997 // const auto& decayDistributions = TPS()->ResonanceFinalStatesDistributions()[rr];
998 // for (size_t br = 0; br < decayDistributions.size(); ++br) {
999 // nij += decayDistributions[br].first * decayDistributions[br].second[i] * decayDistributions[br].second[i];
1000 // ni += decayDistributions[br].first * decayDistributions[br].second[i];
1001 // nj += decayDistributions[br].first * decayDistributions[br].second[i];
1002 // }
1003 // dnij = nij - ni * nj;
1004 // m_TotalCorrel[i][i] += m_densities[rr] / Parameters().T * dnij;
1005 //}
1006 }
1007 }
1008
1009
1010 // Correlations only for stable
1011 for (int i = 0; i < NN; ++i) {
1012 if (m_TPS->Particles()[i].IsStable()) {
1013 for (int j = 0; j < NN; ++j) {
1014 if (j != i && m_TPS->Particles()[j].IsStable()) {
1015 m_TotalCorrel[i][j] = m_PrimCorrel[i][j];
1016
1017 const ThermalParticleSystem::DecayContributionsToParticle& decayContributionsI = m_TPS->DecayContributionsByFeeddown()[Feeddown::StabilityFlag][i];
1018 const ThermalParticleSystem::DecayContributionsToParticle& decayContributionsJ = m_TPS->DecayContributionsByFeeddown()[Feeddown::StabilityFlag][j];
1019
1020 for (size_t r = 0; r < decayContributionsJ.size(); ++r) {
1021 int rr = decayContributionsJ[r].second;
1022 m_TotalCorrel[i][j] += m_PrimCorrel[i][rr] * decayContributionsJ[r].first;
1023 }
1024
1025 for (size_t r = 0; r < decayContributionsI.size(); ++r) {
1026 int rr = decayContributionsI[r].second;
1027 m_TotalCorrel[i][j] += m_PrimCorrel[j][rr] * decayContributionsI[r].first;
1028 }
1029
1030 for (size_t r = 0; r < decayContributionsI.size(); ++r) {
1031 int rr = decayContributionsI[r].second;
1032
1033 for (size_t r2 = 0; r2 < decayContributionsJ.size(); ++r2) {
1034 int rr2 = decayContributionsJ[r2].second;
1035 m_TotalCorrel[i][j] += m_PrimCorrel[rr][rr2] * decayContributionsI[r].first * decayContributionsJ[r2].first;
1036 }
1037 }
1038
1039 // Contribution from probabilistic decays
1040 for (int r = 0; r < m_TPS->ComponentsNumber(); ++r) {
1041 if (r != i && r != j) { // && !m_TPS->Particles()[r].IsStable()) {
1042 double nij = 0., ni = 0., nj = 0., dnij = 0.;
1043 //const ThermalParticle &tpart = m_TPS->Particle(r);
1044 const ThermalParticleSystem::ResonanceFinalStatesDistribution &decayDistributions = m_TPS->ResonanceFinalStatesDistributions()[r];
1045 for (size_t br = 0; br < decayDistributions.size(); ++br) {
1046 nij += decayDistributions[br].first * decayDistributions[br].second[i] * decayDistributions[br].second[j];
1047 ni += decayDistributions[br].first * decayDistributions[br].second[i];
1048 nj += decayDistributions[br].first * decayDistributions[br].second[j];
1049 }
1050 dnij = nij - ni * nj;
1051 m_TotalCorrel[i][j] += m_densities[r] / m_Parameters.T * dnij;
1052 }
1053 }
1054 }
1055 }
1056 }
1057 }
1058
1059 for (int i = 0; i < NN; ++i) {
1060 m_wtot[i] = m_TotalCorrel[i][i];
1061 if (m_densitiestotal[i] > 0.) m_wtot[i] *= m_Parameters.T / m_densitiestotal[i];
1062 else m_wtot[i] = 1.;
1063 }
1064 }
1065
1067 {
1068 if (!IsFluctuationsCalculated()) {
1069 throw std::runtime_error("ThermalModelBase::TwoParticleSusceptibilityPrimordial: fluctuations were not computed beforehand!");
1070 }
1071
1072 return m_PrimCorrel[i][j] / m_Parameters.T / m_Parameters.T / xMath::GeVtoifm() / xMath::GeVtoifm() / xMath::GeVtoifm();
1073 }
1074
1076 {
1077 int i = TPS()->PdgToId(id1);
1078 int j = TPS()->PdgToId(id2);
1079
1080 if (i == -1) {
1081 std::cerr << "**WARNING** ThermalModelBase::TwoParticleSusceptibilityPrimordialByPdg: unknown pdg code " << id1 << std::endl;
1082 return 0.;
1083 }
1084 if (j == -1) {
1085 std::cerr << "**WARNING** ThermalModelBase::TwoParticleSusceptibilityPrimordialByPdg: unknown pdg code " << id2 << std::endl;
1086 return 0.;
1087 }
1088
1090 }
1091
1093 if (!IsFluctuationsCalculated() || !IsTemperatureDerivativesCalculated()) {
1094 throw std::runtime_error("ThermalModelBase::TwoParticleSusceptibilityPrimordial: temperature derivatives of fluctuations were not computed beforehand!");
1095 }
1096
1097 return m_PrimChi2sdT[i][j];
1098 }
1099
1101 {
1102 int i = TPS()->PdgToId(id1);
1103 int j = TPS()->PdgToId(id2);
1104
1105 if (i == -1) {
1106 std::cerr << "**WARNING** ThermalModelBase::TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg: unknown pdg code " << id1 << std::endl;
1107 return 0.;
1108 }
1109 if (j == -1) {
1110 std::cerr << "**WARNING** ThermalModelBase::TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg: unknown pdg code " << id2 << std::endl;
1111 return 0.;
1112 }
1113
1115 }
1116
1118 {
1119 int i1 = TPS()->PdgToId(id1);
1120 int j1 = TPS()->PdgToId(id2);
1121
1122 if (i1 == -1) {
1123 std::cerr << "**WARNING** ThermalModelBase::NetParticleSusceptibilityPrimordialByPdg: unknown pdg code " << id1 << std::endl;
1124 return 0.;
1125 }
1126 if (j1 == -1) {
1127 std::cerr << "**WARNING** ThermalModelBase::NetParticleSusceptibilityPrimordialByPdg: unknown pdg code " << id2 << std::endl;
1128 return 0.;
1129 }
1130
1131 int i2 = TPS()->PdgToId(-id1);
1132 int j2 = TPS()->PdgToId(-id2);
1133
1134 // Both particles are neutral
1135 if (i2 == -1 && j2 == -1) {
1137 }
1138
1139 // First particle species is neutral
1140 if (i2 == -1) {
1142 }
1143
1144 // Second particle species is neutral
1145 if (j2 == -1) {
1147 }
1148
1149 // Both particles have antiparticles
1151 }
1152
1154 {
1155 if (!IsFluctuationsCalculated()) {
1156 throw std::runtime_error("ThermalModelBase::TwoParticleSusceptibilityFinal: fluctuations were not computed beforehand!");
1157 }
1158
1159 if (!m_TPS->Particle(i).IsStable() || !m_TPS->Particle(j).IsStable()) {
1160 int tid = i;
1161 if (!m_TPS->Particle(j).IsStable())
1162 tid = j;
1163 throw std::runtime_error("ThermalModelBase::TwoParticleSusceptibilityFinal: Particle " + std::to_string(tid) + " is not stable! Final correlations not computed for unstable particles.");
1164 }
1165
1166 return m_TotalCorrel[i][j] / m_Parameters.T / m_Parameters.T / xMath::GeVtoifm() / xMath::GeVtoifm() / xMath::GeVtoifm();
1167 }
1168
1170 {
1171 int i = TPS()->PdgToId(id1);
1172 int j = TPS()->PdgToId(id2);
1173
1174 if (i == -1) {
1175 std::cerr << "**WARNING** ThermalModelBase::TwoParticleSusceptibilityFinalByPdg: unknown pdg code " << id1 << std::endl;
1176 return 0.;
1177 }
1178 if (j == -1) {
1179 std::cerr << "**WARNING** ThermalModelBase::TwoParticleSusceptibilityFinalByPdg: unknown pdg code " << id2 << std::endl;
1180 return 0.;
1181 }
1182
1183 return TwoParticleSusceptibilityFinal(i, j);
1184 }
1185
1187 {
1188 int i1 = TPS()->PdgToId(id1);
1189 int j1 = TPS()->PdgToId(id2);
1190
1191 if (i1 == -1) {
1192 std::cerr << "**WARNING** ThermalModelBase::NetParticleSusceptibilityFinalByPdg: unknown pdg code " << id1 << std::endl;
1193 return 0.;
1194 }
1195 if (j1 == -1) {
1196 std::cerr << "**WARNING** ThermalModelBase::NetParticleSusceptibilityFinalByPdg: unknown pdg code " << id2 << std::endl;
1197 return 0.;
1198 }
1199
1200 int i2 = TPS()->PdgToId(-id1);
1201 int j2 = TPS()->PdgToId(-id2);
1202
1203 // Both particles are neutral
1204 if (i2 == -1 && j2 == -1) {
1205 return TwoParticleSusceptibilityFinal(i1, j1);
1206 }
1207
1208 // First particle species is neutral
1209 if (i2 == -1) {
1211 }
1212
1213 // Second particle species is neutral
1214 if (j2 == -1) {
1216 }
1217
1218 // Both particles have antiparticles
1220 }
1221
1223 {
1224 if (!IsFluctuationsCalculated()) {
1225 throw std::runtime_error("ThermalModelBase::PrimordialParticleChargeSusceptibility: fluctuations were not computed beforehand!");
1226 }
1227
1228 return m_PrimChargesCorrel[i][static_cast<int>(chg)] / m_Parameters.T / m_Parameters.T / xMath::GeVtoifm() / xMath::GeVtoifm() / xMath::GeVtoifm();
1229 }
1230
1232 {
1233 int i = TPS()->PdgToId(id1);
1234 if (i == -1) {
1235 std::cerr << "**WARNING** ThermalModelBase::PrimordialParticleChargeSusceptibilityByPdg: unknown pdg code " << id1 << std::endl;
1236 return 0.;
1237 }
1238
1240 }
1241
1243 {
1244 int i1 = TPS()->PdgToId(id1);
1245 if (i1 == -1) {
1246 std::cerr << "**WARNING** ThermalModelBase::PrimordialNetParticleChargeSusceptibilityByPdg: unknown pdg code " << id1 << std::endl;
1247 return 0.;
1248 }
1249
1250 int i2 = TPS()->PdgToId(-id1);
1251 if (i2 == -1)
1253
1255 }
1256
1258 {
1259 if (!IsFluctuationsCalculated()) {
1260 throw std::runtime_error("ThermalModelBase::FinalParticleChargeSusceptibility: fluctuations were not computed beforehand!");
1261 }
1262
1263 return m_FinalChargesCorrel[i][static_cast<int>(chg)] / m_Parameters.T / m_Parameters.T / xMath::GeVtoifm() / xMath::GeVtoifm() / xMath::GeVtoifm();
1264 }
1265
1267 {
1268 int i = TPS()->PdgToId(id1);
1269 if (i == -1) {
1270 std::cerr << "**WARNING** ThermalModelBase::FinalParticleChargeSusceptibilityByPdg: unknown pdg code " << id1 << std::endl;
1271 return 0.;
1272 }
1273
1274 return FinalParticleChargeSusceptibility(i, chg);
1275 }
1276
1278 {
1279 int i1 = TPS()->PdgToId(id1);
1280 if (i1 == -1) {
1281 std::cerr << "**WARNING** ThermalModelBase::FinalNetParticleChargeSusceptibilityByPdg: unknown pdg code " << id1 << std::endl;
1282 return 0.;
1283 }
1284
1285 int i2 = TPS()->PdgToId(-id1);
1286 if (i2 == -1)
1287 return FinalParticleChargeSusceptibility(i1, chg);
1288
1290 }
1291
1293 {
1294 m_Susc.resize(4);
1295 for (int i = 0; i < 4; ++i) m_Susc[i].resize(4);
1296 m_dSuscdT = m_Susc;
1297
1298 for (int i = 0; i < 4; ++i) {
1299 for (int j = 0; j < 4; ++j) {
1300 m_Susc[i][j] = 0.;
1301 m_dSuscdT[i][j] = 0.;
1302 for (size_t k = 0; k < m_PrimCorrel.size(); ++k) {
1303 int c1 = 0;
1304 if (i == 0) c1 = m_TPS->Particles()[k].BaryonCharge();
1305 if (i == 1) c1 = m_TPS->Particles()[k].ElectricCharge();
1306 if (i == 2) c1 = m_TPS->Particles()[k].Strangeness();
1307 if (i == 3) c1 = m_TPS->Particles()[k].Charm();
1308 for (size_t kp = 0; kp < m_PrimCorrel.size(); ++kp) {
1309 int c2 = 0;
1310 if (j == 0) c2 = m_TPS->Particles()[kp].BaryonCharge();
1311 if (j == 1) c2 = m_TPS->Particles()[kp].ElectricCharge();
1312 if (j == 2) c2 = m_TPS->Particles()[kp].Strangeness();
1313 if (j == 3) c2 = m_TPS->Particles()[kp].Charm();
1314 m_Susc[i][j] += c1 * c2 * m_PrimCorrel[k][kp];
1315
1316 if (IsTemperatureDerivativesCalculated())
1317 m_dSuscdT[i][j] += c1 * c2 * m_PrimChi2sdT[k][kp];
1318 }
1319 }
1320 m_Susc[i][j] = m_Susc[i][j] / m_Parameters.T / m_Parameters.T / xMath::GeVtoifm() / xMath::GeVtoifm() / xMath::GeVtoifm();
1321 }
1322 }
1323 m_SusceptibilitiesCalculated = true;
1324 }
1325
1327 {
1328 m_ProxySusc.resize(4);
1329 for (int i = 0; i < 4; ++i)
1330 m_ProxySusc[i].resize(4);
1331
1332 // Up to 3, no charm here yet
1333 for (int i = 0; i < 3; ++i) {
1334 for (int j = 0; j < 3; ++j) {
1335 m_ProxySusc[i][j] = 0.;
1336 for (size_t k = 0; k < m_TotalCorrel.size(); ++k) {
1337 if (m_TPS->Particles()[k].IsStable()) {
1338 int c1 = 0;
1339 //if (i == 0) c1 = m_TPS->Particles()[k].BaryonCharge();
1340 if (i == 0) c1 = 1 * (m_TPS->Particles()[k].PdgId() == 2212) - 1 * (m_TPS->Particles()[k].PdgId() == -2212);
1341 if (i == 1) c1 = m_TPS->Particles()[k].ElectricCharge();
1342 //if (i == 1) c1 = 1 * (m_TPS->Particles()[k].PdgId() == 211) - 1 * (m_TPS->Particles()[k].PdgId() == -211);
1343 if (i == 2) c1 = 1 * (m_TPS->Particles()[k].PdgId() == 321) - 1 * (m_TPS->Particles()[k].PdgId() == -321);
1344 if (i == 3) c1 = m_TPS->Particles()[k].Charm();
1345 for (size_t kp = 0; kp < m_TotalCorrel.size(); ++kp) {
1346 if (m_TPS->Particles()[kp].IsStable()) {
1347 int c2 = 0;
1348 //if (j == 0) c2 = m_TPS->Particles()[kp].BaryonCharge();
1349 if (j == 0) c2 = 1 * (m_TPS->Particles()[kp].PdgId() == 2212) - 1 * (m_TPS->Particles()[kp].PdgId() == -2212);
1350 if (j == 1) c2 = m_TPS->Particles()[kp].ElectricCharge();
1351 //if (j == 1) c2 = 1 * (m_TPS->Particles()[kp].PdgId() == 211) - 1 * (m_TPS->Particles()[kp].PdgId() == -211);
1352 if (j == 2) c2 = 1 * (m_TPS->Particles()[kp].PdgId() == 321) - 1 * (m_TPS->Particles()[kp].PdgId() == -321);
1353 if (j == 3) c2 = m_TPS->Particles()[kp].Charm();
1354 m_ProxySusc[i][j] += c1 * c2 * m_TotalCorrel[k][kp];
1355 }
1356 }
1357 }
1358 }
1359 m_ProxySusc[i][j] = m_ProxySusc[i][j] / m_Parameters.T / m_Parameters.T / xMath::GeVtoifm() / xMath::GeVtoifm() / xMath::GeVtoifm();
1360 }
1361 }
1362 //std::cerr << "chi2netp/chi2skellam = " << m_ProxySusc[0][0] / (m_densitiestotal[m_TPS->PdgToId(2212)] + m_densitiestotal[m_TPS->PdgToId(-2212)]) * pow(m_Parameters.T * xMath::GeVtoifm(), 3) << std::endl;
1363 //std::cerr << "chi2netpi/chi2skellam = " << m_ProxySusc[1][1] / (m_densitiestotal[m_TPS->PdgToId(211)] + m_densitiestotal[m_TPS->PdgToId(-211)]) * pow(m_Parameters.T * xMath::GeVtoifm(), 3) << std::endl;
1364 }
1365
1367 {
1368 m_PrimChargesCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(4, 0.));
1369 m_FinalChargesCorrel = std::vector< std::vector<double> >(TPS()->ComponentsNumber(), std::vector<double>(4, 0.));
1370
1371 for (size_t i = 0; i < TPS()->ComponentsNumber(); ++i) {
1372 int c1 = 1;
1373 for (int chg = 0; chg < 4; ++chg) {
1374 for (size_t j = 0; j < TPS()->ComponentsNumber(); ++j) {
1375 int c2 = TPS()->Particle(j).ConservedCharge((ConservedCharge::Name)chg);
1376 m_PrimChargesCorrel[i][chg] += c1 * c2 * m_PrimCorrel[i][j];
1377 }
1378 }
1379 }
1380
1381 for (size_t i = 0; i < TPS()->ComponentsNumber(); ++i) {
1382 if (m_TPS->Particles()[i].IsStable()) {
1383 int c1 = 1;
1384 for (int chg = 0; chg < 4; ++chg) {
1385 for (size_t j = 0; j < TPS()->ComponentsNumber(); ++j) {
1386 if (m_TPS->Particles()[j].IsStable()) {
1387 int c2 = TPS()->Particle(j).ConservedCharge((ConservedCharge::Name)chg);
1388 m_FinalChargesCorrel[i][chg] += c1 * c2 * m_TotalCorrel[i][j];
1389 }
1390 }
1391 }
1392 }
1393 }
1394 }
1395
1396 std::vector<double> ThermalModelBase::PartialPressures() {
1397 if (!IsCalculated())
1399
1400 std::vector<double> ret(5, 0.);
1401 for (size_t i = 0; i < TPS()->ComponentsNumber(); ++i) {
1402 const ThermalParticle& tpart = TPS()->Particle(i);
1403 int ind = 0;
1404 if (tpart.BaryonCharge() == 1)
1405 ind = 1;
1406 if (tpart.BaryonCharge() == -1)
1407 ind = 2;
1408 if (tpart.BaryonCharge() > 1)
1409 ind = 3;
1410 if (tpart.BaryonCharge() < -1)
1411 ind = 4;
1412 ret[ind] += m_TPS->Particles()[i].Density(m_Parameters, IdealGasFunctions::Pressure, m_UseWidth, FullIdealChemicalPotential(i));
1413 }
1414
1415 return ret;
1416 }
1417
1419 std::cerr << "**WARNING** " << m_TAG << ": Calculation of fluctuations is not implemented" << std::endl;
1420 }
1421
1423 std::cerr << "**WARNING** " << m_TAG << ": Calculation of temperature derivatives is not implemented" << std::endl;
1424 }
1425
1426 std::vector<double> ThermalModelBase::BroydenEquationsChem::Equations(const std::vector<double>& x)
1427 {
1428 std::vector<double> ret(m_N, 0.);
1429
1430 int i1 = 0;
1431 if (m_THM->ConstrainMuB()) { m_THM->SetBaryonChemicalPotential(x[i1]); i1++; }
1432 if (m_THM->ConstrainMuQ()) { m_THM->SetElectricChemicalPotential(x[i1]); i1++; }
1433 if (m_THM->ConstrainMuS()) { m_THM->SetStrangenessChemicalPotential(x[i1]); i1++; }
1434 if (m_THM->ConstrainMuC()) { m_THM->SetCharmChemicalPotential(x[i1]); i1++; }
1435 m_THM->FillChemicalPotentials();
1436 m_THM->CalculatePrimordialDensities();
1437
1438 i1 = 0;
1439
1440 // Baryon charge
1441 if (m_THM->ConstrainMuB()) {
1442 double fBd = m_THM->CalculateBaryonDensity();
1443 double fSd = m_THM->CalculateEntropyDensity();
1444
1445 ret[i1] = (fBd / fSd - 1. / m_THM->SoverB()) * m_THM->SoverB();
1446
1447 i1++;
1448 }
1449
1450 // Electric charge
1451 if (m_THM->ConstrainMuQ()) {
1452 double fBd = m_THM->CalculateBaryonDensity();
1453 double fQd = m_THM->CalculateChargeDensity();
1454
1455 // Update: remove division by Q/B to allow for charge neutrality
1456 ret[i1] = (fQd / fBd - m_THM->QoverB());
1457 //ret[i1] = (fQd / fBd - m_THM->QoverB()) / m_THM->QoverB();
1458
1459 i1++;
1460 }
1461
1462
1463 // Strangeness
1464 if (m_THM->ConstrainMuS()) {
1465 double fSd = m_THM->CalculateStrangenessDensity();
1466 double fASd = m_THM->CalculateAbsoluteStrangenessDensity();
1467 if (fASd < 1.e-25)
1468 fASd = m_THM->CalculateAbsoluteStrangenessDensityModulo();
1469
1470 ret[i1] = fSd / fASd;
1471
1472 i1++;
1473 }
1474
1475
1476 // Charm
1477 if (m_THM->ConstrainMuC()) {
1478 double fCd = m_THM->CalculateCharmDensity();
1479 double fACd = m_THM->CalculateAbsoluteCharmDensity();
1480 if (fACd < 1.e-25)
1481 fACd = m_THM->CalculateAbsoluteCharmDensityModulo();
1482
1483 ret[i1] = fCd / fACd;
1484
1485 i1++;
1486 }
1487
1488 return ret;
1489 }
1490
1491 std::vector<double> ThermalModelBase::BroydenJacobianChem::Jacobian(const std::vector<double>& x)
1492 {
1493 int i1 = 0;
1494 // Analytic calculations of Jacobian not yet supported if entropy per baryon is involved
1495 if (m_THM->ConstrainMuB()) {
1496 throw std::runtime_error("ThermalModelBase::ConstrainChemicalPotentials: analytic calculation of the Jacobian not supported if muB is constrained");
1497 }
1498
1499 if (m_THM->ConstrainMuQ()) { m_THM->SetElectricChemicalPotential(x[i1]); i1++; }
1500 if (m_THM->ConstrainMuS()) { m_THM->SetStrangenessChemicalPotential(x[i1]); i1++; }
1501 if (m_THM->ConstrainMuC()) { m_THM->SetCharmChemicalPotential(x[i1]); i1++; }
1502 m_THM->FillChemicalPotentials();
1503 m_THM->CalculatePrimordialDensities();
1504
1505 double fBd = m_THM->CalculateBaryonDensity();
1506 double fQd = m_THM->CalculateChargeDensity();
1507 double fSd = m_THM->CalculateStrangenessDensity();
1508 double fASd = m_THM->CalculateAbsoluteStrangenessDensity();
1509 if (fASd < 1.e-25) {
1510 fASd = m_THM->CalculateAbsoluteStrangenessDensityModulo();
1511 }
1512 double fCd = m_THM->CalculateCharmDensity();
1513 double fACd = m_THM->CalculateAbsoluteCharmDensity();
1514 if (fACd < 1.e-25) {
1515 fACd = m_THM->CalculateAbsoluteCharmDensityModulo();
1516 }
1517
1518 vector<double> chi2s;
1519 chi2s.resize(m_THM->Densities().size());
1520 for (size_t i = 0; i < chi2s.size(); ++i) {
1521 chi2s[i] = m_THM->TPS()->Particle(i).chiDimensionfull(2, m_THM->Parameters(), m_THM->UseWidth(), m_THM->ChemicalPotential(i) + m_THM->MuShift(i))
1522 * xMath::GeVtoifm3();
1523 }
1524
1525 int NNN = 0;
1526 if (m_THM->ConstrainMuQ()) NNN++;
1527 if (m_THM->ConstrainMuS()) NNN++;
1528 if (m_THM->ConstrainMuC()) NNN++;
1529 MatrixXd ret(NNN, NNN);
1530
1531 i1 = 0;
1532 // Electric charge derivatives
1533 if (m_THM->ConstrainMuQ()) {
1534 int i2 = 0;
1535
1536 double d1 = 0., d2 = 0.;
1537
1538 if (m_THM->ConstrainMuQ()) {
1539 d1 = 0.;
1540 for (size_t i = 0; i < chi2s.size(); ++i)
1541 d1 += m_THM->TPS()->Particle(i).ElectricCharge() * m_THM->TPS()->Particle(i).ElectricCharge() * chi2s[i];
1542
1543 d2 = 0.;
1544 for (size_t i = 0; i < chi2s.size(); ++i)
1545 d2 += m_THM->TPS()->Particle(i).BaryonCharge() * m_THM->TPS()->Particle(i).ElectricCharge() * chi2s[i];
1546
1547 // Update: remove division by Q/B to allow for charge neutrality
1548 ret(i1, i2) = (d1 / fBd - fQd / fBd / fBd * d2);
1549 //ret(i1, i2) = (d1 / fBd - fQd / fBd / fBd * d2) / m_THM->QoverB();
1550
1551 i2++;
1552 }
1553
1554
1555 if (m_THM->ConstrainMuS()) {
1556 d1 = 0.;
1557 for (size_t i = 0; i < chi2s.size(); ++i)
1558 d1 += m_THM->TPS()->Particle(i).ElectricCharge() * m_THM->TPS()->Particle(i).Strangeness() * chi2s[i];
1559
1560 d2 = 0.;
1561 for (size_t i = 0; i < chi2s.size(); ++i)
1562 d2 += m_THM->TPS()->Particle(i).BaryonCharge() * m_THM->TPS()->Particle(i).Strangeness() * chi2s[i];
1563
1564 // Update: remove division by Q/B to allow for charge neutrality
1565 ret(i1, i2) = (d1 / fBd - fQd / fBd / fBd * d2);
1566 //ret(i1, i2) = (d1 / fBd - fQd / fBd / fBd * d2) / m_THM->QoverB();
1567
1568 i2++;
1569 }
1570
1571
1572 if (m_THM->ConstrainMuC()) {
1573 d1 = 0.;
1574 for (size_t i = 0; i < chi2s.size(); ++i)
1575 d1 += m_THM->TPS()->Particle(i).ElectricCharge() * m_THM->TPS()->Particle(i).Charm() * chi2s[i];
1576
1577 d2 = 0.;
1578 for (size_t i = 0; i < chi2s.size(); ++i)
1579 d2 += m_THM->TPS()->Particle(i).BaryonCharge() * m_THM->TPS()->Particle(i).Charm() * chi2s[i];
1580
1581 // Update: remove division by Q/B to allow for charge neutrality
1582 ret(i1, i2) = (d1 / fBd - fQd / fBd / fBd * d2);
1583 //ret(i1, i2) = (d1 / fBd - fQd / fBd / fBd * d2) / m_THM->QoverB();
1584
1585 i2++;
1586 }
1587
1588 i1++;
1589 }
1590
1591
1592 // Strangeness derivatives
1593 if (m_THM->ConstrainMuS()) {
1594 int i2 = 0;
1595
1596 double d1 = 0., d2 = 0.;
1597
1598 if (m_THM->ConstrainMuQ()) {
1599 d1 = 0.;
1600 for (size_t i = 0; i < chi2s.size(); ++i)
1601 d1 += m_THM->TPS()->Particle(i).Strangeness() * m_THM->TPS()->Particle(i).ElectricCharge() * chi2s[i];
1602
1603 d2 = 0.;
1604 for (size_t i = 0; i < chi2s.size(); ++i)
1605 d2 += m_THM->TPS()->Particle(i).AbsoluteStrangeness() * m_THM->TPS()->Particle(i).ElectricCharge() * chi2s[i];
1606
1607 ret(i1, i2) = d1 / fASd - fSd / fASd / fASd * d2;
1608
1609 i2++;
1610 }
1611
1612
1613 if (m_THM->ConstrainMuS()) {
1614 d1 = 0.;
1615 for (size_t i = 0; i < chi2s.size(); ++i)
1616 d1 += m_THM->TPS()->Particle(i).Strangeness() * m_THM->TPS()->Particle(i).Strangeness() * chi2s[i];
1617
1618 d2 = 0.;
1619 for (size_t i = 0; i < chi2s.size(); ++i)
1620 d2 += m_THM->TPS()->Particle(i).AbsoluteStrangeness() * m_THM->TPS()->Particle(i).Strangeness() * chi2s[i];
1621
1622 ret(i1, i2) = d1 / fASd - fSd / fASd / fASd * d2;
1623
1624 i2++;
1625 }
1626
1627
1628 if (m_THM->ConstrainMuC()) {
1629 d1 = 0.;
1630 for (size_t i = 0; i < chi2s.size(); ++i)
1631 d1 += m_THM->TPS()->Particle(i).Strangeness() * m_THM->TPS()->Particle(i).Charm() * chi2s[i];
1632
1633 d2 = 0.;
1634 for (size_t i = 0; i < chi2s.size(); ++i)
1635 d2 += m_THM->TPS()->Particle(i).AbsoluteStrangeness() * m_THM->TPS()->Particle(i).Charm() * chi2s[i];
1636
1637 ret(i1, i2) = d1 / fASd - fSd / fASd / fASd * d2;
1638
1639 i2++;
1640 }
1641
1642 i1++;
1643 }
1644
1645
1646 // Charm derivatives
1647 if (m_THM->ConstrainMuC()) {
1648 int i2 = 0;
1649
1650 double d1 = 0., d2 = 0.;
1651
1652 if (m_THM->ConstrainMuQ()) {
1653 d1 = 0.;
1654 for (size_t i = 0; i < chi2s.size(); ++i)
1655 d1 += m_THM->TPS()->Particle(i).Charm() * m_THM->TPS()->Particle(i).ElectricCharge() * chi2s[i];
1656
1657 d2 = 0.;
1658 for (size_t i = 0; i < chi2s.size(); ++i)
1659 d2 += m_THM->TPS()->Particle(i).AbsoluteCharm() * m_THM->TPS()->Particle(i).ElectricCharge() *chi2s[i];
1660
1661 ret(i1, i2) = d1 / fACd - fCd / fACd / fACd * d2;
1662
1663 i2++;
1664 }
1665
1666
1667 if (m_THM->ConstrainMuS()) {
1668 d1 = 0.;
1669 for (size_t i = 0; i < chi2s.size(); ++i)
1670 d1 += m_THM->TPS()->Particle(i).Charm() * m_THM->TPS()->Particle(i).Strangeness() * chi2s[i];
1671
1672 d2 = 0.;
1673 for (size_t i = 0; i < chi2s.size(); ++i)
1674 d2 += m_THM->TPS()->Particle(i).AbsoluteCharm() * m_THM->TPS()->Particle(i).Strangeness() * chi2s[i];
1675
1676 ret(i1, i2) = d1 / fACd - fCd / fACd / fACd * d2;
1677
1678 i2++;
1679 }
1680
1681
1682 if (m_THM->ConstrainMuC()) {
1683 d1 = 0.;
1684 for (size_t i = 0; i < chi2s.size(); ++i)
1685 d1 += m_THM->TPS()->Particle(i).Charm() * m_THM->TPS()->Particle(i).Charm() * chi2s[i];
1686
1687 d2 = 0.;
1688 for (size_t i = 0; i < chi2s.size(); ++i)
1689 d2 += m_THM->TPS()->Particle(i).AbsoluteCharm() * m_THM->TPS()->Particle(i).Charm() * chi2s[i];
1690
1691 ret(i1, i2) = d1 / fACd - fCd / fACd / fACd * d2;
1692
1693 i2++;
1694 }
1695
1696 i1++;
1697 }
1698
1699 std::vector<double> retVec(NNN*NNN, 0);
1700 for (int i = 0; i < NNN; ++i)
1701 for (int j = 0; j < NNN; ++j)
1702 retVec[i*NNN + j] = ret(i, j);
1703
1704
1705 return retVec;
1706 }
1707
1708 std::vector<double> ThermalModelBase::BroydenChem::Solve(const std::vector<double>& x0, BroydenSolutionCriterium * solcrit, int max_iterations)
1709 {
1710 if (m_Equations == NULL) {
1711 throw std::runtime_error("Broyden::Solve: Equations to solve not specified!");
1712 }
1713
1714 int NNN = 0;
1715 std::vector<double> xcur;
1716 if (m_THM->ConstrainMuB()) { xcur.push_back(x0[0]); NNN++; }
1717 if (m_THM->ConstrainMuQ()) { xcur.push_back(x0[1]); NNN++; }
1718 if (m_THM->ConstrainMuS()) { xcur.push_back(x0[2]); NNN++; }
1719 if (m_THM->ConstrainMuC()) { xcur.push_back(x0[3]); NNN++; }
1720 if (NNN == 0) {
1721 m_THM->FillChemicalPotentials();
1722 return xcur;
1723 }
1724
1725 m_Equations->SetDimension(NNN);
1726
1727 m_MaxIterations = max_iterations;
1728
1729 BroydenSolutionCriterium *SolutionCriterium = solcrit;
1730 bool UseDefaultSolutionCriterium = false;
1731 if (SolutionCriterium == NULL) {
1732 SolutionCriterium = new BroydenSolutionCriterium(TOL);
1733 UseDefaultSolutionCriterium = true;
1734 }
1735
1736 BroydenJacobian *JacobianInUse = m_Jacobian;
1737 bool UseDefaultJacobian = false;
1738 if (JacobianInUse == NULL || m_THM->ConstrainMuB()) {
1739 JacobianInUse = new BroydenJacobian(m_Equations);
1740 UseDefaultJacobian = true;
1741 }
1742 m_Iterations = 0;
1743 double &maxdiff = m_MaxDifference;
1744 int N = m_Equations->Dimension();
1745
1746
1747
1748 std::vector<double> tmpvec, xdeltavec = xcur;
1749 VectorXd xold(N), xnew(N), xdelta(N);
1750 VectorXd fold(N), fnew(N), fdelta(N);
1751
1752 xold = VectorXd::Map(&xcur[0], xcur.size());
1753
1754 MatrixXd Jac = Eigen::Map< Matrix<double, Dynamic, Dynamic, RowMajor> >(&JacobianInUse->Jacobian(xcur)[0], N, N);
1755
1756 bool constrmuB = m_THM->ConstrainMuB();
1757 bool constrmuQ = m_THM->ConstrainMuQ();
1758 bool constrmuS = m_THM->ConstrainMuS();
1759 bool constrmuC = m_THM->ConstrainMuC();
1760 bool repeat = false;
1761 NNN = 0;
1762 if (m_THM->ConstrainMuB()) {
1763 for (int j = 0; j < Jac.rows(); ++j)
1764 if (Jac(NNN, j) > 1.e8) { repeat = true; m_THM->ConstrainMuB(false); }
1765 double S = m_THM->CalculateEntropyDensity();
1766 double nB = m_THM->CalculateBaryonDensity();
1767 if (abs(S) < 1.e-25 || abs(nB) < 1.e-25) { repeat = true; m_THM->ConstrainMuB(false); }
1768 NNN++;
1769 }
1770 if (m_THM->ConstrainMuQ()) {
1771 for (int j = 0; j < Jac.rows(); ++j)
1772 if (Jac(NNN, j) > 1.e8) { repeat = true; m_THM->ConstrainMuQ(false); }
1773 double nQ = m_THM->CalculateChargeDensity();
1774 double nB = m_THM->CalculateBaryonDensity();
1775 if (abs(nQ) < 1.e-25 || abs(nB) < 1.e-25) { repeat = true; m_THM->ConstrainMuQ(false); }
1776 NNN++;
1777 }
1778 if (m_THM->ConstrainMuS()) {
1779 for (int j = 0; j < Jac.rows(); ++j)
1780 if (Jac(NNN, j) > 1.e8) { repeat = true; m_THM->ConstrainMuS(false); }
1781
1782 double nS = m_THM->CalculateAbsoluteStrangenessDensity();
1783 if (abs(nS) < 1.e-25) { nS = m_THM->CalculateAbsoluteStrangenessDensityModulo(); }
1784 if (abs(nS) < 1.e-25) { repeat = true; m_THM->ConstrainMuS(false); }
1785 NNN++;
1786 }
1787 if (m_THM->ConstrainMuC()) {
1788 for (int j = 0; j < Jac.rows(); ++j)
1789 if (Jac(NNN, j) > 1.e8) { repeat = true; m_THM->ConstrainMuC(false); }
1790 double nC = m_THM->CalculateAbsoluteCharmDensity();
1791 if (abs(nC) < 1.e-25) { nC = m_THM->CalculateAbsoluteCharmDensityModulo(); }
1792 if (abs(nC) < 1.e-25) { repeat = true; m_THM->ConstrainMuC(false); }
1793 NNN++;
1794 }
1795 if (repeat) {
1796 std::vector<double> ret = Solve(x0, solcrit, max_iterations);
1797 m_THM->ConstrainMuQ(constrmuB);
1798 m_THM->ConstrainMuQ(constrmuQ);
1799 m_THM->ConstrainMuS(constrmuS);
1800 m_THM->ConstrainMuC(constrmuC);
1801 return ret;
1802 }
1803
1804 if (Jac.determinant() == 0.0)
1805 {
1806 std::cerr << "**WARNING** Singular Jacobian in Broyden::Solve" << std::endl;
1807 return xcur;
1808 }
1809
1810 MatrixXd Jinv = Jac.inverse();
1811 tmpvec = m_Equations->Equations(xcur);
1812 fold = VectorXd::Map(&tmpvec[0], tmpvec.size());
1813
1814 for (m_Iterations = 1; m_Iterations < max_iterations; ++m_Iterations) {
1815 xnew = xold - Jinv * fold;
1816
1817 VectorXd::Map(&xcur[0], xcur.size()) = xnew;
1818
1819 tmpvec = m_Equations->Equations(xcur);
1820 fnew = VectorXd::Map(&tmpvec[0], tmpvec.size());
1821
1822
1823 maxdiff = 0.;
1824 for (size_t i = 0; i < xcur.size(); ++i) {
1825 maxdiff = std::max(maxdiff, fabs(fnew[i]));
1826 }
1827
1828 xdelta = xnew - xold;
1829 fdelta = fnew - fold;
1830
1831 VectorXd::Map(&xdeltavec[0], xdeltavec.size()) = xdelta;
1832
1833 if (SolutionCriterium->IsSolved(xcur, tmpvec, xdeltavec))
1834 break;
1835
1836 if (!m_UseNewton) // Use Broyden's method
1837 {
1838
1839 double norm = 0.;
1840 for (int i = 0; i < N; ++i)
1841 for (int j = 0; j < N; ++j)
1842 norm += xdelta[i] * Jinv(i, j) * fdelta[j];
1843 VectorXd p1(N);
1844 p1 = (xdelta - Jinv * fdelta);
1845 for (int i = 0; i < N; ++i) p1[i] *= 1. / norm;
1846 Jinv = Jinv + (p1 * xdelta.transpose()) * Jinv;
1847 }
1848 else // Use Newton's method
1849 {
1850 Jac = Eigen::Map< Matrix<double, Dynamic, Dynamic, RowMajor> >(&JacobianInUse->Jacobian(xcur)[0], N, N);
1851 Jinv = Jac.inverse();
1852 }
1853
1854 xold = xnew;
1855 fold = fnew;
1856 }
1857
1858 if (m_Iterations == max_iterations) {
1859 std::cerr << "**WARNING** Reached maximum number of iterations in Broyden procedure" << std::endl;
1860 }
1861
1862 if (UseDefaultSolutionCriterium) {
1863 delete SolutionCriterium;
1864 SolutionCriterium = NULL;
1865 }
1866 if (UseDefaultJacobian) {
1867 delete JacobianInUse;
1868 JacobianInUse = NULL;
1869 }
1870 return xcur;
1871 }
1872
1873
1874 ThermalModelBase::BroydenEquationsChemTotals::BroydenEquationsChemTotals(const std::vector<int>& vConstr, const std::vector<int>& vType, const std::vector<double>& vTotals, ThermalModelBase * model) :
1875 BroydenEquations(), m_Constr(vConstr), m_Type(vType), m_Totals(vTotals), m_THM(model)
1876 {
1877 m_N = 0;
1878 for (size_t i = 0; i < m_Constr.size(); ++i)
1879 m_N += m_Constr[i];
1880 }
1881
1882 std::vector<double> ThermalModelBase::BroydenEquationsChemTotals::Equations(const std::vector<double>& x)
1883 {
1884 std::vector<double> ret(m_N, 0.);
1885
1886 int i1 = 0;
1887 for (int i = 0; i < 4; ++i) {
1888 if (m_Constr[i]) {
1889 if (i == 0) m_THM->SetBaryonChemicalPotential(x[i1]);
1890 if (i == 1) m_THM->SetElectricChemicalPotential(x[i1]);
1891 if (i == 2) m_THM->SetStrangenessChemicalPotential(x[i1]);
1892 if (i == 3) m_THM->SetCharmChemicalPotential(x[i1]);
1893 i1++;
1894 }
1895 }
1896 m_THM->FillChemicalPotentials();
1897 m_THM->CalculatePrimordialDensities();
1898
1899 vector<double> dens(4, 0.), absdens(4, 0.);
1900 if (m_Constr[0]) {
1901 dens[0] = m_THM->CalculateBaryonDensity();
1902 absdens[0] = m_THM->CalculateAbsoluteBaryonDensity();
1903 }
1904 if (m_Constr[1]) {
1905 dens[1] = m_THM->CalculateChargeDensity();
1906 absdens[1] = m_THM->CalculateAbsoluteChargeDensity();
1907 }
1908 if (m_Constr[2]) {
1909 dens[2] = m_THM->CalculateStrangenessDensity();
1910 absdens[2] = m_THM->CalculateAbsoluteStrangenessDensity();
1911 if (absdens[2] < 1.e-25) {
1912 absdens[2] = m_THM->CalculateAbsoluteStrangenessDensityModulo();
1913 }
1914 }
1915 if (m_Constr[3]) {
1916 dens[3] = m_THM->CalculateCharmDensity();
1917 absdens[3] = m_THM->CalculateAbsoluteCharmDensity();
1918 if (absdens[3] < 1.e-25) {
1919 absdens[3] = m_THM->CalculateAbsoluteCharmDensityModulo();
1920 }
1921 }
1922
1923 i1 = 0;
1924
1925 for (int i = 0; i < 4; ++i) {
1926 if (m_Constr[i]) {
1927 if (m_Type[i] == 0)
1928 ret[i1] = (dens[i] * m_THM->Parameters().V - m_Totals[i]) / m_Totals[i];
1929 else
1930 ret[i1] = dens[i] / absdens[i];
1931 i1++;
1932 }
1933 }
1934
1935 return ret;
1936 }
1937
1938 std::vector<double> ThermalModelBase::BroydenJacobianChemTotals::Jacobian(const std::vector<double>& x)
1939 {
1940 int NNN = 0;
1941 for (int i = 0; i < 4; ++i) NNN += m_Constr[i];
1942
1943 int i1 = 0;
1944
1945 for (int i = 0; i < 4; ++i) {
1946 if (m_Constr[i]) {
1947 if (i == 0) m_THM->SetBaryonChemicalPotential(x[i1]);
1948 if (i == 1) m_THM->SetElectricChemicalPotential(x[i1]);
1949 if (i == 2) m_THM->SetStrangenessChemicalPotential(x[i1]);
1950 if (i == 3) m_THM->SetCharmChemicalPotential(x[i1]);
1951 i1++;
1952 }
1953 }
1954
1955 m_THM->FillChemicalPotentials();
1956 m_THM->CalculatePrimordialDensities();
1957
1958 vector<double> dens(4, 0.), absdens(4, 0.);
1959 if (m_Constr[0]) {
1960 dens[0] = m_THM->CalculateBaryonDensity();
1961 absdens[0] = m_THM->CalculateAbsoluteBaryonDensity();
1962 }
1963 if (m_Constr[1]) {
1964 dens[1] = m_THM->CalculateChargeDensity();
1965 absdens[1] = m_THM->CalculateAbsoluteChargeDensity();
1966 }
1967 if (m_Constr[2]) {
1968 dens[2] = m_THM->CalculateStrangenessDensity();
1969 absdens[2] = m_THM->CalculateAbsoluteStrangenessDensity();
1970 if (absdens[2] < 1.e-25) {
1971 absdens[2] = m_THM->CalculateAbsoluteStrangenessDensityModulo();
1972 }
1973 }
1974 if (m_Constr[3]) {
1975 dens[3] = m_THM->CalculateCharmDensity();
1976 absdens[3] = m_THM->CalculateAbsoluteCharmDensity();
1977 if (absdens[3] < 1.e-25) {
1978 absdens[3] = m_THM->CalculateAbsoluteCharmDensityModulo();
1979 }
1980 }
1981
1982 vector<double> chi2s;
1983 chi2s.resize(m_THM->Densities().size());
1984 for (size_t i = 0; i < chi2s.size(); ++i) {
1985 chi2s[i] = m_THM->TPS()->Particle(i).chiDimensionfull(2, m_THM->Parameters(), m_THM->UseWidth(), m_THM->ChemicalPotential(i) + m_THM->MuShift(i))
1986 * xMath::GeVtoifm3();
1987 }
1988
1989 vector< vector<double> > deriv(4, vector<double>(4)), derivabs(4, vector<double>(4));
1990 for (int i = 0; i < 4; ++i)
1991 for (int j = 0; j < 4; ++j) {
1992 deriv[i][j] = 0.;
1993 for (size_t part = 0; part < chi2s.size(); ++part)
1994 deriv[i][j] += m_THM->TPS()->Particles()[part].GetCharge(i) * m_THM->TPS()->Particles()[part].GetCharge(j) * chi2s[part];
1995
1996 derivabs[i][j] = 0.;
1997 for (size_t part = 0; part < chi2s.size(); ++part)
1998 derivabs[i][j] += m_THM->TPS()->Particles()[part].GetAbsCharge(i) * m_THM->TPS()->Particles()[part].GetCharge(j) * chi2s[part];
1999 }
2000
2001
2002 MatrixXd ret(NNN, NNN);
2003
2004 i1 = 0;
2005
2006 for (int i = 0; i < 4; ++i) {
2007 if (m_Constr[i]) {
2008 int i2 = 0;
2009 for (int j = 0; j < 4; ++j)
2010 if (m_Constr[j]) {
2011 ret(i1, i2) = 0.;
2012 if (m_Type[i] == 0)
2013 ret(i1, i2) = deriv[i][j] * m_THM->Parameters().V / m_Totals[i];
2014 else
2015 ret(i1, i2) = deriv[i][j] / absdens[i] - dens[i] / absdens[i] / absdens[i] * derivabs[i][j];
2016 i2++;
2017 }
2018 i1++;
2019 }
2020 }
2021
2022 std::vector<double> retVec(NNN*NNN, 0);
2023 for (int i = 0; i < NNN; ++i)
2024 for (int j = 0; j < NNN; ++j)
2025 retVec[i*NNN + j] = ret(i, j);
2026
2027
2028 return retVec;
2029 }
2030
2031 void ThermalModelBase::SetDensityModelForParticleSpecies(int i, GeneralizedDensity *density_model)
2032 {
2033 if (i >= 0 && i < ComponentsNumber()) {
2034 TPS()->Particle(i).SetGeneralizedDensity(density_model);
2035 } else {
2036 std::cerr << "**WARNING** ThermalModelBase::SetDensityModelForParticleSpecies(): Particle id " << i << " is outside the range!" << std::endl;
2037 }
2038 }
2039
2040 void ThermalModelBase::SetDensityModelForParticleSpeciesByPdg(long long PDGID, GeneralizedDensity *density_model)
2041 {
2042 int id = TPS()->PdgToId(PDGID);
2043 if (id != -1) {
2044 TPS()->Particle(id).SetGeneralizedDensity(density_model);
2045 } else {
2046 std::cerr << "**WARNING** ThermalModelBase::SetDensityModelForParticleSpeciesByPdg(): Pdg code " << PDGID << " is outside the range!" << std::endl;
2047 }
2048 }
2049
2050 void ThermalModelBase::ClearDensityModels() {
2051 for (auto& particle : TPS()->Particles())
2052 particle.ClearGeneralizedDensity();
2053 }
2054
2055 void ThermalModelBase::SetMagneticField(double B, int lmax) {
2056 m_IGFExtraConfig.MagneticField.B = B;
2057 if (lmax >= 0)
2058 m_IGFExtraConfig.MagneticField.lmax = lmax;
2059 for (auto& particle : TPS()->Particles())
2060 particle.SetMagneticField(B, m_IGFExtraConfig.MagneticField.lmax);
2061 RecomputeThresholdsDueToMagneticField();
2062 ResetCalculatedFlags();
2063 }
2064
2065 void ThermalModelBase::RecomputeThresholdsDueToMagneticField() {
2066 const double &B = m_IGFExtraConfig.MagneticField.B;
2067 for (int ipart = 0; ipart < TPS()->ComponentsNumber(); ++ipart) {
2068 ThermalParticle &part = TPS()->Particle(ipart);
2069 if (!part.ZeroWidthEnforced() || part.ElectricCharge()==0) {
2070 if (B == 0.) {
2071 part.CalculateAndSetDynamicalThreshold();
2072 part.FillCoefficientsDynamical();
2073 } else if (!part.ZeroWidthEnforced()) {
2074 double szmax = (part.Degeneracy() - 1.) / 2.;
2075 if (szmax > 0.5) {
2076 double mmin = sqrt(2. * abs(part.ElectricCharge()) * B * (szmax - 0.5));
2077 part.CalculateAndSetDynamicalThreshold();
2078 if (mmin > part.DecayThresholdMassDynamical()) {
2079 part.SetDecayThresholdMassDynamical(mmin);
2080 }
2081 part.FillCoefficientsDynamical();
2082 }
2083 }
2084 }
2085 }
2086 }
2087
2088 void ThermalModelBase::ClearMagneticField() {
2089 m_IGFExtraConfig.MagneticField.B = 0.;
2090 for (auto& particle : TPS()->Particles())
2091 particle.ClearMagneticField();
2092 ResetCalculatedFlags();
2093 }
2094
2095 double ThermalModelBase::Susc(ConservedCharge::Name i, ConservedCharge::Name j) const {
2096 assert(IsSusceptibilitiesCalculated());
2097 return m_Susc[i][j];
2098 }
2099
2100 double ThermalModelBase::dSuscdT(ConservedCharge::Name i, ConservedCharge::Name j) const {
2101 assert(IsSusceptibilitiesCalculated() && IsTemperatureDerivativesCalculated());
2102 return m_dSuscdT[i][j];
2103 }
2104
2105 double ThermalModelBase::ProxySusc(ConservedCharge::Name i, ConservedCharge::Name j) const {
2106 assert(IsFluctuationsCalculated());
2107 return m_ProxySusc[i][j];
2108 }
2109
2110 double ThermalModelBase::GetdndT(int i) const {
2111 if (IsTemperatureDerivativesCalculated() && i >= 0 && i < ComponentsNumber())
2112 return m_dndT[i];
2113 else {
2114 std::cerr << "**WARNING** ThermalModelBase::GetdndT(): Temperature derivatives are not calculated or particle id " << i << " is outside the range!" << std::endl;
2115 return 0.;
2116 }
2117 return 0.;
2118 }
2119
2121 double ret = 0.;
2122 for (size_t k = 0; k < m_PrimCorrel.size(); ++k) {
2123 int c1 = m_TPS->Particles()[k].ConservedCharge(i);
2124 for (size_t kp = 0; kp < m_PrimCorrel.size(); ++kp) {
2125 int c2 = m_TPS->Particles()[kp].ConservedCharge(j);
2126 ret += c1 * c2 * m_PrimCorrel[k][kp];
2127 }
2128 }
2129 return ret / xMath::GeVtoifm() / xMath::GeVtoifm() / xMath::GeVtoifm();
2130 }
2131
2132 double ThermalModelBase::CalculateSpecificHeatChem()
2133 {
2134 double dedT = CalculateEnergyDensityDerivativeT(); // fm^-3
2135
2136 // T * ds/dT = (de/dT - mu_i d n_i / d T)
2137 double ret = dedT;
2138 for(int i = 0; i < TPS()->ComponentsNumber(); ++i) {
2139 ret -= m_Chem[i] * m_dndT[i];
2140 }
2141 return ret;
2142 }
2143
2144 // Auxiliary function to compute the susceptibility matrix
2145 MatrixXd GetSusceptibilityMatrix(ThermalModelBase* model, const vector<int>& ConservedDensities) {
2146 int N = 0;
2147 for (int i = 0; i < 4; ++i)
2148 if (ConservedDensities[i] == 1)
2149 N++;
2150
2151 assert(N > 0);
2152
2153 MatrixXd ret(N, N);
2154
2155 // Pi(i,j) = d P / dmui dmuj
2156 {
2157 int i1 = 0;
2158 for(int i = 0; i < 4; ++i) {
2159 if (ConservedDensities[i] == 0)
2160 continue;
2161 int i2 = 0;
2162 for(int j = 0; j < 4; ++j) {
2163 if (ConservedDensities[j] == 0)
2164 continue;
2166 i2++;
2167 }
2168 i1++;
2169 }
2170 }
2171
2172 return ret;
2173 }
2174
2175 // Auxiliary function to compute the Hessian matrix of pressure
2176 MatrixXd GetPressureHessianMatrix(ThermalModelBase* model, const vector<int>& ConservedDensities) {
2177 int N = 1;
2178 for (int i = 0; i < 4; ++i)
2179 if (ConservedDensities[i] == 1)
2180 N++;
2181
2182 MatrixXd ret(N, N);
2183
2184 // Compute the Hessian matrix of pressure
2185 // Pi(0,0) = d^2 P / d T^2 = ds/dT
2186 ret(0, 0) = model->SpecificHeatChem() / model->Parameters().T; // fm^-3 GeV^-1
2187
2188 // Pi(0,i) = Pi(i,0) = d^2 P / dT dmui = d ni / d T
2189 {
2190 int i1 = 0;
2191 for(int i = 0; i < 4; ++i) {
2192 if (ConservedDensities[i] == 0)
2193 continue;
2194 ret(0, i1 + 1) = ret(i1 + 1, 0) = model->ConservedChargeDensitydT((ConservedCharge::Name)i);
2195 i1++;
2196 }
2197 }
2198
2199 // Pi(i,j) = d^2 P / dmui dmuj
2200 {
2201 int i1 = 0;
2202 for(int i = 0; i < 4; ++i) {
2203 if (ConservedDensities[i] == 0)
2204 continue;
2205 int i2 = 0;
2206 for(int j = 0; j < 4; ++j) {
2207 if (ConservedDensities[j] == 0)
2208 continue;
2210 i2++;
2211 }
2212 i1++;
2213 }
2214 }
2215
2216 return ret;
2217 }
2218
2219 double ThermalModelBase::CalculateAdiabaticSpeedOfSoundSquared(bool rhoBconst, bool rhoQconst, bool rhoSconst, bool rhoCconst) {
2220 double ret = 0.;
2221
2222 double T = Parameters().T;
2223
2224 std::vector<double> mus = {Parameters().muB, Parameters().muQ, Parameters().muS, Parameters().muC};
2225 std::vector<int> ConservedCharges = {static_cast<int>(rhoBconst), static_cast<int>(rhoQconst), static_cast<int>(rhoSconst), static_cast<int>(rhoCconst)};
2226
2227 if (!TPS()->hasBaryons())
2228 ConservedCharges[0] = 0;
2229 if (!TPS()->hasCharged())
2230 ConservedCharges[1] = 0;
2231 if (!TPS()->hasStrange())
2232 ConservedCharges[2] = 0;
2233 if (!TPS()->hasCharmed())
2234 ConservedCharges[3] = 0;
2235
2236
2237 // Cross-check (if charge is not conserved, its chemical potential should be zero)
2238 for(int i = 0; i < ConservedCharges.size(); ++i)
2239 assert(ConservedCharges[i] == 1 || mus[i] == 0.0);
2240
2241 // Zero chemical potentials
2242 {
2243 bool AllMuZero = true;
2244 for(int i = 0; i < ConservedCharges.size(); ++i)
2245 AllMuZero &= (ConservedCharges[i] == 0 || mus[i] == 0.0);
2246 if (AllMuZero) {
2247 // At mu = 0 cs2 = s(T) / e'(T)
2248 // Calculate temperature derivatives if not already
2249 if (!IsTemperatureDerivativesCalculated())
2251 return EntropyDensity() / CalculateEnergyDensityDerivativeT();
2252 //return EntropyDensity() / SpecificHeatChem();
2253 }
2254 }
2255
2256 int Ndens = 0;
2257 for(int i = 0; i < 4; ++i)
2258 if (ConservedCharges[i] == 1)
2259 Ndens++;
2260
2261 // Compute fluctuations if not already
2262 if (!IsFluctuationsCalculated())
2264
2265 // Zero temperature
2266 if (T == 0.) {
2267
2268 MatrixXd chi2matr = GetSusceptibilityMatrix(this, ConservedCharges);
2269 MatrixXd chi2inv = chi2matr.inverse();
2270 ret = 0.;
2271 int i1 = 0, i2 = 0;
2272 for(int i = 0; i < 4; ++i) {
2273 if (ConservedCharges[i] == 0)
2274 continue;
2275 i2 = 0;
2276 for(int j = 0; j < 4; ++j) {
2277 if (ConservedCharges[j] == 0)
2278 continue;
2279
2280 auto chg1 = (ConservedCharge::Name)i;
2281 auto chg2 = (ConservedCharge::Name)j;
2282 ret += ConservedChargeDensity(chg1) * ConservedChargeDensity(chg2) * chi2inv(i1, i2);
2283 i2++;
2284 }
2285 i1++;
2286 }
2287
2288 ret /= (EnergyDensity() + Pressure());
2289 return ret;
2290 }
2291
2292 // General case
2293
2294 // Calculate temperature derivatives if not already
2295 if (!IsTemperatureDerivativesCalculated())
2297
2298 // double dedT = SpecificHeatChem(); // fm^-3
2299 // Compute the Hessian matrix of pressure
2300 MatrixXd PiMatr = GetPressureHessianMatrix(this, ConservedCharges);
2301
2302 // Compute the inverse of the Hessian matrix
2303 MatrixXd PiMatrInv = PiMatr.inverse();
2304 // Prepare the vector x
2305 double s = EntropyDensity();
2306 vector<double> x = {s};
2307 for(int i = 0; i < ConservedCharges.size(); ++i)
2308 if (ConservedCharges[i] == 1)
2309 x.push_back(ConservedChargeDensity((ConservedCharge::Name)i));
2310
2311 ret = 0.;
2312 for(int i = 0; i < x.size(); ++i)
2313 for(int j = 0; j < x.size(); ++j)
2314 ret += x[i] * x[j] * PiMatrInv(i,j);
2315 ret /= (EnergyDensity() + Pressure());
2316
2317 return ret;
2318 }
2319
2320 double ThermalModelBase::CalculateIsothermalSpeedOfSoundSquared(bool rhoBconst, bool rhoQconst, bool rhoSconst, bool rhoCconst) {
2321 double ret = 0.;
2322
2323 double T = Parameters().T;
2324
2325 std::vector<double> mus = {Parameters().muB, Parameters().muQ, Parameters().muS, Parameters().muC};
2326 std::vector<int> ConservedCharges = {static_cast<int>(rhoBconst), static_cast<int>(rhoQconst), static_cast<int>(rhoSconst), static_cast<int>(rhoCconst)};
2327
2328 if (!TPS()->hasBaryons())
2329 ConservedCharges[0] = 0;
2330 if (!TPS()->hasCharged())
2331 ConservedCharges[1] = 0;
2332 if (!TPS()->hasStrange())
2333 ConservedCharges[2] = 0;
2334 if (!TPS()->hasCharmed())
2335 ConservedCharges[3] = 0;
2336
2337 // If all chemical potentials are zero, c_T is undefined (return 0.)
2338 {
2339 bool AllMuZero = true;
2340 for(int i = 0; i < ConservedCharges.size(); ++i)
2341 AllMuZero &= (ConservedCharges[i] == 0 || mus[i] == 0.0);
2342 if (AllMuZero) {
2343 return 0.;
2344 }
2345 }
2346
2347
2348 // Cross-check (if charge is not conserved, its chemical potential should be zero)
2349 for(int i = 0; i < ConservedCharges.size(); ++i)
2350 assert(ConservedCharges[i] == 1 || mus[i] == 0.0);
2351
2352 auto chi2Matr = GetSusceptibilityMatrix(this, ConservedCharges);
2353 // Compute the inverse of the susceptibility matrix
2354 MatrixXd chi2inv= chi2Matr.inverse();
2355
2356 double numerator = 0.;
2357 double denominator = 0.;
2358 int i1 = 0, i2 = 0;
2359 for(int i = 0; i < 4; ++i) {
2360 \
2361 if (ConservedCharges[i] == 0)
2362 continue;
2363 denominator += mus[i] * ConservedChargeDensity((ConservedCharge::Name)i);
2364 i2 = 0;
2365 for(int j = 0; j < 4; ++j) {
2366 if (ConservedCharges[j] == 0)
2367 continue;
2368
2369 auto chg1 = (ConservedCharge::Name)i;
2370 auto chg2 = (ConservedCharge::Name)j;
2371
2372 numerator += ConservedChargeDensity(chg1) * ConservedChargeDensity(chg2) * chi2inv(i1, i2);
2373 denominator += T * ConservedChargeDensitydT(chg1) * chi2inv(i1, i2) * ConservedChargeDensity(chg2);
2374
2375 ret += ConservedChargeDensity(chg1) * ConservedChargeDensity(chg2) * chi2inv(i1, i2);
2376 i2++;
2377 }
2378 i1++;
2379 }
2380
2381 return numerator / denominator;
2382 }
2383
2384 double ThermalModelBase::CalculateHeatCapacity(bool rhoBconst, bool rhoQconst, bool rhoSconst, bool rhoCconst) {
2385 double T = Parameters().T;
2386
2387 assert(T > 0.);
2388
2389 std::vector<double> mus = {Parameters().muB, Parameters().muQ, Parameters().muS, Parameters().muC};
2390 std::vector<int> ConservedCharges = {static_cast<int>(rhoBconst), static_cast<int>(rhoQconst), static_cast<int>(rhoSconst), static_cast<int>(rhoCconst)};
2391
2392 if (!TPS()->hasBaryons())
2393 ConservedCharges[0] = 0;
2394 if (!TPS()->hasCharged())
2395 ConservedCharges[1] = 0;
2396 if (!TPS()->hasStrange())
2397 ConservedCharges[2] = 0;
2398 if (!TPS()->hasCharmed())
2399 ConservedCharges[3] = 0;
2400
2401 // Calculate the numerator
2402 double TdsdT = SpecificHeatChem();
2403
2404 // If all chemical potentials are zero, retur TdsT
2405 {
2406 bool AllMuZero = true;
2407 for(int i = 0; i < ConservedCharges.size(); ++i)
2408 AllMuZero &= (ConservedCharges[i] == 0 || mus[i] == 0.0);
2409 if (AllMuZero) {
2410 return TdsdT;
2411 }
2412 }
2413
2414 if (!IsFluctuationsCalculated())
2416
2417 // Calculate the denominator
2418 double denominator = 1.;
2419 auto PiMatr = GetPressureHessianMatrix(this, ConservedCharges);
2420 auto PiMatrInv = PiMatr.inverse();
2421 int N = PiMatr.rows();
2422 for(int i = 1; i < N; ++i) {
2423 denominator -= PiMatr(0, i) * PiMatrInv(i, 0);
2424 }
2425
2426 double heatCapacity = TdsdT / denominator;
2427 return heatCapacity;
2428 }
2429
2430} // namespace thermalfist
map< string, double > params
Contains some helper functions.
Sub-class where it is determined whether the required accuracy is achieved in the Broyden's method.
Definition Broyden.h:144
Abstract class which defines the system of non-linear equations to be solved by the Broyden's method.
Definition Broyden.h:32
Class implementing the Broyden method to solve a system of non-linear equations.
Definition Broyden.h:132
virtual std::vector< double > Solve(const std::vector< double > &x0, BroydenSolutionCriterium *solcrit=NULL, int max_iterations=MAX_ITERS)
Definition Broyden.cpp:83
int Iterations() const
Definition Broyden.h:229
int MaxIterations() const
Definition Broyden.h:234
Class which implements calculation of the Jacobian needed for the Broyden's method.
Definition Broyden.h:77
static bool PrintDisclaimer()
Definition Utility.cpp:47
static bool DisclaimerPrinted
Definition Utility.h:27
Implements the possibility of a generalized calculation of the densities. For example,...
Abstract base class for an HRG model implementation.
virtual void CalculateSusceptibilityMatrix()
Calculates the conserved charges susceptibility matrix.
virtual double FullIdealChemicalPotential(int i) const
Chemical potential entering the ideal gas expressions of particle species i.
virtual void FixParametersNoReset()
Method which actually implements ConstrainChemicalPotentialsNoReset() (for backward compatibility).
virtual void ChangeTPS(ThermalParticleSystem *TPS)
Change the particle list.
virtual void SetGammaC(double gammaC)
Set the charm quark fugacity factor.
virtual void CalculateParticleChargeCorrelationMatrix()
Calculates the matrix of correlators between primordial (and also final) particle numbers and conserv...
virtual void SetTemperature(double T)
Set the temperature.
virtual void CalculateFluctuations()
Computes the fluctuation observables.
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 double TwoParticleSusceptibilityPrimordialByPdg(long long id1, long long id2)
Returns the computed primordial particle number (cross-)susceptibility for particles with pdg codes ...
virtual void SetVirial(int, int, double)
Set the excluded volume coefficient .
virtual void FillChemicalPotentials()
Sets the chemical potentials of all particles.
virtual void CalculateTemperatureDerivatives()
Computes the temperature derivatives of densities, shifted chemical potentials, and primordial hadron...
void ConstrainChemicalPotentials(bool resetInitialValues=true)
Constrains the chemical potentials by the conservation laws imposed.
virtual void FixParameters()
Method which actually implements ConstrainChemicalPotentials() (for backward compatibility).
virtual double TwoParticleSusceptibilityFinalByPdg(long long id1, long long id2)
Returns the computed final particle number (cross-)susceptibility for particles with pdg codes id1 a...
virtual double TwoParticleSusceptibilityFinal(int i, int j) const
Returns the computed final particle number (cross-)susceptibility for particles with ids i and j....
virtual void SetGammaS(double gammaS)
Set the strange quark fugacity factor.
virtual void SetCharm(int C)
Set the total charm (for canonical ensemble only)
virtual double FinalNetParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed final (after decays) net particle vs conserved charge (cross-)susceptibility fo...
void SetResonanceWidthIntegrationType(ThermalParticle::ResonanceWidthIntegration type)
Set the ThermalParticle::ResonanceWidthIntegration scheme for all particles. Calls the corresponding ...
virtual void CalculateFeeddown()
Calculates the total densities which include feeddown contributions.
virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordial(int i, int j) const
Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility ...
void SetUseWidth(bool useWidth)
Sets whether finite resonance widths are used. Deprecated.
virtual double FinalParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility for pa...
virtual void DisableBaryonAntiBaryonAttraction()
virtual double FinalParticleChargeSusceptibility(int i, ConservedCharge::Name chg) const
Returns the computed final (after decays) particle vs conserved charge (cross-)susceptibility for pa...
virtual void SetChemicalPotentials(const std::vector< double > &chem=std::vector< double >(0))
Sets the chemical potentials of all particles.
double ChemicalPotential(int i) const
Chemical potential of particle species i.
virtual void SetElectricCharge(int Q)
Set the total electric charge (for canonical ensemble only)
virtual std::vector< double > CalculateChargeFluctuations(const std::vector< double > &chgs, int order=4)
Calculates fluctuations (diagonal susceptibilities) of an arbitrary "conserved" charge.
virtual void SetBaryonChemicalPotential(double muB)
Set the baryon chemical potential.
virtual void CalculateTwoParticleCorrelations()
Computes the fluctuations and correlations of the primordial particle numbers.
virtual double SusceptibilityDimensionfull(ConservedCharge::Name i, ConservedCharge::Name j) const
A 2nd order susceptibility of conserved charges.
void SetNormBratio(bool normBratio)
Whether branching ratios are renormalized to 100%.
virtual void SetParameters(const ThermalModelParameters &params)
The thermal parameters.
virtual void SetElectricChemicalPotential(double muQ)
Set the electric chemical potential.
virtual double PrimordialNetParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed primordial net particle vs conserved charge (cross-)susceptibility for particle...
virtual double TwoParticleSusceptibilityTemperatureDerivativePrimordialByPdg(long long id1, long long id2)
Returns the computed temperature derivative of the primordial particle number (cross-)susceptibility ...
virtual double NetParticleSusceptibilityFinalByPdg(long long id1, long long id2)
Returns the computed (cross-)susceptibility between final net-particle numbers for pdg codes id1 and...
virtual void SetChemicalPotential(int i, double chem)
Sets the chemical potential of particle species i.
virtual bool FixChemicalPotentialsThroughDensities(double rhoB=0., double rhoQ=0., double rhoS=0., double rhoC=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 baryon,...
virtual void ValidateCalculation()
Checks whether issues have occured during the calculation of particle densities in the CalculateDensi...
int ComponentsNumber() const
Number of different particle species in the list.
virtual void SetBaryonCharge(int B)
Set the total baryon number (for canonical ensemble only)
virtual void CalculateProxySusceptibilityMatrix()
Calculates the susceptibility matrix of conserved charges proxies.
ThermalModelBase(ThermalParticleSystem *TPS, const ThermalModelParameters &params=ThermalModelParameters())
Construct a new ThermalModelBase object.
virtual double PrimordialParticleChargeSusceptibility(int i, ConservedCharge::Name chg) const
Returns the computed primordial particle vs conserved charge (cross-)susceptibility for particle wit...
virtual void CalculatePrimordialDensities()=0
Calculates the primordial densities of all species.
virtual void SetStrangenessChemicalPotential(double muS)
Set the strangeness chemical potential.
virtual double TwoParticleSusceptibilityPrimordial(int i, int j) const
Returns the computed primordial particle number (cross-)susceptibility for particles with ids i and ...
virtual void CalculateDensities()
Calculates the primordial and total (after decays) densities of all species.
virtual void FillVirial(const std::vector< double > &ri=std::vector< double >(0))
Fills the excluded volume coefficients based on the provided radii parameters for all species.
virtual void SetGammaq(double gammaq)
Set the light quark fugacity factor.
virtual void SetStrangeness(int S)
Set the total strangeness (for canonical ensemble only)
virtual void SetAttraction(int, int, double)
Set the vdW mean field attraction coefficient .
virtual double PrimordialParticleChargeSusceptibilityByPdg(long long id1, ConservedCharge::Name chg)
Returns the computed primordial particle vs conserved charge (cross-)susceptibility for particle wit...
const ThermalModelParameters & Parameters() const
virtual double NetParticleSusceptibilityPrimordialByPdg(long long id1, long long id2)
Returns the computed (cross-)susceptibility between primordial net-particle numbers for pdg codes id...
virtual void SetCharmChemicalPotential(double muC)
Set the charm chemical potential.
virtual void SetStatistics(bool stats)
double Volume() const
System volume (fm )
@ GCE
Grand canonical ensemble.
virtual void CalculateTwoParticleFluctuationsDecays()
Computes particle number correlations and fluctuations for all final-state particles which are marked...
Class containing all information about a particle specie.
double AbsoluteCharm() const
Absolute charm quark content |s|.
int BaryonCharge() const
Particle's baryon number.
double Density(const ThermalModelParameters &params, IdealGasFunctions::Quantity type=IdealGasFunctions::ParticleDensity, bool useWidth=0, double mu=0.) const
Computes a specified ideal gas thermodynamic function.
const ParticleDecaysVector & Decays() const
A vector of particle's decays.
ResonanceWidthIntegration
Treatment of finite resonance widths.
@ eBW
Energy-dependent Breit-Wigner scheme (eBW)
@ BWTwoGamma
Energy-independent Breit-Wigner in +-2\Gamma interval.
@ ZeroWidth
Zero-width approximation.
double AbsoluteStrangeness() const
Absolute strange quark content |s|.
double AbsoluteQuark() const
Absolute light quark content |u,d|.
const ParticleDecaysVector & DecaysOriginal() const
A backup copy of particle's decays.
Class containing the particle list.
std::vector< std::pair< double, std::vector< int > > > ResonanceFinalStatesDistribution
std::vector< SingleDecayContribution > DecayContributionsToParticle
A vector of SingleDecayContribution where each element corresponds to a certain resonance species.
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
MatrixXd GetSusceptibilityMatrix(ThermalModelBase *model, const vector< int > &ConservedDensities)
MatrixXd GetPressureHessianMatrix(ThermalModelBase *model, const vector< int > &ConservedDensities)
Name
Set of all conserved charges considered.
@ ElectricCharge
Electric charge.
@ Strong
Feeddown from strong decays.
@ Weak
Feeddown from strong, electromagnetic, and weak decays.
@ StabilityFlag
Feeddown from all particles marked as unstable.
@ Primordial
No feeddown.
static const int NumberOfDecayTypes
Structure containing all thermal parameters of the model.