1 // FragmentationSystems.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2008 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
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.
12 #ifndef Pythia8_FragmentationSystems_H
13 #define Pythia8_FragmentationSystems_H
17 #include "FragmentationFlavZpT.h"
18 #include "ParticleData.h"
19 #include "PythiaStdlib.h"
24 //**************************************************************************
26 // The ColSinglet class contains info on an individual singlet.
27 // Only to be used inside ColConfig, so no private members.
34 ColSinglet() : pSum(0., 0., 0., 0.), mass(0.), massExcess(0.),
35 hasJunction(false), isClosed(false), isCollected(false) {}
36 ColSinglet(vector<int>& iPartonIn, Vec4 pSumIn, double massIn,
37 double massExcessIn, bool hasJunctionIn = false,
38 bool isClosedIn = false, bool isCollectedIn = false)
39 : iParton(iPartonIn), pSum(pSumIn), mass(massIn),
40 massExcess(massExcessIn), hasJunction(hasJunctionIn),
41 isClosed(isClosedIn), isCollected(isCollectedIn) {}
43 // Size of iParton array.
44 int size() const { return iParton.size();}
49 double mass, massExcess;
50 bool hasJunction, isClosed, isCollected;
54 //**************************************************************************
56 // The ColConfig class describes the colour configuration of the whole event.
65 // Initialize and save pointers.
66 void init(StringFlav* flavSelPtrIn);
68 // Number of colour singlets.
69 int size() const {return singlets.size();}
71 // Overload index operator to access separate colour singlets.
72 ColSinglet& operator[](int iSub) {return singlets[iSub];}
75 void clear() {singlets.resize(0);}
77 // Insert a new colour singlet system in ascending mass order.
78 // Calculate its properties. Join nearby partons.
79 void insert( vector<int>& iPartonIn, Event& event);
81 // Collect all partons of singlet to be consecutively ordered.
82 void collect(int iSub, Event& event);
84 // List all currently identified singlets.
85 void list(ostream& os = cout);
89 // Pointer to class for flavour generation.
90 StringFlav* flavSelPtr;
92 // Initialization data, to be read from Settings.
93 double mJoin, mJoinJunction, mStringMin;
95 // List of all separate colour singlets.
96 vector<ColSinglet> singlets;
98 // Join two legs of junction to a diquark for small invariant masses.
99 bool joinJunction( vector<int>& iPartonIn, Event& event,
100 double massExcessIn);
104 //**************************************************************************
106 // The StringRegion class contains the information related to
107 // one string section in the evolution of a multiparton system.
108 // Only to be used inside StringFragmentation and MiniStringFragmentation,
109 // so no private members.
116 StringRegion() : isSetUp(false), isEmpty(true) {}
118 // Constants: could only be changed in the code itself.
119 static const double MJOIN, TINY;
122 bool isSetUp, isEmpty;
123 Vec4 pPos, pNeg, eX, eY;
124 double w2, xPosProj, xNegProj, pxProj, pyProj;
126 // Set up four-vectors for longitudinal and transverse directions.
127 void setUp(Vec4 p1, Vec4 p2, bool isMassless = false);
129 // Construct a four-momentum from (x+, x-, px, py).
130 Vec4 pHad( double xPosIn, double xNegIn, double pxIn, double pyIn)
131 { return xPosIn * pPos + xNegIn * pNeg + pxIn * eX + pyIn * eY; }
133 // Project a four-momentum onto (x+, x-, px, py). Read out projection.
134 void project(Vec4 pIn);
135 void project( double pxIn, double pyIn, double pzIn, double eIn)
136 { project( Vec4( pxIn, pyIn, pzIn, eIn) ); }
137 double xPos() const {return xPosProj;}
138 double xNeg() const {return xNegProj;}
139 double px() const {return pxProj;}
140 double py() const {return pyProj;}
144 //**************************************************************************
146 // The StringSystem class contains the complete set of all string regions.
147 // Only to be used inside StringFragmentation, so no private members.
156 // Set up system from parton list.
157 void setUp(vector<int>& iSys, Event& event);
159 // Calculate string region from (iPos, iNeg) pair.
160 int iReg( int iPos, int iNeg) const
161 {return (iPos * (indxReg - iPos)) / 2 + iNeg;}
163 // Reference to string region specified by (iPos, iNeg) pair.
164 StringRegion& region(int iPos, int iNeg) {return system[iReg(iPos, iNeg)];}
166 // Reference to low string region specified either by iPos or iNeg.
167 StringRegion& regionLowPos(int iPos) {return system[iReg(iPos, iMax - iPos)];}
168 StringRegion& regionLowNeg(int iNeg) {return system[iReg(iMax - iNeg, iNeg)];}
170 // Main content: a vector with all the string regions of the system.
171 vector<StringRegion> system;
173 // Other data members.
174 int sizePartons, sizeStrings, sizeRegions, indxReg, iMax;
175 double mJoin, m2Join;
179 //**************************************************************************
181 } // end namespace Pythia8
183 #endif // Pythia8_FragmentationSystems_H