PYTHIA  8.313
StringFragmentation.h
1 // StringFragmentation.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 the classes for string fragmentation.
7 // StringEnd: keeps track of the fragmentation step.
8 // StringFragmentation: is the top-level class.
9 
10 #ifndef Pythia8_StringFragmentation_H
11 #define Pythia8_StringFragmentation_H
12 
13 #include "Pythia8/FragmentationModel.h"
14 
15 namespace Pythia8 {
16 
17 //==========================================================================
18 
19 // The StringEnd class contains the information related to
20 // one of the current endpoints of the string system.
21 // Only to be used inside StringFragmentation, so no private members.
22 
23 class StringEnd {
24 
25 public:
26 
27  // Constructor.
28  StringEnd() : particleDataPtr(), flavSelPtr(), pTSelPtr(), zSelPtr(),
29  fromPos(), thermalModel(), mT2suppression(), iEnd(), iMax(), idHad(),
30  iPosOld(), iNegOld(), iPosNew(), iNegNew(), hadSoFar(), colOld(), colNew(),
31  pxOld(), pyOld(), pxNew(), pyNew(), pxHad(), pyHad(), mHad(), mT2Had(),
32  GammaOld(), GammaNew(), xPosOld(), xPosNew(), xPosHad(), xNegOld(),
33  xNegNew(), xNegHad(), aLund(), bLund(), iPosOldPrev(), iNegOldPrev(),
34  colOldPrev(), pxOldPrev(), pyOldPrev(), GammaOldPrev(), xPosOldPrev(),
35  xNegOldPrev() {}
36 
37  // Save pointers.
38  void init( ParticleData* particleDataPtrIn, StringFlav* flavSelPtrIn,
39  StringPT* pTSelPtrIn, StringZ* zSelPtrIn, Settings& settings) {
40  particleDataPtr = particleDataPtrIn; flavSelPtr = flavSelPtrIn;
42  pTSelPtr = pTSelPtrIn; zSelPtr = zSelPtrIn;
43  bLund = zSelPtr->bAreaLund(); aLund = zSelPtr->aAreaLund();
44  thermalModel = settings.flag("StringPT:thermalModel");
45  mT2suppression = settings.flag("StringPT:mT2suppression");
46  closePacking = settings.flag("ClosePacking:doClosePacking"); }
47 
48  // Set up initial endpoint values from input.
49  void setUp(bool fromPosIn, int iEndIn, int idOldIn, int iMaxIn,
50  double pxIn, double pyIn, double GammaIn, double xPosIn,
51  double xNegIn, int colIn);
52 
53  // Fragment off one hadron from the string system, in flavour and pT.
54  void newHadron(double kappaModifier, bool forbidPopcornNow = false,
55  double strangeJunc = 0., double probQQmod = 1.);
56 
57  // Fragment off one hadron from the string system, in momentum space,
58  // by taking steps either from positive or from negative end.
59  Vec4 kinematicsHadron(StringSystem& system, StringVertex& newVertex,
60  double zHad);
61 
62  // Generate momentum for some possible next hadron, based on mean values
63  // to get an estimate for rapidity and pT.
64  Vec4 kinematicsHadronTmp(StringSystem system, Vec4 pRem, double phi,
65  double mult);
66 
67  // Update string end information after a hadron has been removed.
68  void update();
69  void storePrev();
70  void updateToPrev();
71 
72  // Constants: could only be changed in the code itself.
73  static const double TINY, PT2SAME, MEANMMIN, MEANM, MEANPT;
74 
75  // Pointer to the particle data table.
77 
78  // Pointers to classes for flavour, pT and z generation.
80  StringPT* pTSelPtr;
81  StringZ* zSelPtr;
82 
83  // Local copy of flavSelPtr for modified flavour selection.
85 
86  // Data members.
87  bool fromPos, thermalModel, mT2suppression, closePacking;
88  int iEnd, iMax, idHad, iPosOld, iNegOld, iPosNew, iNegNew, hadSoFar,
89  colOld, colNew;
90  double pxOld, pyOld, pxNew, pyNew, pxHad, pyHad, mHad, mT2Had,
91  GammaOld, GammaNew, xPosOld, xPosNew, xPosHad, xNegOld, xNegNew,
92  xNegHad, aLund, bLund;
93  int iPosOldPrev, iNegOldPrev, colOldPrev;
94  double pxOldPrev, pyOldPrev, GammaOldPrev, xPosOldPrev, xNegOldPrev;
95  FlavContainer flavOld, flavNew, flavOldPrev;
96  Vec4 pHad, pSoFar;
97 
98 };
99 
100 //==========================================================================
101 
102 // The StringFragmentation class contains the top-level routines
103 // to fragment a colour singlet partonic system.
104 
106 
107 public:
108 
109  // Constructor.
111  FragmentationModel(), flavRopePtr(),
112  closePacking(), setVertices(), constantTau(), smearOn(),
113  traceColours(false), hadronVertex(), stopMass(), stopNewFlav(),
114  stopSmear(), pNormJunction(), pMaxJunction(), eBothLeftJunction(),
115  eMaxLeftJunction(), eMinLeftJunction(), mJoin(), bLund(),
116  closePackingFluxRatio(1.), closePackingPT20(1.), pT20(), xySmear(),
117  maxSmear(), maxTau(), kappaVtx(), mc(), mb(), hasJunction(), isClosed(),
118  iPos(), iNeg(), nExtraJoin(), w2Rem(), stopMassNow(), idDiquark(),
119  legMin(), legMid() {}
120 
121  // Initialize and save pointers.
122  bool init(StringFlav* flavSelPtrIn = nullptr, StringPT* pTSelPtrIn = nullptr,
123  StringZ* zSelPtrIn = nullptr, FragModPtr fragModPtrIn = nullptr) override;
124 
125  // Do the fragmentation: driver routine.
126  bool fragment(int iSub, ColConfig& colConfig, Event& event,
127  bool isDiff = false, bool systemRecoil = true) override;
128 
129  // Local copy of flavSelPtr for modified flavour selection.
131 
132  // Find the boost matrix to the rest frame of a junction.
133  Vec4 junctionRestFrame(const Vec4& p0, const Vec4& p1, const Vec4& p2,
134  const bool angleCheck = true) const;
135 
136 private:
137 
138  // Constants: could only be changed in the code itself.
139  static const int NTRYFLAV, NTRYJOIN, NSTOPMASS,
140  NTRYJNMATCH, NTRYJRFEQ, NTRYSMEAR, MAXVETOFINTWO;
141  static const double FACSTOPMASS, CLOSEDM2MAX, CLOSEDM2FRAC, EXPMAX,
142  MATCHPOSNEG, M2MINJRF, EMINJRF, EEXTRAJNMATCH,
143  MDIQUARKMIN, CONVJRFEQ, CHECKPOS;
144 
145  // Pointer to flavour-composition-changing ropes.
146  FragModPtr flavRopePtr;
147 
148  // Initialization data, read from Settings.
149  bool closePacking, setVertices, constantTau, smearOn,
150  traceColours, hardRemn, doStrangeJunc;
151  int hadronVertex;
152  double stopMass, stopNewFlav, stopSmear, pNormJunction, pMaxJunction,
153  eBothLeftJunction, eMaxLeftJunction, eMinLeftJunction,
154  mJoin, bLund, closePackingFluxRatio, closePackingPT20,
155  qqSupPar, qqSupAnti, pT20, xySmear, maxSmear, maxTau,
156  kappaVtx, mc, mb, dampPopcorn, aRemn, bRemn, strangeJuncParm;
157 
158  // Data members.
159  bool hasJunction, isClosed;
160  int iPos, iNeg, nExtraJoin;
161  double w2Rem, stopMassNow, kappaModifier, probQQmod;
162  Vec4 pSum, pRem, pJunctionHadrons;
163 
164  // UserHooks flags.
165  bool doChangeFragPar = false, doVetoFrag = false;
166 
167  // List of partons in string system.
168  vector<int> iParton, iPartonMinLeg, iPartonMidLeg, iPartonMax;
169 
170  // Vertex information from the fragmentation process.
171  vector<StringVertex> stringVertices, legMinVertices, legMidVertices;
172  StringVertex newVertex;
173 
174  // Variables used for JRF finding.
175  vector<Vec4> listJRF;
176  vector<double> weightJRF;
177  int iLeg[3], idLeg[3], legEnd[3];
178  double weightSum, pSumJRF, m2Leg[3];
179  Vec4 pLeg[3];
180  bool lastJRF, endpoint[3];
181 
182  // Boost from/to rest frame of a junction to original frame.
183  RotBstMatrix MfromJRF, MtoJRF;
184 
185  // Information on diquark created at the junction.
186  int idDiquark;
187 
188  // Fictitious opposing partons in JRF: string ends for vertex location.
189  Vec4 pMinEnd, pMidEnd;
190 
191  // Temporary event record for the produced particles.
192  Event hadrons;
193 
194  // Information on the system of string regions.
195  StringSystem system, systemMin, systemMid;
196 
197  // Information on the two current endpoints of the fragmenting
198  // system, with backup copies if there is a veto.
199  StringEnd posEnd, negEnd, posEndSave, negEndSave;
200 
201  // Find region where to put first string break for closed gluon loop.
202  vector<int> findFirstRegion(int iSub, const ColConfig& colConfig,
203  const Event& event) const;
204 
205  // Set flavours and momentum position for initial string endpoints.
206  void setStartEnds(int idPos, int idNeg, const StringSystem& systemNow,
207  int legNow = 3);
208 
209  // Check remaining energy-momentum whether it is OK to continue.
210  bool energyUsedUp(bool fromPos);
211 
212  // Produce the final two partons to complete the system.
213  bool finalTwo(bool fromPos, const Event& event, bool usedPosJun,
214  bool usedNegJun);
215 
216  // Final region information.
217  Vec4 pPosFinalReg, pNegFinalReg, eXFinalReg, eYFinalReg;
218 
219  // Set hadron production points in space-time picture.
220  bool setHadronVertices(Event& event);
221 
222  // Construct a special joining region for the final two hadrons.
223  StringRegion finalRegion();
224 
225  // Store the hadrons in the normal event record, ordered from one end.
226  void store(Event& event);
227 
228  // Fragment off two of the string legs in to a junction.
229  bool fragmentToJunction(Event& event,
230  vector< vector< pair<double,double> > >& rapPairs);
231 
232  // Initially considered legs from the junction.
233  int legMin, legMid;
234 
235  // Functions used in JRF iterative procedure.
236  bool collinearPair(Event& event);
237  bool perturbedJRF(Event& event);
238  int updateLegs(Event& event, Vec4 vJunIn, bool juncCoM = false);
239  double updateWeights(double pSmall, Vec4 vJunIn);
240  void nextParton(Event& event, int leg);
241 
242  // Join extra nearby partons when stuck.
243  int extraJoin(double facExtra, Event& event);
244 
245  // Get the number of nearby strings given the energies.
246  void kappaEffModifier(StringSystem& systemNow,
247  StringEnd end, bool fromPos, vector<int> partonList,
248  vector< vector< pair<double,double> > >& rapPairs,
249  double mRem, Event& event);
250 
251  double yMax(Particle pIn, double mTiny) {
252  double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, pIn.mT()) );
253  return (pIn.pz() > 0) ? temp : -temp; }
254 
255 };
256 
257 //==========================================================================
258 
259 } // end namespace Pythia8
260 
261 #endif // Pythia8_StringFragmentation_H
Definition: FragmentationSystems.h:223
bool flag(string keyIn)
Give back current value, with check that key exists.
Definition: Settings.cc:1599
virtual double aAreaLund()
a and b fragmentation parameters needed in some operations.
Definition: FragmentationFlavZpT.h:299
bool fromPos
Data members.
Definition: StringFragmentation.h:87
The Event class holds all info on the generated event.
Definition: Event.h:408
void newHadron(double kappaModifier, bool forbidPopcornNow=false, double strangeJunc=0., double probQQmod=1.)
Fragment off one hadron from the string system, in flavour and pT.
Definition: StringFragmentation.cc:64
The StringPT class is used to select select transverse momenta.
Definition: FragmentationFlavZpT.h:326
static const double PT2SAME
Assume two (eX, eY) regions are related if pT2 differs by less.
Definition: StringFragmentation.h:73
StringFlav * flavSelPtr
Pointers to classes for flavour, pT and z generation.
Definition: StringFragmentation.h:79
ParticleData * particleDataPtr
Pointer to the particle data table.
Definition: StringFragmentation.h:76
Definition: StringFragmentation.h:105
Definition: FragmentationSystems.h:185
StringFragmentation()
Constructor.
Definition: StringFragmentation.h:110
Definition: FragmentationSystems.h:136
StringEnd()
Constructor.
Definition: StringFragmentation.h:28
The StringZ class is used to sample the fragmentation function f(z).
Definition: FragmentationFlavZpT.h:265
void updateToPrev()
Update string end information to previous string break.
Definition: StringFragmentation.cc:568
Vec4 kinematicsHadron(StringSystem &system, StringVertex &newVertex, double zHad)
Definition: StringFragmentation.cc:129
Definition: StringFragmentation.h:23
Vec4 kinematicsHadronTmp(StringSystem system, Vec4 pRem, double phi, double mult)
Definition: StringFragmentation.cc:317
void storePrev()
Store old string end information.
Definition: StringFragmentation.cc:550
Definition: Event.h:32
StringFlav flavSelNow
Local copy of flavSelPtr for modified flavour selection.
Definition: StringFragmentation.h:84
The ColConfig class describes the colour configuration of the whole event.
Definition: FragmentationSystems.h:60
void update()
Update string end information after a hadron has been removed.
Definition: StringFragmentation.cc:532
Definition: Basics.h:254
Definition: FragmentationFlavZpT.h:41
FragmentationModel is the base class for handling fragmentation algorithms.
Definition: FragmentationModel.h:28
double phi(const Vec4 &v1, const Vec4 &v2)
phi is azimuthal angle between v1 and v2 around z axis.
Definition: Basics.cc:704
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:84
StringFlav flavSelNow
Local copy of flavSelPtr for modified flavour selection.
Definition: StringFragmentation.h:130
static const double MEANMMIN
Definition: StringFragmentation.h:73
void init(ParticleData *particleDataPtrIn, StringFlav *flavSelPtrIn, StringPT *pTSelPtrIn, StringZ *zSelPtrIn, Settings &settings)
Save pointers.
Definition: StringFragmentation.h:38
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
static const double TINY
Constants: could only be changed in the code itself.
Definition: StringFragmentation.h:73
void setUp(bool fromPosIn, int iEndIn, int idOldIn, int iMaxIn, double pxIn, double pyIn, double GammaIn, double xPosIn, double xNegIn, int colIn)
Set up initial endpoint values from input.
Definition: StringFragmentation.cc:39
Definition: Basics.h:32
Definition: Settings.h:196