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) 2018 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() {singlets.resize(0);}
66 
67  // Initialize and save pointers.
68  void init(Info* infoPtrIn, Settings& settings, StringFlav* flavSelPtrIn);
69 
70  // Number of colour singlets.
71  int size() const {return singlets.size();}
72 
73  // Overload index operator to access separate colour singlets.
74  ColSinglet& operator[](int iSub) {return singlets[iSub];}
75 
76  // Clear contents.
77  void clear() {singlets.resize(0);}
78 
79  // Insert a new colour singlet system in ascending mass order.
80  // Calculate its properties. Join nearby partons.
81  bool insert( vector<int>& iPartonIn, Event& event);
82 
83  // Erase a colour singlet system. (Rare operation.)
84  void erase(int iSub) {singlets.erase(singlets.begin() + iSub);}
85 
86  // Collect all partons of singlet to be consecutively ordered.
87  void collect(int iSub, Event& event, bool skipTrivial = true);
88 
89  // Find to which singlet system a particle belongs.
90  int findSinglet(int i);
91 
92  // List all currently identified singlets.
93  void list() const;
94 
95  // Rapidity range [y_min, y_max] of all string pieces in all singlets.
96  // Only used when stringPT:closePacking is on.
97  vector< vector< pair<double,double> > > rapPairs;
98 
99 private:
100 
101  // Constants: could only be changed in the code itself.
102  static const double CONSTITUENTMASS;
103 
104  // Pointer to various information on the generation.
105  Info* infoPtr;
106 
107  // Pointer to class for flavour generation.
108  StringFlav* flavSelPtr;
109 
110  // Initialization data, to be read from Settings.
111  double mJoin, mJoinJunction, mStringMin;
112 
113  // List of all separate colour singlets.
114  vector<ColSinglet> singlets;
115 
116  // Join two legs of junction to a diquark for small invariant masses.
117  bool joinJunction( vector<int>& iPartonIn, Event& event,
118  double massExcessIn);
119 
120 };
121 
122 //==========================================================================
123 
124 // The StringRegion class contains the information related to
125 // one string section in the evolution of a multiparton system.
126 // Only to be used inside StringFragmentation and MiniStringFragmentation,
127 // so no private members.
128 
129 class StringRegion {
130 
131 public:
132 
133  // Constructor.
134  StringRegion() : isSetUp(false), isEmpty(true), w2(0.), xPosProj(0.),
135  xNegProj(0.), pxProj(0.), pyProj(0.) {}
136 
137  // Constants: could only be changed in the code itself.
138  static const double MJOIN, TINY;
139 
140  // Data members.
141  bool isSetUp, isEmpty;
142  Vec4 pPos, pNeg, eX, eY, pPosMass, pNegMass, massOffset;
143  double w2, xPosProj, xNegProj, pxProj, pyProj;
144 
145  // Calculate offset of the region from parton list. Special junction case.
146  Vec4 gluonOffset(vector<int>& iSys, Event& event, int iPos, int iNeg);
147  Vec4 gluonOffsetJRF(vector<int>& iSys, Event& event, int iPos, int iNeg,
148  RotBstMatrix MtoJRF);
149 
150  // If massive case, the offset of the initial regions is calculated.
151  bool massiveOffset(int iPos, int iNeg, int iMax, int id1, int id2,
152  double mc, double mb);
153 
154  // Set up four-vectors for longitudinal and transverse directions.
155  void setUp(Vec4 p1, Vec4 p2, bool isMassless = false);
156 
157  // Construct a four-momentum from (x+, x-, px, py).
158  Vec4 pHad( double xPosIn, double xNegIn, double pxIn, double pyIn)
159  { return xPosIn * pPos + xNegIn * pNeg + pxIn * eX + pyIn * eY; }
160 
161  // Project a four-momentum onto (x+, x-, px, py). Read out projection.
162  void project(Vec4 pIn);
163  void project( double pxIn, double pyIn, double pzIn, double eIn)
164  { project( Vec4( pxIn, pyIn, pzIn, eIn) ); }
165  double xPos() const {return xPosProj;}
166  double xNeg() const {return xNegProj;}
167  double px() const {return pxProj;}
168  double py() const {return pyProj;}
169 
170 };
171 
172 //==========================================================================
173 
174 // The StringSystem class contains the complete set of all string regions.
175 // Only to be used inside StringFragmentation, so no private members.
176 
177 class StringSystem {
178 
179 public:
180 
181  // Constructor.
182  StringSystem() {}
183 
184  // Set up system from parton list.
185  void setUp(vector<int>& iSys, Event& event);
186 
187  // Calculate string region from (iPos, iNeg) pair.
188  int iReg( int iPos, int iNeg) const
189  {return (iPos * (indxReg - iPos)) / 2 + iNeg;}
190 
191  // Reference to string region specified by (iPos, iNeg) pair.
192  StringRegion& region(int iPos, int iNeg) {return system[iReg(iPos, iNeg)];}
193 
194  // Reference to low string region specified either by iPos or iNeg.
195  StringRegion& regionLowPos(int iPos) {
196  return system[iReg(iPos, iMax - iPos)]; }
197  StringRegion& regionLowNeg(int iNeg) {
198  return system[iReg(iMax - iNeg, iNeg)]; }
199 
200  // Main content: a vector with all the string regions of the system.
201  vector<StringRegion> system;
202 
203  // Other data members.
204  int sizePartons, sizeStrings, sizeRegions, indxReg, iMax;
205  double mJoin, m2Join;
206 
207 };
208 
209 //==========================================================================
210 
211 // The StringVertex class contains the space-time vertex location information
212 // stored during the fragmentation process. No private members.
213 
215 
216 public:
217 
218  // Constructors.
219  StringVertex(bool fromPosIn = true, int iRegPosIn = 0,
220  int iRegNegIn = 0, double xRegPosIn = 0., double xRegNegIn = 0.)
221  : fromPos(fromPosIn), iRegPos(iRegPosIn), iRegNeg(iRegNegIn),
222  xRegPos(xRegPosIn), xRegNeg(xRegNegIn) { }
223 
224  StringVertex(const StringVertex& v): fromPos(v.fromPos),
225  iRegPos(v.iRegPos), iRegNeg(v.iRegNeg),
226  xRegPos(v.xRegPos), xRegNeg(v.xRegNeg) { }
227 
228  StringVertex& operator = (const StringVertex& v) {if (this != &v)
229  {fromPos = v.fromPos; iRegPos = v.iRegPos; iRegNeg = v.iRegNeg;
230  xRegPos = v.xRegPos; xRegNeg = v.xRegNeg;} return *this; }
231 
232  // Variable members.
233  bool fromPos;
234  int iRegPos, iRegNeg;
235  double xRegPos, xRegNeg;
236 
237 };
238 
239 //==========================================================================
240 
241 } // end namespace Pythia8
242 
243 #endif // Pythia8_FragmentationSystems_H
Definition: AgUStep.h:26