StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
FragmentationSystems.h
1 // FragmentationSystems.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 auxiliary classes in the fragmentation process.
7 // ColSinglet contains info on an individual singlet.
8 // ColConfig describes the colour configuration of the whole event.
9 // StringRegion keeps track on string momenta and directions.
10 // StringSystem contains all the StringRegions of the colour singlet.
11 // StringVertex contains information on space-time location of string breaks.
12 
13 #ifndef Pythia8_FragmentationSystems_H
14 #define Pythia8_FragmentationSystems_H
15 
16 #include "Pythia8/Basics.h"
17 #include "Pythia8/Event.h"
18 #include "Pythia8/FragmentationFlavZpT.h"
19 #include "Pythia8/Info.h"
20 #include "Pythia8/ParticleData.h"
21 #include "Pythia8/PythiaStdlib.h"
22 #include "Pythia8/Settings.h"
23 
24 namespace Pythia8 {
25 
26 //==========================================================================
27 
28 // The ColSinglet class contains info on an individual singlet.
29 // Only to be used inside ColConfig, so no private members.
30 
31 class ColSinglet {
32 
33 public:
34 
35  // Constructors.
36  ColSinglet() : pSum(0., 0., 0., 0.), mass(0.), massExcess(0.),
37  hasJunction(false), isClosed(false), isCollected(false) {}
38  ColSinglet(vector<int>& iPartonIn, Vec4 pSumIn, double massIn,
39  double massExcessIn, bool hasJunctionIn = false,
40  bool isClosedIn = false, bool isCollectedIn = false)
41  : iParton(iPartonIn), pSum(pSumIn), mass(massIn),
42  massExcess(massExcessIn), hasJunction(hasJunctionIn),
43  isClosed(isClosedIn), isCollected(isCollectedIn) {}
44 
45  // Size of iParton array.
46  int size() const { return iParton.size();}
47 
48  // Stored quantities.
49  vector<int> iParton;
50  Vec4 pSum;
51  double mass, massExcess;
52  bool hasJunction, isClosed, isCollected;
53 
54 };
55 
56 //==========================================================================
57 
58 // The ColConfig class describes the colour configuration of the whole event.
59 
60 class ColConfig {
61 
62 public:
63 
64  // Constructor.
65  ColConfig() : infoPtr(), flavSelPtr(), mJoin(), mJoinJunction(),
66  mStringMin() {singlets.resize(0);}
67 
68  // Initialize and save pointers.
69  void init(Info* infoPtrIn, StringFlav* flavSelPtrIn);
70 
71  // Number of colour singlets.
72  int size() const {return singlets.size();}
73 
74  // Overload index operator to access separate colour singlets.
75  ColSinglet& operator[](int iSub) {return singlets[iSub];}
76 
77  // Clear contents.
78  void clear() {singlets.resize(0);}
79 
80  // Insert a new colour singlet system in ascending mass order.
81  // Calculate its properties. Join nearby partons.
82  bool insert( vector<int>& iPartonIn, Event& event);
83 
84  // Insert a new qqbar colour singlet system in ascending mass order.
85  // Calculate its properties.
86  bool simpleInsert( vector<int>& iPartonIn, Event& event,
87  bool fixOrder = false);
88 
89  // Erase a colour singlet system. (Rare operation.)
90  void erase(int iSub) {singlets.erase(singlets.begin() + iSub);}
91 
92  // Collect all partons of singlet to be consecutively ordered.
93  void collect(int iSub, Event& event, bool skipTrivial = true);
94 
95  // Find to which singlet system a particle belongs.
96  int findSinglet(int i);
97 
98  // List all currently identified singlets.
99  void list() const;
100 
101  // Rapidity range [y_min, y_max] of all string pieces in all singlets.
102  // Only used when stringPT:closePacking is on.
103  vector< vector< pair<double,double> > > rapPairs;
104 
105 private:
106 
107  // Constants: could only be changed in the code itself.
108  static const double CONSTITUENTMASS;
109 
110  // Pointer to various information on the generation.
111  Info* infoPtr;
112 
113  // Pointer to class for flavour generation.
114  StringFlav* flavSelPtr;
115 
116  // Initialization data, to be read from Settings.
117  double mJoin, mJoinJunction, mStringMin;
118 
119  // List of all separate colour singlets.
120  vector<ColSinglet> singlets;
121 
122  // Join two legs of junction to a diquark for small invariant masses.
123  bool joinJunction( vector<int>& iPartonIn, Event& event,
124  double massExcessIn);
125 
126 };
127 
128 //==========================================================================
129 
130 // The StringRegion class contains the information related to
131 // one string section in the evolution of a multiparton system.
132 // Only to be used inside StringFragmentation and MiniStringFragmentation,
133 // so no private members.
134 
135 class StringRegion {
136 
137 public:
138 
139  // Constructor.
140  StringRegion() : isSetUp(false), isEmpty(true), w2(0.), xPosProj(0.),
141  xNegProj(0.), pxProj(0.), pyProj(0.), colPos(0), colNeg(0) {}
142 
143  // Constants: could only be changed in the code itself.
144  static const double MJOIN, TINY;
145 
146  // Data members.
147  bool isSetUp, isEmpty;
148  Vec4 pPos, pNeg, eX, eY, pPosMass, pNegMass, massOffset;
149  double w2, xPosProj, xNegProj, pxProj, pyProj;
150  int colPos, colNeg;
151 
152  // Calculate offset of the region from parton list. Special junction case.
153  Vec4 gluonOffset(vector<int>& iSys, Event& event, int iPos, int iNeg);
154  Vec4 gluonOffsetJRF(vector<int>& iSys, Event& event, int iPos, int iNeg,
155  RotBstMatrix MtoJRF);
156 
157  // If massive case, the offset of the initial regions is calculated.
158  bool massiveOffset(int iPos, int iNeg, int iMax, int id1, int id2,
159  double mc, double mb);
160 
161  // Set up four-vectors for longitudinal and transverse directions.
162  void setUp(Vec4 p1, Vec4 p2, int col1, int col2, bool isMassless = false);
163 
164  // Construct a four-momentum from (x+, x-, px, py).
165  Vec4 pHad( double xPosIn, double xNegIn, double pxIn, double pyIn)
166  { return xPosIn * pPos + xNegIn * pNeg + pxIn * eX + pyIn * eY; }
167 
168  // Project a four-momentum onto (x+, x-, px, py). Read out projection.
169  void project(Vec4 pIn);
170  void project( double pxIn, double pyIn, double pzIn, double eIn)
171  { project( Vec4( pxIn, pyIn, pzIn, eIn) ); }
172  double xPos() const {return xPosProj;}
173  double xNeg() const {return xNegProj;}
174  double px() const {return pxProj;}
175  double py() const {return pyProj;}
176 
177 };
178 
179 //==========================================================================
180 
181 // The StringSystem class contains the complete set of all string regions.
182 // Only to be used inside StringFragmentation, so no private members.
183 
184 class StringSystem {
185 
186 public:
187 
188  // Constructor.
189  StringSystem() : sizePartons(), sizeStrings(), sizeRegions(), indxReg(),
190  iMax(), mJoin(), m2Join() {}
191 
192  // Set up system from parton list.
193  void setUp(vector<int>& iSys, Event& event);
194 
195  // Calculate string region from (iPos, iNeg) pair.
196  int iReg( int iPos, int iNeg) const
197  {return (iPos * (indxReg - iPos)) / 2 + iNeg;}
198 
199  // Reference to string region specified by (iPos, iNeg) pair.
200  StringRegion& region(int iPos, int iNeg) {return system[iReg(iPos, iNeg)];}
201 
202  // Reference to low string region specified either by iPos or iNeg.
203  StringRegion& regionLowPos(int iPos) {
204  return system[iReg(iPos, iMax - iPos)]; }
205  StringRegion& regionLowNeg(int iNeg) {
206  return system[iReg(iMax - iNeg, iNeg)]; }
207 
208  // Main content: a vector with all the string regions of the system.
209  vector<StringRegion> system;
210 
211  // Other data members.
212  int sizePartons, sizeStrings, sizeRegions, indxReg, iMax;
213  double mJoin, m2Join;
214 
215 };
216 
217 //==========================================================================
218 
219 // The StringVertex class contains the space-time vertex location information
220 // stored during the fragmentation process. No private members.
221 
222 class StringVertex {
223 
224 public:
225 
226  // Constructors.
227  StringVertex(bool fromPosIn = true, int iRegPosIn = 0,
228  int iRegNegIn = 0, double xRegPosIn = 0., double xRegNegIn = 0.)
229  : fromPos(fromPosIn), iRegPos(iRegPosIn), iRegNeg(iRegNegIn),
230  xRegPos(xRegPosIn), xRegNeg(xRegNegIn) { }
231 
232  StringVertex(const StringVertex& v): fromPos(v.fromPos),
233  iRegPos(v.iRegPos), iRegNeg(v.iRegNeg),
234  xRegPos(v.xRegPos), xRegNeg(v.xRegNeg) { }
235 
236  StringVertex& operator = (const StringVertex& v) {if (this != &v)
237  {fromPos = v.fromPos; iRegPos = v.iRegPos; iRegNeg = v.iRegNeg;
238  xRegPos = v.xRegPos; xRegNeg = v.xRegNeg;} return *this; }
239 
240  // Variable members.
241  bool fromPos;
242  int iRegPos, iRegNeg;
243  double xRegPos, xRegNeg;
244 
245 };
246 
247 //==========================================================================
248 
249 } // end namespace Pythia8
250 
251 #endif // Pythia8_FragmentationSystems_H
Definition: AgUStep.h:26