// FragmentationSystems.h is a part of the PYTHIA event generator. // Copyright (C) 2008 Torbjorn Sjostrand. // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details. // Please respect the MCnet Guidelines, see GUIDELINES for details. // This file contains auxiliary classes in the fragmentation process. // ColSinglet contains info on an individual singlet. // ColConfig describes the colour configuration of the whole event. // StringRegion keeps track on string momenta and directions. // StringSystem contains all the StringRegions of the colour singlet. #ifndef Pythia8_FragmentationSystems_H #define Pythia8_FragmentationSystems_H #include "Basics.h" #include "Event.h" #include "FragmentationFlavZpT.h" #include "ParticleData.h" #include "PythiaStdlib.h" #include "Settings.h" namespace Pythia8 { //************************************************************************** // The ColSinglet class contains info on an individual singlet. // Only to be used inside ColConfig, so no private members. class ColSinglet { public: // Constructors. ColSinglet() : pSum(0., 0., 0., 0.), mass(0.), massExcess(0.), hasJunction(false), isClosed(false), isCollected(false) {} ColSinglet(vector& iPartonIn, Vec4 pSumIn, double massIn, double massExcessIn, bool hasJunctionIn = false, bool isClosedIn = false, bool isCollectedIn = false) : iParton(iPartonIn), pSum(pSumIn), mass(massIn), massExcess(massExcessIn), hasJunction(hasJunctionIn), isClosed(isClosedIn), isCollected(isCollectedIn) {} // Size of iParton array. int size() const { return iParton.size();} // Stored quantities. vector iParton; Vec4 pSum; double mass, massExcess; bool hasJunction, isClosed, isCollected; }; //************************************************************************** // The ColConfig class describes the colour configuration of the whole event. class ColConfig { public: // Constructor. ColConfig() {} // Initialize and save pointers. void init(StringFlav* flavSelPtrIn); // Number of colour singlets. int size() const {return singlets.size();} // Overload index operator to access separate colour singlets. ColSinglet& operator[](int iSub) {return singlets[iSub];} // Clear contents. void clear() {singlets.resize(0);} // Insert a new colour singlet system in ascending mass order. // Calculate its properties. Join nearby partons. void insert( vector& iPartonIn, Event& event); // Collect all partons of singlet to be consecutively ordered. void collect(int iSub, Event& event); // List all currently identified singlets. void list(ostream& os = cout); private: // Pointer to class for flavour generation. StringFlav* flavSelPtr; // Initialization data, to be read from Settings. double mJoin, mJoinJunction, mStringMin; // List of all separate colour singlets. vector singlets; // Join two legs of junction to a diquark for small invariant masses. bool joinJunction( vector& iPartonIn, Event& event, double massExcessIn); }; //************************************************************************** // The StringRegion class contains the information related to // one string section in the evolution of a multiparton system. // Only to be used inside StringFragmentation and MiniStringFragmentation, // so no private members. class StringRegion { public: // Constructor. StringRegion() : isSetUp(false), isEmpty(true) {} // Constants: could only be changed in the code itself. static const double MJOIN, TINY; // Data members. bool isSetUp, isEmpty; Vec4 pPos, pNeg, eX, eY; double w2, xPosProj, xNegProj, pxProj, pyProj; // Set up four-vectors for longitudinal and transverse directions. void setUp(Vec4 p1, Vec4 p2, bool isMassless = false); // Construct a four-momentum from (x+, x-, px, py). Vec4 pHad( double xPosIn, double xNegIn, double pxIn, double pyIn) { return xPosIn * pPos + xNegIn * pNeg + pxIn * eX + pyIn * eY; } // Project a four-momentum onto (x+, x-, px, py). Read out projection. void project(Vec4 pIn); void project( double pxIn, double pyIn, double pzIn, double eIn) { project( Vec4( pxIn, pyIn, pzIn, eIn) ); } double xPos() const {return xPosProj;} double xNeg() const {return xNegProj;} double px() const {return pxProj;} double py() const {return pyProj;} }; //************************************************************************** // The StringSystem class contains the complete set of all string regions. // Only to be used inside StringFragmentation, so no private members. class StringSystem { public: // Constructor. StringSystem() {} // Set up system from parton list. void setUp(vector& iSys, Event& event); // Calculate string region from (iPos, iNeg) pair. int iReg( int iPos, int iNeg) const {return (iPos * (indxReg - iPos)) / 2 + iNeg;} // Reference to string region specified by (iPos, iNeg) pair. StringRegion& region(int iPos, int iNeg) {return system[iReg(iPos, iNeg)];} // Reference to low string region specified either by iPos or iNeg. StringRegion& regionLowPos(int iPos) {return system[iReg(iPos, iMax - iPos)];} StringRegion& regionLowNeg(int iNeg) {return system[iReg(iMax - iNeg, iNeg)];} // Main content: a vector with all the string regions of the system. vector system; // Other data members. int sizePartons, sizeStrings, sizeRegions, indxReg, iMax; double mJoin, m2Join; }; //************************************************************************** } // end namespace Pythia8 #endif // Pythia8_FragmentationSystems_H