Line data Source code
1 : // HadronLevel.h is a part of the PYTHIA event generator.
2 : // Copyright (C) 2015 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.
5 :
6 : // This file contains the main class for hadron-level generation.
7 : // HadronLevel: handles administration of fragmentation and decay.
8 :
9 : #ifndef Pythia8_HadronLevel_H
10 : #define Pythia8_HadronLevel_H
11 :
12 : #include "Pythia8/Basics.h"
13 : #include "Pythia8/BoseEinstein.h"
14 : #include "Pythia8/ColourTracing.h"
15 : #include "Pythia8/Event.h"
16 : #include "Pythia8/FragmentationFlavZpT.h"
17 : #include "Pythia8/FragmentationSystems.h"
18 : #include "Pythia8/HadronScatter.h"
19 : #include "Pythia8/HiddenValleyFragmentation.h"
20 : #include "Pythia8/Info.h"
21 : #include "Pythia8/JunctionSplitting.h"
22 : #include "Pythia8/MiniStringFragmentation.h"
23 : #include "Pythia8/ParticleData.h"
24 : #include "Pythia8/ParticleDecays.h"
25 : #include "Pythia8/PythiaStdlib.h"
26 : #include "Pythia8/RHadrons.h"
27 : #include "Pythia8/Settings.h"
28 : #include "Pythia8/StringFragmentation.h"
29 : #include "Pythia8/TimeShower.h"
30 :
31 : namespace Pythia8 {
32 :
33 : //==========================================================================
34 :
35 : // The HadronLevel class contains the top-level routines to generate
36 : // the transition from the partonic to the hadronic stage of an event.
37 :
38 0 : class HadronLevel {
39 :
40 : public:
41 :
42 : // Constructor.
43 0 : HadronLevel() {}
44 :
45 : // Initialize HadronLevel classes as required.
46 : bool init(Info* infoPtrIn, Settings& settings,
47 : ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
48 : Couplings* couplingsPtrIn, TimeShower* timesDecPtr,
49 : RHadrons* rHadronsPtrIn, DecayHandler* decayHandlePtr,
50 : vector<int> handledParticles);
51 :
52 : // Get pointer to StringFlav instance (needed by BeamParticle).
53 0 : StringFlav* getStringFlavPtr() {return &flavSel;}
54 :
55 : // Generate the next event.
56 : bool next(Event& event);
57 :
58 : // Special routine to allow more decays if on/off switches changed.
59 : bool moreDecays(Event& event);
60 :
61 : private:
62 :
63 : // Initialization data, read from Settings.
64 : bool doHadronize, doDecay, doBoseEinstein, allowRH;
65 : double mStringMin, eNormJunction, widthSepBE;
66 :
67 : // Settings for hadron scattering --rjc
68 : bool doHadronScatter, hsAfterDecay;
69 :
70 : // Pointer to various information on the generation.
71 : Info* infoPtr;
72 :
73 : // Pointer to the particle data table.
74 : ParticleData* particleDataPtr;
75 :
76 : // Pointer to the random number generator.
77 : Rndm* rndmPtr;
78 :
79 : // Pointers to Standard Model couplings.
80 : Couplings* couplingsPtr;
81 :
82 : // Configuration of colour-singlet systems.
83 : ColConfig colConfig;
84 :
85 : // Colour and mass information.
86 : vector<int> iParton, iJunLegA, iJunLegB, iJunLegC,
87 : iAntiLegA, iAntiLegB, iAntiLegC, iGluLeg;
88 : vector<double> m2Pair;
89 :
90 : // The generator class for normal string fragmentation.
91 : StringFragmentation stringFrag;
92 :
93 : // The generator class for special low-mass string fragmentation.
94 : MiniStringFragmentation ministringFrag;
95 :
96 : // The generator class for normal decays.
97 : ParticleDecays decays;
98 :
99 : // The generator class for hadron scattering --rjc
100 : HadronScatter hadronScatter;
101 :
102 : // The generator class for Bose-Einstein effects.
103 : BoseEinstein boseEinstein;
104 :
105 : // Classes for flavour, pT and z generation.
106 : StringFlav flavSel;
107 : StringPT pTSel;
108 : StringZ zSel;
109 :
110 : // Class for colour tracing.
111 : ColourTracing colTrace;
112 :
113 : // Junction splitting class.
114 : JunctionSplitting junctionSplitting;
115 :
116 : // The RHadrons class is used to fragment off and decay R-hadrons.
117 : RHadrons* rHadronsPtr;
118 :
119 : // Special class for Hidden-Valley hadronization. Not always used.
120 : HiddenValleyFragmentation hiddenvalleyFrag;
121 : bool useHiddenValley;
122 :
123 : // Special case: colour-octet onium decays, to be done initially.
124 : bool decayOctetOnia(Event& event);
125 :
126 : // Trace colour flow in the event to form colour singlet subsystems.
127 : bool findSinglets(Event& event);
128 :
129 : };
130 :
131 : //==========================================================================
132 :
133 : } // end namespace Pythia8
134 :
135 : #endif // Pythia8_HadronLevel_H
|