File indexing completed on 2026-04-18 08:20:44
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012 #ifndef Pythia8_ThermalFragmentation_H
0013 #define Pythia8_ThermalFragmentation_H
0014
0015 #include "Pythia8/FragmentationModel.h"
0016 #include "Pythia8/MiniStringFragmentation.h"
0017 #include "Pythia8/StringFragmentation.h"
0018
0019 namespace Pythia8 {
0020
0021
0022
0023
0024
0025 class ThermalStringFlav : public StringFlav {
0026
0027 public:
0028
0029
0030 ThermalStringFlav() : mesonNonetL1(), temperature(), tempPreFactor(),
0031 nNewQuark(), mesMixRate1(), mesMixRate2(), mesMixRate3(),
0032 baryonOctWeight(), baryonDecWeight(), hadronConstIDs(), possibleHadrons(),
0033 possibleRatePrefacs(), possibleHadronsLast(), possibleRatePrefacsLast(),
0034 hadronIDwin(), quarkIDwin(), hadronMassWin() {}
0035
0036
0037 ~ThermalStringFlav() {}
0038
0039
0040 void init() override;
0041
0042
0043 FlavContainer pick(FlavContainer& flavOld,
0044 double pT, double kappaModifier, bool allowPop) override;
0045
0046
0047 virtual int getHadronIDwin() { return hadronIDwin; }
0048
0049
0050 double getHadronMassWin(int idHad) override { return
0051 ((hadronMassWin < 0.0) ? particleDataPtr->mSel(idHad) : hadronMassWin); }
0052
0053
0054
0055 int combineLastThermal(FlavContainer& flav1, FlavContainer& flav2,
0056 double pT, double kappaModifier);
0057
0058
0059 int getHadronID(FlavContainer& flav1, FlavContainer& flav2, double pT = -1.0,
0060 double kappaModifier = -1.0, bool finalTwo = false) override {
0061 if (finalTwo) return combineLastThermal(flav1, flav2, pT, kappaModifier);
0062 if ( (hadronIDwin != 0) && (quarkIDwin != 0)) return getHadronIDwin();
0063 return combine(flav1, flav2); }
0064
0065
0066 void addQuarkDiquark(vector< pair<int,int> >& quarkCombis,
0067 int qID, int diqID, int hadronID) {
0068 bool allowed = true;
0069 for (int iCombi = 0; iCombi < int(quarkCombis.size()); iCombi++)
0070 if ( (qID == quarkCombis[iCombi].first ) &&
0071 (diqID == quarkCombis[iCombi].second) ) allowed = false;
0072 if (allowed) quarkCombis.push_back( (hadronID > 0) ?
0073 make_pair( qID, diqID) : make_pair(-qID, -diqID) ); }
0074
0075 protected:
0076
0077
0078 bool mesonNonetL1;
0079 double temperature, tempPreFactor;
0080 int nNewQuark;
0081 double mesMixRate1[2][6], mesMixRate2[2][6], mesMixRate3[2][6];
0082 double baryonOctWeight[6][6][6][2], baryonDecWeight[6][6][6][2];
0083
0084
0085 map< int, vector< pair<int,int> > > hadronConstIDs;
0086
0087
0088 map< int, vector< pair<int,int> > > possibleHadrons;
0089
0090 map< int, vector<double> > possibleRatePrefacs;
0091
0092 map< pair<int,int>, vector< pair<int,int> > > possibleHadronsLast;
0093 map< pair<int,int>, vector<double> > possibleRatePrefacsLast;
0094
0095
0096 int hadronIDwin, quarkIDwin;
0097 double hadronMassWin;
0098
0099 };
0100
0101
0102
0103
0104
0105 class ThermalStringPT : public StringPT {
0106
0107 public:
0108
0109
0110 ThermalStringPT() : temperature(), tempPreFactor(), fracSmallX() {}
0111
0112
0113 ~ThermalStringPT() {}
0114
0115
0116 void init() override;
0117
0118
0119 pair<double, double> pxy(int idIn = 0, double kappaModifier = -1.0) override;
0120
0121
0122 double suppressPT2(double pT2) override {
0123 return exp(-sqrt(pT2)/temperature); }
0124
0125 protected:
0126
0127
0128 double temperature, tempPreFactor, fracSmallX;
0129
0130 private:
0131
0132
0133 double BesselK14(double x);
0134
0135 };
0136
0137
0138
0139
0140
0141
0142 class ThermalFragmentation : public FragmentationModel {
0143
0144 public:
0145
0146
0147 ThermalFragmentation();
0148
0149
0150 ~ThermalFragmentation() override;
0151
0152
0153 bool init(StringFlav* flavSelPtrIn = nullptr,
0154 StringPT* pTSelPtrIn = nullptr, StringZ* zSelPtrIn = nullptr,
0155 FragModPtr fragModPtrIn = nullptr) override;
0156
0157
0158 bool fragment(int iSub, ColConfig& colConfig, Event& event,
0159 bool isDiff = false, bool systemRecoil = true) override;
0160
0161
0162 ThermalStringFlav* thermalFlavSelPtr{};
0163 ThermalStringPT* thermalPTSelPtr{};
0164 StringZ* zSelPtr{};
0165
0166
0167 StringFragmentation* stringFragPtr{};
0168 MiniStringFragmentation* ministringFragPtr{};
0169
0170 private:
0171
0172
0173 double mStringMin{};
0174 bool tryMiniAfterFailedFrag{};
0175
0176 };
0177
0178
0179
0180 }
0181
0182 #endif