StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HadronLevel.h
1 // HadronLevel.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 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 #include "Pythia8/UserHooks.h"
31 
32 namespace Pythia8 {
33 
34 //==========================================================================
35 
36 // The HadronLevel class contains the top-level routines to generate
37 // the transition from the partonic to the hadronic stage of an event.
38 
39 class HadronLevel {
40 
41 public:
42 
43  // Constructor.
44  HadronLevel() {}
45 
46  // Initialize HadronLevel classes as required.
47  bool init(Info* infoPtrIn, Settings& settings,
48  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
49  Couplings* couplingsPtrIn, TimeShower* timesDecPtr,
50  RHadrons* rHadronsPtrIn, DecayHandler* decayHandlePtr,
51  vector<int> handledParticles, UserHooks* userHooksPtrIn);
52 
53  // Get pointer to StringFlav instance (needed by BeamParticle).
54  StringFlav* getStringFlavPtr() {return &flavSel;}
55 
56  // Generate the next event.
57  bool next(Event& event);
58 
59  // Special routine to allow more decays if on/off switches changed.
60  bool moreDecays(Event& event);
61 
62 private:
63 
64  // Constants: could only be changed in the code itself.
65  static const double MTINY;
66 
67  // Initialization data, read from Settings.
68  bool doHadronize, doDecay, doBoseEinstein, allowRH, closePacking;
69  double mStringMin, eNormJunction, widthSepBE;
70 
71  // Settings for hadron scattering.
72  bool doHadronScatter, hsAfterDecay;
73  int hadronScatMode;
74 
75  // Pointer to various information on the generation.
76  Info* infoPtr;
77 
78  // Pointer to the particle data table.
79  ParticleData* particleDataPtr;
80 
81  // Pointer to the random number generator.
82  Rndm* rndmPtr;
83 
84  // Pointer to the user hooks object.
85  UserHooks* userHooksPtr;
86 
87  // Pointers to Standard Model couplings.
88  Couplings* couplingsPtr;
89 
90  // Configuration of colour-singlet systems.
91  ColConfig colConfig;
92 
93  // Colour and mass information.
94  vector<int> iParton, iJunLegA, iJunLegB, iJunLegC,
95  iAntiLegA, iAntiLegB, iAntiLegC, iGluLeg;
96  vector<double> m2Pair;
97 
98  // The generator class for normal string fragmentation.
99  StringFragmentation stringFrag;
100 
101  // The generator class for special low-mass string fragmentation.
102  MiniStringFragmentation ministringFrag;
103 
104  // The generator class for normal decays.
105  ParticleDecays decays;
106 
107  // The generator class for hadron scattering.
108  HadronScatter hadronScatter;
109 
110  // Class for event geometry for Rope Hadronization. Production vertices.
111  Ropewalk ropewalk;
112  bool doRopes, doShoving, doFlavour, doVertex, doBuffon;
113 
114  // Flavour change with Rope Hadronization.
115  FlavourRope flavourRope;
116 
117  // The generator class for Bose-Einstein effects.
118  BoseEinstein boseEinstein;
119 
120  // Classes for flavour, pT and z generation.
121  StringFlav flavSel;
122  StringPT pTSel;
123  StringZ zSel;
124 
125  // Class for colour tracing.
126  ColourTracing colTrace;
127 
128  // Junction splitting class.
129  JunctionSplitting junctionSplitting;
130 
131  // The RHadrons class is used to fragment off and decay R-hadrons.
132  RHadrons* rHadronsPtr;
133 
134  // Special class for Hidden-Valley hadronization. Not always used.
135  HiddenValleyFragmentation hiddenvalleyFrag;
136  bool useHiddenValley;
137 
138  // Special case: colour-octet onium decays, to be done initially.
139  bool decayOctetOnia(Event& event);
140 
141  // Trace colour flow in the event to form colour singlet subsystems.
142  // Option to keep junctions, needed for rope hadronization.
143  bool findSinglets(Event& event, bool keepJunctions = false);
144 
145  // Extract rapidity pairs.
146  vector< vector< pair<double,double> > > rapidityPairs(Event& event);
147 
148  // Calculate the rapidity for string ends, protected against too large y.
149  double yMax(Particle pIn, double mTiny) {
150  double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, pIn.mT()) );
151  return (pIn.pz() > 0) ? temp : -temp; }
152 
153 };
154 
155 //==========================================================================
156 
157 } // end namespace Pythia8
158 
159 #endif // Pythia8_HadronLevel_H
Definition: AgUStep.h:26