1 // FragmentationSystems.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2010 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"
19 #include "ParticleData.h"
20 #include "PythiaStdlib.h"
25 //==========================================================================
27 // The ColSinglet class contains info on an individual singlet.
28 // Only to be used inside ColConfig, so no private members.
35 ColSinglet() : pSum(0., 0., 0., 0.), mass(0.), massExcess(0.),
36 hasJunction(false), isClosed(false), isCollected(false) {}
37 ColSinglet(vector<int>& iPartonIn, Vec4 pSumIn, double massIn,
38 double massExcessIn, bool hasJunctionIn = false,
39 bool isClosedIn = false, bool isCollectedIn = false)
40 : iParton(iPartonIn), pSum(pSumIn), mass(massIn),
41 massExcess(massExcessIn), hasJunction(hasJunctionIn),
42 isClosed(isClosedIn), isCollected(isCollectedIn) {}
44 // Size of iParton array.
45 int size() const { return iParton.size();}
50 double mass, massExcess;
51 bool hasJunction, isClosed, isCollected;
55 //==========================================================================
57 // The ColConfig class describes the colour configuration of the whole event.
64 ColConfig() {singlets.resize(0);}
66 // Initialize and save pointers.
67 void init(Info* infoPtrIn, Settings& settings, StringFlav* flavSelPtrIn);
69 // Number of colour singlets.
70 int size() const {return singlets.size();}
72 // Overload index operator to access separate colour singlets.
73 ColSinglet& operator[](int iSub) {return singlets[iSub];}
76 void clear() {singlets.resize(0);}
78 // Insert a new colour singlet system in ascending mass order.
79 // Calculate its properties. Join nearby partons.
80 bool insert( vector<int>& iPartonIn, Event& event);
82 // Collect all partons of singlet to be consecutively ordered.
83 void collect(int iSub, Event& event);
85 // List all currently identified singlets.
86 void list(ostream& os = cout) const;
90 // Constants: could only be changed in the code itself.
91 static const double CONSTITUENTMASS;
93 // Pointer to various information on the generation.
96 // Pointer to class for flavour generation.
97 StringFlav* flavSelPtr;
99 // Initialization data, to be read from Settings.
100 double mJoin, mJoinJunction, mStringMin;
102 // List of all separate colour singlets.
103 vector<ColSinglet> singlets;
105 // Join two legs of junction to a diquark for small invariant masses.
106 bool joinJunction( vector<int>& iPartonIn, Event& event,
107 double massExcessIn);
111 //==========================================================================
113 // The StringRegion class contains the information related to
114 // one string section in the evolution of a multiparton system.
115 // Only to be used inside StringFragmentation and MiniStringFragmentation,
116 // so no private members.
123 StringRegion() : isSetUp(false), isEmpty(true) {}
125 // Constants: could only be changed in the code itself.
126 static const double MJOIN, TINY;
129 bool isSetUp, isEmpty;
130 Vec4 pPos, pNeg, eX, eY;
131 double w2, xPosProj, xNegProj, pxProj, pyProj;
133 // Set up four-vectors for longitudinal and transverse directions.
134 void setUp(Vec4 p1, Vec4 p2, bool isMassless = false);
136 // Construct a four-momentum from (x+, x-, px, py).
137 Vec4 pHad( double xPosIn, double xNegIn, double pxIn, double pyIn)
138 { return xPosIn * pPos + xNegIn * pNeg + pxIn * eX + pyIn * eY; }
140 // Project a four-momentum onto (x+, x-, px, py). Read out projection.
141 void project(Vec4 pIn);
142 void project( double pxIn, double pyIn, double pzIn, double eIn)
143 { project( Vec4( pxIn, pyIn, pzIn, eIn) ); }
144 double xPos() const {return xPosProj;}
145 double xNeg() const {return xNegProj;}
146 double px() const {return pxProj;}
147 double py() const {return pyProj;}
151 //==========================================================================
153 // The StringSystem class contains the complete set of all string regions.
154 // Only to be used inside StringFragmentation, so no private members.
163 // Set up system from parton list.
164 void setUp(vector<int>& iSys, Event& event);
166 // Calculate string region from (iPos, iNeg) pair.
167 int iReg( int iPos, int iNeg) const
168 {return (iPos * (indxReg - iPos)) / 2 + iNeg;}
170 // Reference to string region specified by (iPos, iNeg) pair.
171 StringRegion& region(int iPos, int iNeg) {return system[iReg(iPos, iNeg)];}
173 // Reference to low string region specified either by iPos or iNeg.
174 StringRegion& regionLowPos(int iPos) {return system[iReg(iPos, iMax - iPos)];}
175 StringRegion& regionLowNeg(int iNeg) {return system[iReg(iMax - iNeg, iNeg)];}
177 // Main content: a vector with all the string regions of the system.
178 vector<StringRegion> system;
180 // Other data members.
181 int sizePartons, sizeStrings, sizeRegions, indxReg, iMax;
182 double mJoin, m2Join;
186 //==========================================================================
188 } // end namespace Pythia8
190 #endif // Pythia8_FragmentationSystems_H