PYTHIA  8.314
FragmentationFlavZpT.h
1 // FragmentationFlavZpT.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2025 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // This file contains helper classes for fragmentation.
7 // StringFlav is used to select quark and hadron flavours.
8 // StringPT is used to select transverse momenta.
9 // StringZ is used to sample the fragmentation function f(z).
10 
11 #ifndef Pythia8_FragmentationFlavZpT_H
12 #define Pythia8_FragmentationFlavZpT_H
13 
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"
20 
21 namespace Pythia8 {
22 
23 //==========================================================================
24 
25 // Functions for unnormalised, <z>, and RMSD(z) of Lund FF. The two latter
26 // return negative values in case of failure.
27 
28 double LundFFRaw(double z, double a, double b, double c, double mT2);
29 
30 double LundFFAvg(double a, double b, double mT2, double tol);
31 
32 double LundFFRms(double a, double b, double mT2, double tol);
33 
34 //==========================================================================
35 
36 // The FlavContainer class is a simple container for flavour,
37 // including the extra properties needed for popcorn baryon handling.
38 // id = current flavour.
39 // rank = current rank; 0 for endpoint flavour and then increase by 1.
40 // nPop = number of popcorn mesons yet to be produced (1 or 0).
41 // idPop = (absolute sign of) popcorn quark, shared between B and Bbar.
42 // idVtx = (absolute sign of) vertex (= non-shared) quark in diquark.
43 
45 
46 public:
47 
48  // Constructor.
49  FlavContainer(int idIn = 0, int rankIn = 0, int nPopIn = 0,
50  int idPopIn = 0, int idVtxIn = 0) : id(idIn), rank(rankIn),
51  nPop(nPopIn), idPop(idPopIn), idVtx(idVtxIn) {}
52 
53  // Copy constructor.
55  id = flav.id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
56  idVtx = flav.idVtx;}
57 
58  // Overloaded equal operator.
59  FlavContainer& operator=(const FlavContainer& flav) { if (this != &flav) {
60  id = flav.id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
61  idVtx = flav.idVtx; } return *this; }
62 
63  // Invert flavour.
64  FlavContainer& anti() {id = -id; return *this;}
65 
66  // Read in a container into another, without/with id sign flip.
67  FlavContainer& copy(const FlavContainer& flav) { if (this != &flav) {
68  id = flav.id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
69  idVtx = flav.idVtx; } return *this; }
70  FlavContainer& anti(const FlavContainer& flav) { if (this != &flav) {
71  id = -flav.id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
72  idVtx = flav.idVtx; } return *this; }
73 
74  // Check whether is diquark.
75  bool isDiquark() {int idAbs = abs(id);
76  return (idAbs > 1000 && idAbs < 10000 && (idAbs/10)%10 == 0);}
77 
78  // Stored properties.
79  int id, rank, nPop, idPop, idVtx;
80 
81 };
82 
83 //==========================================================================
84 
85 // The StringFlav class is used to select quark and hadron flavours.
86 
87 class StringFlav : public PhysicsBase {
88 
89 public:
90 
91  // Constructor.
93  suppressLeadingB(),
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) {}
108 
109  // Destructor.
110  virtual ~StringFlav() {}
111 
112  // Initialize data members.
113  virtual void init();
114 
115  // Initialise parameters when using close packing.
116  virtual void init(double kappaModifier, double strangeJunc,
117  double probQQmod);
118 
119  // Pick a light d, u or s quark according to fixed ratios.
120  int pickLightQ() { double rndmFlav = probQandS * rndmPtr->flat();
121  if (rndmFlav < 1.) return 1;
122  if (rndmFlav < 2.) return 2;
123  return 3; }
124 
125  // Pick a new flavour (including diquarks) given an incoming one,
126  // either by old standard Gaussian or new alternative exponential.
127  virtual FlavContainer pick(FlavContainer& flavOld, double pT = -1.0,
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); }
133  virtual FlavContainer pickGauss(FlavContainer& flavOld,
134  bool allowPop = true);
135  virtual FlavContainer pickThermal(FlavContainer& flavOld,
136  double pT, double kappaModifier);
137 
138  // Combine two flavours (including diquarks) to produce a hadron.
139  virtual int combine(FlavContainer& flav1, FlavContainer& flav2);
140 
141  // Ditto, simplified input argument for simple configurations.
142  virtual int combineId( int id1, int id2, bool keepTrying = true) {
143  FlavContainer flag1(id1); FlavContainer flag2(id2);
144  for (int i = 0; i < 100; ++i) { int idNew = combine( flag1, flag2);
145  if (idNew != 0 || !keepTrying) return idNew;} return 0;}
146 
147  // Combine three (di-) quark flavours into two hadrons.
148  virtual pair<int,int> combineDiquarkJunction(int id1, int id2, int id3);
149 
150  // Combine two flavours to produce a hadron with lowest possible mass.
151  virtual int combineToLightest( int id1, int id2);
152 
153  // Lightest flavour-neutral meson.
154  virtual int idLightestNeutralMeson() { return 111; }
155 
156  // Return chosen hadron in case of thermal model.
157  virtual int getHadronIDwin() { return hadronIDwin; }
158 
159  // Combine two flavours into hadron for last two remaining flavours
160  // for thermal model.
161  virtual int combineLastThermal(FlavContainer& flav1, FlavContainer& flav2,
162  double pT, double kappaModifier);
163 
164  // General function, decides whether to just return the hadron id
165  // if thermal model was use or whether to combine the two flavours.
166  virtual int getHadronID(FlavContainer& flav1, FlavContainer& flav2,
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); }
174 
175  // Return hadron mass. Used one if present, pick otherwise.
176  virtual double getHadronMassWin(int idHad) { return
177  ((hadronMassWin < 0.0) ? particleDataPtr->mSel(idHad) : hadronMassWin); }
178 
179  // Assign popcorn quark inside an original (= rank 0) diquark.
180  void assignPopQ(FlavContainer& flav);
181 
182  // Combine two quarks to produce a diquark.
183  int makeDiquark(int id1, int id2, int idHad = 0);
184 
185  // Check if quark-diquark combination should be added. If so add.
186  void addQuarkDiquark(vector< pair<int,int> >& quarkCombis,
187  int qID, int diqID, int hadronID) {
188  bool allowed = true;
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) ); }
194 
195  // Get spin counter for mesons.
196  int getMesonSpinCounter(int hadronID) { hadronID = abs(hadronID);
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;
201  return -1; }
202 
203  // Get the flavour and spin ratios calculated from the diquark weights.
204  // i: (0) q -> B B, (1) q -> B M B, (2) qq -> M B
205  // j: (0) s/u popcorn ratio, (1/2) s/u ratio for vertex quark if popcorn
206  // quark is u/d or s, (3) q/q' vertex quark ratio if popcorn quark is
207  // light and = q, (4/5/6) (spin 1)/(spin 0) ratio for su, us and ud
208  double getFlavourSpinRatios(int i, int j) {
209  return (i < 3 && j < 7) ? dWT[i][j] : -1.0;}
210 
211 protected:
212 
213  // Initialise derived parameters.
214  virtual void initDerived();
215 
216  // Constants: could only be changed in the code itself.
217  static const int mesonMultipletCode[6];
218  static const double baryonCGOct[6], baryonCGDec[6];
219 
220  // Settings for Gaussian model.
221  bool suppressLeadingB, mT2suppression, useWidthPre;
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,
228  heavyLeadingBSup;
229  bool qqKappa;
230  double probStoUDSav, probQQtoQSav, probSQtoQQSav, probQQ1toQQ0Sav,
231  alphaQQSav, sigmaHad, widthPreStrange, widthPreDiquark;
232 
233  // Settings for thermal model.
234  bool thermalModel, mesonNonetL1;
235  double temperature, tempPreFactor;
236  int nNewQuark;
237  double mesMixRate1[2][6], mesMixRate2[2][6], mesMixRate3[2][6];
238  double baryonOctWeight[6][6][6][2], baryonDecWeight[6][6][6][2];
239 
240  // Settings used by both models.
241  bool closePacking, doEnhanceDiquark;
242  double enhanceStrange, enhancePT, enhanceDiquark, exponentMPI, exponentNSP;
243 
244  // Key = hadron id, value = list of constituent ids.
245  map< int, vector< pair<int,int> > > hadronConstIDs;
246  // Key = initial (di)quark id, value = list of possible hadron ids
247  // + nr in hadronConstIDs.
248  map< int, vector< pair<int,int> > > possibleHadrons;
249  // Key = initial (di)quark id, value = prefactor to multiply rate.
250  map< int, vector<double> > possibleRatePrefacs;
251  // Similar, but for combining the last two (di)quarks. Key = (di)quark pair.
252  map< pair<int,int>, vector< pair<int,int> > > possibleHadronsLast;
253  map< pair<int,int>, vector<double> > possibleRatePrefacsLast;
254 
255  // Selection in thermal model.
256  int hadronIDwin, idNewWin;
257  double hadronMassWin;
258 
259  // Fragmentation weights container.
261 
262 };
263 
264 //==========================================================================
265 
266 // The StringZ class is used to sample the fragmentation function f(z).
267 
268 class StringZ : public PhysicsBase {
269 
270 public:
271 
272  // Constructor.
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() {}
278 
279  // Destructor.
280  virtual ~StringZ() {}
281 
282  // Initialize data members.
283  virtual void init();
284 
285  // Fragmentation function: top-level to determine parameters.
286  virtual double zFrag( int idOld, int idNew = 0, double mT2 = 1.);
287 
288  // Fragmentation function: select z according to provided parameters.
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.);
295 
296  // Parameters for stopping in the middle; overloaded for Hidden Valley.
297  virtual double stopMass() {return stopM;}
298  virtual double stopNewFlav() {return stopNF;}
299  virtual double stopSmear() {return stopS;}
300 
301  // a and b fragmentation parameters needed in some operations.
302  virtual double aAreaLund() {return aLund;}
303  virtual double bAreaLund() {return bLund;}
304 
305  // Method to derive both a and b parameters (from <z> and RMSD(z)).
306  bool deriveABLund( bool derivaA = false, bool deriveAExtraDiquark = false,
307  bool deriveAExtraSQuark = false);
308  // Method to derive bLund from <z> (for fixed a and reference mT2).
309  double deriveBLund( double avgZ, double a, double mT2ref);
310 
311  protected:
312 
313  // Constants: could only be changed in the code itself.
314  static const double CFROMUNITY, AFROMZERO, AFROMC, EXPMAX;
315 
316  // Initialization data, to be read from Settings.
317  bool useNonStandC, useNonStandB, useNonStandH,
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;
322 
323  // Fragmentation weights container.
325 
326 };
327 
328 //==========================================================================
329 
330 // The StringPT class is used to select select transverse momenta.
331 
332 class StringPT : public PhysicsBase {
333 
334 public:
335 
336  // Constructor.
337  StringPT() : useWidthPre(), sigmaQ(), enhancedFraction(), enhancedWidth(),
338  sigma2Had(), widthPreStrange(), widthPreDiquark(),
339  thermalModel(), temperature(), tempPreFactor(), fracSmallX(),
340  closePacking(), enhancePT(), exponentMPI(), exponentNSP() {}
341 
342  // Destructor.
343  virtual ~StringPT() {}
344 
345  // Initialize data members.
346  virtual void init();
347 
348  // General function, return px and py as a pair in the same call
349  // in either model.
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);
355 
356  // Gaussian suppression of given pT2; used in MiniStringFragmentation.
357  double suppressPT2(double pT2) { return (thermalModel ?
358  exp(-sqrt(pT2)/temperature) : exp(-pT2/sigma2Had)); }
359 
360 protected:
361 
362  // Constants: could only be changed in the code itself.
363  static const double SIGMAMIN;
364 
365  // Initialization data, to be read from Settings.
366  // Gaussian model.
368  double sigmaQ, enhancedFraction, enhancedWidth, sigma2Had,
369  widthPreStrange, widthPreDiquark;
370  // Thermal model.
372  double temperature, tempPreFactor, fracSmallX;
373  // Both.
375  double enhancePT, exponentMPI, exponentNSP;
376 
377 private:
378 
379  // Evaluate Bessel function K_{1/4}(x).
380  double BesselK14(double x);
381 
382  // Fragmentation weights container.
383  WeightsFragmentation* wgtsPtr{};
384 
385 };
386 
387 //==========================================================================
388 
389 } // end namespace Pythia8
390 
391 #endif // Pythia8_FragmentationFlavZpT_H
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