11 #ifndef Pythia8_FragmentationFlavZpT_H 12 #define Pythia8_FragmentationFlavZpT_H 14 #include "Pythia8/Basics.h" 15 #include "Pythia8/MathTools.h" 16 #include "Pythia8/ParticleData.h" 17 #include "Pythia8/PhysicsBase.h" 18 #include "Pythia8/PythiaStdlib.h" 19 #include "Pythia8/Settings.h" 28 double LundFFRaw(
double z,
double a,
double b,
double c,
double mT2);
30 double LundFFAvg(
double a,
double b,
double mT2,
double tol);
32 double LundFFRms(
double a,
double b,
double mT2,
double tol);
50 int idPopIn = 0,
int idVtxIn = 0) :
id(idIn), rank(rankIn),
51 nPop(nPopIn), idPop(idPopIn), idVtx(idVtxIn) {}
55 id = flav.
id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
60 id = flav.
id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
61 idVtx = flav.idVtx; }
return *
this; }
68 id = flav.
id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
69 idVtx = flav.idVtx; }
return *
this; }
71 id = -flav.
id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
72 idVtx = flav.idVtx; }
return *
this; }
76 return (idAbs > 1000 && idAbs < 10000 && (idAbs/10)%10 == 0);}
79 int id, rank, nPop, idPop, idVtx;
94 mT2suppression(), useWidthPre(), probQQtoQ(), probStoUD(), probSQtoQQ(),
95 probQQ1toQQ0(), probQandQQ(), probQandS(), probQandSinQQ(), probQQ1corr(),
96 probQQ1corrInv(), probQQ1norm(), probQQ1join(), mesonRate(),
97 mesonRateSum(), mesonMix1(), mesonMix2(), etaSup(), etaPrimeSup(),
98 decupletSup(), baryonCGSum(), baryonCGMax(), popcornRate(), popcornSpair(),
99 popcornSmeson(), barCGMax(), scbBM(), popFrac(), popS(), dWT(),
100 lightLeadingBSup(), heavyLeadingBSup(), probStoUDSav(), probQQtoQSav(),
101 probSQtoQQSav(), probQQ1toQQ0Sav(), alphaQQSav(), sigmaHad(),
102 widthPreStrange(), widthPreDiquark(), thermalModel(), mesonNonetL1(),
103 temperature(), tempPreFactor(), nNewQuark(), mesMixRate1(), mesMixRate2(),
104 mesMixRate3(), baryonOctWeight(), baryonDecWeight(), closePacking(),
105 doEnhanceDiquark(), enhanceStrange(), enhancePT(), enhanceDiquark(),
106 exponentMPI(), exponentNSP(), hadronIDwin(0), idNewWin(0),
107 hadronMassWin(-1.0) {}
116 virtual void init(
double kappaModifier,
double strangeJunc,
120 int pickLightQ() {
double rndmFlav = probQandS * rndmPtr->flat();
121 if (rndmFlav < 1.)
return 1;
122 if (rndmFlav < 2.)
return 2;
128 double kappaModifier = -1.0,
bool allowPop =
true) {
129 hadronIDwin = 0; idNewWin = 0; hadronMassWin = -1.0;
130 if ( (thermalModel || mT2suppression) && (pT >= 0.0) )
131 return pickThermal(flavOld, pT, kappaModifier);
132 return pickGauss(flavOld, allowPop); }
134 bool allowPop =
true);
136 double pT,
double kappaModifier);
142 virtual int combineId(
int id1,
int id2,
bool keepTrying =
true) {
144 for (
int i = 0; i < 100; ++i) {
int idNew = combine( flag1, flag2);
145 if (idNew != 0 || !keepTrying)
return idNew;}
return 0;}
148 virtual pair<int,int> combineDiquarkJunction(
int id1,
int id2,
int id3);
151 virtual int combineToLightest(
int id1,
int id2);
162 double pT,
double kappaModifier);
167 double pT = -1.0,
double kappaModifier = -1.0,
bool finalTwo =
false) {
168 if (finalTwo)
return ((thermalModel || mT2suppression) ?
169 combineLastThermal(flav1, flav2, pT, kappaModifier)
170 : combine(flav1, flav2));
171 if ((thermalModel || mT2suppression)&& (hadronIDwin != 0)
172 && (idNewWin != 0))
return getHadronIDwin();
173 return combine(flav1, flav2); }
177 ((hadronMassWin < 0.0) ? particleDataPtr->mSel(idHad) : hadronMassWin); }
183 int makeDiquark(
int id1,
int id2,
int idHad = 0);
187 int qID,
int diqID,
int hadronID) {
189 for (
int iCombi = 0; iCombi < int(quarkCombis.size()); iCombi++)
190 if ( (qID == quarkCombis[iCombi].first ) &&
191 (diqID == quarkCombis[iCombi].second) ) allowed =
false;
192 if (allowed) quarkCombis.push_back( (hadronID > 0) ?
193 make_pair( qID, diqID) : make_pair(-qID, -diqID) ); }
197 int j = (hadronID % 10);
198 if (hadronID < 1000)
return ((j==1) ? 0 : ( (j==3) ? 1 : 5 ));
199 if (hadronID < 20000)
return ((j==1) ? 3 : 2);
200 if (hadronID > 20000)
return 4;
209 return (i < 3 && j < 7) ? dWT[i][j] : -1.0;}
214 virtual void initDerived();
217 static const int mesonMultipletCode[6];
218 static const double baryonCGOct[6], baryonCGDec[6];
222 double probQQtoQ, probStoUD, probSQtoQQ, probQQ1toQQ0, probQandQQ,
223 probQandS, probQandSinQQ, probQQ1corr, probQQ1corrInv, probQQ1norm,
224 probQQ1join[4], mesonRate[4][6], mesonRateSum[4], mesonMix1[2][6],
225 mesonMix2[2][6], etaSup, etaPrimeSup, decupletSup, baryonCGSum[6],
226 baryonCGMax[6], popcornRate, popcornSpair, popcornSmeson, barCGMax[8],
227 scbBM[3], popFrac, popS[3], dWT[3][7], lightLeadingBSup,
230 double probStoUDSav, probQQtoQSav, probSQtoQQSav, probQQ1toQQ0Sav,
231 alphaQQSav, sigmaHad, widthPreStrange, widthPreDiquark;
235 double temperature, tempPreFactor;
237 double mesMixRate1[2][6], mesMixRate2[2][6], mesMixRate3[2][6];
238 double baryonOctWeight[6][6][6][2], baryonDecWeight[6][6][6][2];
242 double enhanceStrange, enhancePT, enhanceDiquark, exponentMPI, exponentNSP;
253 map< pair<int,int>, vector<double> > possibleRatePrefacsLast;
257 double hadronMassWin;
273 StringZ() : useNonStandC(), useNonStandB(), useNonStandH(), usePetersonC(),
274 usePetersonB(), usePetersonH(), mc2(), mb2(), aLund(), bLund(),
275 aExtraSQuark(), aExtraDiquark(), rFactC(), rFactB(), rFactH(), aNonC(),
276 aNonB(), aNonH(), bNonC(), bNonB(), bNonH(), epsilonC(), epsilonB(),
277 epsilonH(), stopM(), stopNF(), stopS() {}
286 virtual double zFrag(
int idOld,
int idNew = 0,
double mT2 = 1.);
289 virtual double zLund(
double a,
double b,
double c = 1.,
290 double head = 1.,
double bNow = 0.,
int idFrag = 0,
291 bool isOldSQuark =
false,
bool isNewSQuark =
false,
292 bool isOldDiquark =
false,
bool isNewDiquark =
false);
293 virtual double zPeterson(
double epsilon);
294 virtual double zLundMax(
double a,
double b,
double c = 1.);
298 virtual double stopNewFlav() {
return stopNF;}
299 virtual double stopSmear() {
return stopS;}
303 virtual double bAreaLund() {
return bLund;}
306 bool deriveABLund(
bool derivaA =
false,
bool deriveAExtraDiquark =
false,
307 bool deriveAExtraSQuark =
false);
309 double deriveBLund(
double avgZ,
double a,
double mT2ref);
314 static const double CFROMUNITY, AFROMZERO, AFROMC,
EXPMAX;
318 usePetersonC, usePetersonB, usePetersonH;
319 double mc2, mb2, aLund, bLund, aExtraSQuark, aExtraDiquark, rFactC,
320 rFactB, rFactH, aNonC, aNonB, aNonH, bNonC, bNonB, bNonH,
321 epsilonC, epsilonB, epsilonH, stopM, stopNF, stopS;
337 StringPT() : useWidthPre(), sigmaQ(), enhancedFraction(), enhancedWidth(),
338 sigma2Had(), widthPreStrange(), widthPreDiquark(),
339 thermalModel(), temperature(), tempPreFactor(), fracSmallX(),
340 closePacking(), enhancePT(), exponentMPI(), exponentNSP() {}
350 pair<double, double>
pxy(
int idIn,
double kappaModifier = -1.0) {
351 return (thermalModel ? pxyThermal(idIn, kappaModifier) :
352 pxyGauss(idIn, kappaModifier)); }
353 pair<double, double> pxyGauss(
int idIn = 0,
double kappaModifier = -1.0);
354 pair<double, double> pxyThermal(
int idIn,
double kappaModifier = -1.0);
358 exp(-sqrt(pT2)/temperature) : exp(-pT2/sigma2Had)); }
368 double sigmaQ, enhancedFraction, enhancedWidth, sigma2Had,
369 widthPreStrange, widthPreDiquark;
372 double temperature, tempPreFactor, fracSmallX;
375 double enhancePT, exponentMPI, exponentNSP;
380 double BesselK14(
double x);
int hadronIDwin
Selection in thermal model.
Definition: FragmentationFlavZpT.h:256
pair< double, double > pxy(int idIn, double kappaModifier=-1.0)
Definition: FragmentationFlavZpT.h:350
double LundFFRms(double a, double b, double mT2, double tol)
Definition: FragmentationFlavZpT.cc:70
FlavContainer & operator=(const FlavContainer &flav)
Overloaded equal operator.
Definition: FragmentationFlavZpT.h:59
virtual double aAreaLund()
a and b fragmentation parameters needed in some operations.
Definition: FragmentationFlavZpT.h:302
Definition: PhysicsBase.h:27
virtual double getHadronMassWin(int idHad)
Return hadron mass. Used one if present, pick otherwise.
Definition: FragmentationFlavZpT.h:176
int getMesonSpinCounter(int hadronID)
Get spin counter for mesons.
Definition: FragmentationFlavZpT.h:196
map< int, vector< pair< int, int > > > possibleHadrons
Definition: FragmentationFlavZpT.h:248
static const double SIGMAMIN
Constants: could only be changed in the code itself.
Definition: FragmentationFlavZpT.h:363
bool thermalModel
Settings for thermal model.
Definition: FragmentationFlavZpT.h:234
The StringPT class is used to select select transverse momenta.
Definition: FragmentationFlavZpT.h:332
virtual ~StringPT()
Destructor.
Definition: FragmentationFlavZpT.h:343
double LundFFAvg(double a, double b, double mT2, double tol)
Definition: FragmentationFlavZpT.cc:36
bool thermalModel
Thermal model.
Definition: FragmentationFlavZpT.h:371
bool isDiquark()
Check whether is diquark.
Definition: FragmentationFlavZpT.h:75
Wave4 epsilon(Wave4 w1, Wave4 w2, Wave4 w3)
Permutation operator.
Definition: HelicityBasics.cc:56
virtual ~StringZ()
Destructor.
Definition: FragmentationFlavZpT.h:280
virtual ~StringFlav()
Destructor.
Definition: FragmentationFlavZpT.h:110
void addQuarkDiquark(vector< pair< int, int > > &quarkCombis, int qID, int diqID, int hadronID)
Check if quark-diquark combination should be added. If so add.
Definition: FragmentationFlavZpT.h:186
virtual int getHadronID(FlavContainer &flav1, FlavContainer &flav2, double pT=-1.0, double kappaModifier=-1.0, bool finalTwo=false)
Definition: FragmentationFlavZpT.h:166
The StringZ class is used to sample the fragmentation function f(z).
Definition: FragmentationFlavZpT.h:268
FlavContainer & anti()
Invert flavour.
Definition: FragmentationFlavZpT.h:64
FlavContainer(int idIn=0, int rankIn=0, int nPopIn=0, int idPopIn=0, int idVtxIn=0)
Constructor.
Definition: FragmentationFlavZpT.h:49
int id
Stored properties.
Definition: FragmentationFlavZpT.h:79
double getFlavourSpinRatios(int i, int j)
Definition: FragmentationFlavZpT.h:208
bool suppressLeadingB
Settings for Gaussian model.
Definition: FragmentationFlavZpT.h:221
bool useWidthPre
Definition: FragmentationFlavZpT.h:367
virtual int combineId(int id1, int id2, bool keepTrying=true)
Ditto, simplified input argument for simple configurations.
Definition: FragmentationFlavZpT.h:142
double suppressPT2(double pT2)
Gaussian suppression of given pT2; used in MiniStringFragmentation.
Definition: FragmentationFlavZpT.h:357
FlavContainer & copy(const FlavContainer &flav)
Read in a container into another, without/with id sign flip.
Definition: FragmentationFlavZpT.h:67
virtual int idLightestNeutralMeson()
Lightest flavour-neutral meson.
Definition: FragmentationFlavZpT.h:154
FlavContainer(const FlavContainer &flav)
Copy constructor.
Definition: FragmentationFlavZpT.h:54
double LundFFRaw(double z, double a, double b, double c, double mT2)
Functions for unnormalised and average Lund FF.
Definition: FragmentationFlavZpT.cc:21
StringZ()
Constructor.
Definition: FragmentationFlavZpT.h:273
map< int, vector< pair< int, int > > > hadronConstIDs
Key = hadron id, value = list of constituent ids.
Definition: FragmentationFlavZpT.h:245
virtual double stopMass()
Parameters for stopping in the middle; overloaded for Hidden Valley.
Definition: FragmentationFlavZpT.h:297
Definition: FragmentationFlavZpT.h:44
Definition: Weights.h:354
map< pair< int, int >, vector< pair< int, int > > > possibleHadronsLast
Similar, but for combining the last two (di)quarks. Key = (di)quark pair.
Definition: FragmentationFlavZpT.h:252
static const double EXPMAX
Do not take exponent of too large or small number.
Definition: FragmentationFlavZpT.h:314
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
bool closePacking
Both.
Definition: FragmentationFlavZpT.h:374
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:87
virtual FlavContainer pick(FlavContainer &flavOld, double pT=-1.0, double kappaModifier=-1.0, bool allowPop=true)
Definition: FragmentationFlavZpT.h:127
bool closePacking
Settings used by both models.
Definition: FragmentationFlavZpT.h:241
bool useNonStandC
Initialization data, to be read from Settings.
Definition: FragmentationFlavZpT.h:317
int pickLightQ()
Pick a light d, u or s quark according to fixed ratios.
Definition: FragmentationFlavZpT.h:120
StringPT()
Constructor.
Definition: FragmentationFlavZpT.h:337
StringFlav()
Constructor.
Definition: FragmentationFlavZpT.h:92
map< int, vector< double > > possibleRatePrefacs
Key = initial (di)quark id, value = prefactor to multiply rate.
Definition: FragmentationFlavZpT.h:250
virtual int getHadronIDwin()
Return chosen hadron in case of thermal model.
Definition: FragmentationFlavZpT.h:157