StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
FragmentationFlavZpT.h
1 // FragmentationFlavZpT.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 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 and average Lund FF.
26 
27 double LundFFRaw(double z, double a, double b, double c, double mT2);
28 
29 double LundFFAvg(double a, double b, double c, double mT2, double tol);
30 
31 //==========================================================================
32 
33 // The FlavContainer class is a simple container for flavour,
34 // including the extra properties needed for popcorn baryon handling.
35 // id = current flavour.
36 // rank = current rank; 0 for endpoint flavour and then increase by 1.
37 // nPop = number of popcorn mesons yet to be produced (1 or 0).
38 // idPop = (absolute sign of) popcorn quark, shared between B and Bbar.
39 // idVtx = (absolute sign of) vertex (= non-shared) quark in diquark.
40 
41 class FlavContainer {
42 
43 public:
44 
45  // Constructor.
46  FlavContainer(int idIn = 0, int rankIn = 0, int nPopIn = 0,
47  int idPopIn = 0, int idVtxIn = 0) : id(idIn), rank(rankIn),
48  nPop(nPopIn), idPop(idPopIn), idVtx(idVtxIn) {}
49 
50  // Copy constructor.
51  FlavContainer(const FlavContainer& flav) {
52  id = flav.id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
53  idVtx = flav.idVtx;}
54 
55  // Overloaded equal operator.
56  FlavContainer& operator=(const FlavContainer& flav) { if (this != &flav) {
57  id = flav.id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
58  idVtx = flav.idVtx; } return *this; }
59 
60  // Invert flavour.
61  FlavContainer& anti() {id = -id; return *this;}
62 
63  // Read in a container into another, without/with id sign flip.
64  FlavContainer& copy(const FlavContainer& flav) { if (this != &flav) {
65  id = flav.id; rank = flav.rank; nPop = flav.nPop; idPop = flav.idPop;
66  idVtx = flav.idVtx; } return *this; }
67  FlavContainer& anti(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 
71  // Check whether is diquark.
72  bool isDiquark() {int idAbs = abs(id);
73  return (idAbs > 1000 && idAbs < 10000 && (idAbs/10)%10 == 0);}
74 
75  // Stored properties.
76  int id, rank, nPop, idPop, idVtx;
77 
78 };
79 
80 //==========================================================================
81 
82 // The StringFlav class is used to select quark and hadron flavours.
83 
84 class StringFlav : public PhysicsBase {
85 
86 public:
87 
88  // Constructor.
89  StringFlav() :
90  suppressLeadingB(),
91  mT2suppression(), useWidthPre(), probQQtoQ(), probStoUD(), probSQtoQQ(),
92  probQQ1toQQ0(), probQandQQ(), probQandS(), probQandSinQQ(), probQQ1corr(),
93  probQQ1corrInv(), probQQ1norm(), probQQ1join(), mesonRate(),
94  mesonRateSum(), mesonMix1(), mesonMix2(), etaSup(), etaPrimeSup(),
95  decupletSup(), baryonCGSum(), baryonCGMax(), popcornRate(), popcornSpair(),
96  popcornSmeson(), scbBM(), popFrac(), popS(), dWT(), lightLeadingBSup(),
97  heavyLeadingBSup(), sigmaHad(), widthPreStrange(), widthPreDiquark(),
98  thermalModel(), mesonNonetL1(), temperature(), tempPreFactor(),
99  nNewQuark(), mesMixRate1(), mesMixRate2(), mesMixRate3(),
100  baryonOctWeight(), baryonDecWeight(), closePacking(), exponentMPI(),
101  exponentNSP(), hadronIDwin(0), idNewWin(0), hadronMassWin(-1.0) {}
102 
103  // Destructor.
104  virtual ~StringFlav() {}
105 
106  // Initialize data members.
107  virtual void init();
108 
109  // Pick a light d, u or s quark according to fixed ratios.
110  int pickLightQ() { double rndmFlav = probQandS * rndmPtr->flat();
111  if (rndmFlav < 1.) return 1;
112  if (rndmFlav < 2.) return 2;
113  return 3; }
114 
115  // Pick a new flavour (including diquarks) given an incoming one,
116  // either by old standard Gaussian or new alternative exponential.
117  virtual FlavContainer pick(FlavContainer& flavOld, double pT = -1.0,
118  double nNSP = 0.0, bool allowPop = true) {
119  hadronIDwin = 0; idNewWin = 0; hadronMassWin = -1.0;
120  if ( (thermalModel || mT2suppression) && (pT >= 0.0) )
121  return pickThermal(flavOld, pT, nNSP);
122  return pickGauss(flavOld, allowPop); }
123  virtual FlavContainer pickGauss(FlavContainer& flavOld,
124  bool allowPop = true);
125  virtual FlavContainer pickThermal(FlavContainer& flavOld,
126  double pT, double nNSP);
127 
128  // Combine two flavours (including diquarks) to produce a hadron.
129  virtual int combine(FlavContainer& flav1, FlavContainer& flav2);
130 
131  // Ditto, simplified input argument for simple configurations.
132  virtual int combineId( int id1, int id2, bool keepTrying = true) {
133  FlavContainer flag1(id1); FlavContainer flag2(id2);
134  for (int i = 0; i < 100; ++i) { int idNew = combine( flag1, flag2);
135  if (idNew != 0 || !keepTrying) return idNew;} return 0;}
136 
137  // Combine two flavours to produce a hadron with lowest possible mass.
138  virtual int combineToLightest( int id1, int id2);
139 
140  // Return chosen hadron in case of thermal model.
141  virtual int getHadronIDwin() { return hadronIDwin; }
142 
143  // Combine two flavours into hadron for last two remaining flavours
144  // for thermal model.
145  virtual int combineLastThermal(FlavContainer& flav1, FlavContainer& flav2,
146  double pT, double nNSP);
147 
148  // General function, decides whether to just return the hadron id
149  // if thermal model was use or whether to combine the two flavours.
150  virtual int getHadronID(FlavContainer& flav1, FlavContainer& flav2,
151  double pT = -1.0, double nNSP = 0, bool finalTwo = false) {
152  if (finalTwo) return ((thermalModel || mT2suppression) ?
153  combineLastThermal(flav1, flav2, pT, nNSP) : combine(flav1, flav2));
154  if ((thermalModel || mT2suppression)&& (hadronIDwin != 0)
155  && (idNewWin != 0)) return getHadronIDwin();
156  return combine(flav1, flav2); }
157 
158  // Return hadron mass. Used one if present, pick otherwise.
159  virtual double getHadronMassWin(int idHad) { return
160  ((hadronMassWin < 0.0) ? particleDataPtr->mSel(idHad) : hadronMassWin); }
161 
162  // Assign popcorn quark inside an original (= rank 0) diquark.
163  void assignPopQ(FlavContainer& flav);
164 
165  // Combine two quarks to produce a diquark.
166  int makeDiquark(int id1, int id2, int idHad = 0);
167 
168  // Check if quark-diquark combination should be added. If so add.
169  void addQuarkDiquark(vector< pair<int,int> >& quarkCombis,
170  int qID, int diqID, int hadronID) {
171  bool allowed = true;
172  for (int iCombi = 0; iCombi < int(quarkCombis.size()); iCombi++)
173  if ( (qID == quarkCombis[iCombi].first ) &&
174  (diqID == quarkCombis[iCombi].second) ) allowed = false;
175  if (allowed) quarkCombis.push_back( (hadronID > 0) ?
176  make_pair( qID, diqID) : make_pair(-qID, -diqID) ); }
177 
178  // Get spin counter for mesons.
179  int getMesonSpinCounter(int hadronID) { hadronID = abs(hadronID);
180  int j = (hadronID % 10);
181  if (hadronID < 1000) return ((j==1) ? 0 : ( (j==3) ? 1 : 5 ));
182  if (hadronID < 20000) return ((j==1) ? 3 : 2);
183  if (hadronID > 20000) return 4;
184  return -1; }
185 
186 protected:
187 
188  // Constants: could only be changed in the code itself.
189  static const int mesonMultipletCode[6];
190  static const double baryonCGOct[6], baryonCGDec[6];
191 
192  // Settings for Gaussian model.
193  bool suppressLeadingB, mT2suppression, useWidthPre;
194  double probQQtoQ, probStoUD, probSQtoQQ, probQQ1toQQ0, probQandQQ,
195  probQandS, probQandSinQQ, probQQ1corr, probQQ1corrInv, probQQ1norm,
196  probQQ1join[4], mesonRate[4][6], mesonRateSum[4], mesonMix1[2][6],
197  mesonMix2[2][6], etaSup, etaPrimeSup, decupletSup, baryonCGSum[6],
198  baryonCGMax[6], popcornRate, popcornSpair, popcornSmeson, scbBM[3],
199  popFrac, popS[3], dWT[3][7], lightLeadingBSup, heavyLeadingBSup;
200  double sigmaHad, widthPreStrange, widthPreDiquark;
201 
202  // Settings for thermal model.
203  bool thermalModel, mesonNonetL1;
204  double temperature, tempPreFactor;
205  int nNewQuark;
206  double mesMixRate1[2][6], mesMixRate2[2][6], mesMixRate3[2][6];
207  double baryonOctWeight[6][6][6][2], baryonDecWeight[6][6][6][2];
208 
209  // Settings used by both models.
210  bool closePacking;
211  double exponentMPI, exponentNSP;
212 
213  // Key = hadron id, value = list of constituent ids.
214  map< int, vector< pair<int,int> > > hadronConstIDs;
215  // Key = initial (di)quark id, value = list of possible hadron ids
216  // + nr in hadronConstIDs.
217  map< int, vector< pair<int,int> > > possibleHadrons;
218  // Key = initial (di)quark id, value = prefactor to multiply rate.
219  map< int, vector<double> > possibleRatePrefacs;
220  // Similar, but for combining the last two (di)quarks. Key = (di)quark pair.
221  map< pair<int,int>, vector< pair<int,int> > > possibleHadronsLast;
222  map< pair<int,int>, vector<double> > possibleRatePrefacsLast;
223 
224  // Selection in thermal model.
225  int hadronIDwin, idNewWin;
226  double hadronMassWin;
227 
228 };
229 
230 //==========================================================================
231 
232 // The StringZ class is used to sample the fragmentation function f(z).
233 
234 class StringZ : public PhysicsBase {
235 
236 public:
237 
238  // Constructor.
239  StringZ() : useNonStandC(), useNonStandB(), useNonStandH(), usePetersonC(),
240  usePetersonB(), usePetersonH(), mc2(), mb2(), aLund(), bLund(),
241  aExtraSQuark(), aExtraDiquark(), rFactC(), rFactB(), rFactH(), aNonC(),
242  aNonB(), aNonH(), bNonC(), bNonB(), bNonH(), epsilonC(), epsilonB(),
243  epsilonH(), stopM(), stopNF(), stopS() {}
244 
245  // Destructor.
246  virtual ~StringZ() {}
247 
248  // Initialize data members.
249  virtual void init();
250 
251  // Fragmentation function: top-level to determine parameters.
252  virtual double zFrag( int idOld, int idNew = 0, double mT2 = 1.);
253 
254  // Parameters for stopping in the middle; overloaded for Hidden Valley.
255  virtual double stopMass() {return stopM;}
256  virtual double stopNewFlav() {return stopNF;}
257  virtual double stopSmear() {return stopS;}
258 
259  // a and b fragmentation parameters needed in some operations.
260  virtual double aAreaLund() {return aLund;}
261  virtual double bAreaLund() {return bLund;}
262 
263  // Method to derive bLund from <z> (for fixed a and reference mT2).
264  bool deriveBLund();
265 
266 protected:
267 
268  // Constants: could only be changed in the code itself.
269  static const double CFROMUNITY, AFROMZERO, AFROMC, EXPMAX;
270 
271  // Initialization data, to be read from Settings.
272  bool useNonStandC, useNonStandB, useNonStandH,
273  usePetersonC, usePetersonB, usePetersonH;
274  double mc2, mb2, aLund, bLund, aExtraSQuark, aExtraDiquark, rFactC,
275  rFactB, rFactH, aNonC, aNonB, aNonH, bNonC, bNonB, bNonH,
276  epsilonC, epsilonB, epsilonH, stopM, stopNF, stopS;
277 
278  // Fragmentation function: select z according to provided parameters.
279  double zLund( double a, double b, double c = 1.);
280  double zPeterson( double epsilon);
281 
282 };
283 
284 //==========================================================================
285 
286 // The StringPT class is used to select select transverse momenta.
287 
288 class StringPT : public PhysicsBase {
289 
290 public:
291 
292  // Constructor.
293  StringPT() : useWidthPre(), sigmaQ(), enhancedFraction(), enhancedWidth(),
294  sigma2Had(), widthPreStrange(), widthPreDiquark(), thermalModel(),
295  temperature(), tempPreFactor(), fracSmallX(), closePacking(),
296  exponentMPI(), exponentNSP() {}
297 
298  // Destructor.
299  virtual ~StringPT() {}
300 
301  // Initialize data members.
302  virtual void init();
303 
304  // General function, return px and py as a pair in the same call
305  // in either model.
306  pair<double, double> pxy(int idIn, double nNSP = 0.0) {
307  return (thermalModel ? pxyThermal(idIn, nNSP) :
308  pxyGauss(idIn, nNSP)); }
309  pair<double, double> pxyGauss(int idIn = 0, double nNSP = 0.0);
310  pair<double, double> pxyThermal(int idIn, double nNSP = 0.0);
311 
312  // Gaussian suppression of given pT2; used in MiniStringFragmentation.
313  double suppressPT2(double pT2) { return (thermalModel ?
314  exp(-sqrt(pT2)/temperature) : exp(-pT2/sigma2Had)); }
315 
316 protected:
317 
318  // Constants: could only be changed in the code itself.
319  static const double SIGMAMIN;
320 
321  // Initialization data, to be read from Settings.
322  // Gaussian model.
323  bool useWidthPre;
324  double sigmaQ, enhancedFraction, enhancedWidth, sigma2Had,
325  widthPreStrange, widthPreDiquark;
326  // Thermal model.
327  bool thermalModel;
328  double temperature, tempPreFactor, fracSmallX;
329  // Both.
330  bool closePacking;
331  double exponentMPI, exponentNSP;
332 
333 private:
334 
335  // Evaluate Bessel function K_{1/4}(x).
336  double BesselK14(double x);
337 
338 };
339 
340 //==========================================================================
341 
342 } // end namespace Pythia8
343 
344 #endif // Pythia8_FragmentationFlavZpT_H