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) 2020 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/DeuteronProduction.h"
16 #include "Pythia8/Event.h"
17 #include "Pythia8/FragmentationFlavZpT.h"
18 #include "Pythia8/FragmentationSystems.h"
19 #include "Pythia8/HadronWidths.h"
20 #include "Pythia8/HiddenValleyFragmentation.h"
21 #include "Pythia8/Info.h"
22 #include "Pythia8/JunctionSplitting.h"
23 #include "Pythia8/LowEnergyProcess.h"
24 #include "Pythia8/LowEnergySigma.h"
25 #include "Pythia8/MiniStringFragmentation.h"
26 #include "Pythia8/NucleonExcitations.h"
27 #include "Pythia8/ParticleData.h"
28 #include "Pythia8/ParticleDecays.h"
29 #include "Pythia8/PartonVertex.h"
30 #include "Pythia8/PhysicsBase.h"
31 #include "Pythia8/PythiaStdlib.h"
32 #include "Pythia8/RHadrons.h"
33 #include "Pythia8/Settings.h"
34 #include "Pythia8/StringFragmentation.h"
35 #include "Pythia8/TimeShower.h"
36 #include "Pythia8/UserHooks.h"
37 
38 namespace Pythia8 {
39 
40 //==========================================================================
41 
42 // The HadronLevel class contains the top-level routines to generate
43 // the transition from the partonic to the hadronic stage of an event.
44 
45 class HadronLevel : public PhysicsBase {
46 
47 public:
48 
49  // Constructor.
50  HadronLevel() = default;
51 
52  // Initialize HadronLevel classes as required.
53  bool init( TimeShowerPtr timesDecPtr, RHadrons* rHadronsPtrIn,
54  DecayHandlerPtr decayHandlePtr, vector<int> handledParticles,
55  StringIntPtr stringInteractionsPtrIn, PartonVertexPtr partonVertexPtrIn);
56 
57  // Get pointer to StringFlav instance (needed by BeamParticle).
58  StringFlav* getStringFlavPtr() {return &flavSel;}
59 
60  // Generate the next event.
61  bool next(Event& event);
62 
63  // Special routine to allow more decays if on/off switches changed.
64  bool moreDecays(Event& event);
65 
66  // Prepare and pick process for a low-energy hadron-hadron scattering.
67  bool initLowEnergyProcesses();
68  int pickLowEnergyProcess(int idA, int idB, double eCM, double mA, double mB);
69 
70  // Special routine to do a low-energy hadron-hadron scattering.
71  bool doLowEnergyProcess(int i1, int i2, int type, Event& event) {
72  if (!lowEnergyProcess.collide( i1, i2, type, event)) {
73  infoPtr->errorMsg("Error in HadronLevel::doLowEnergyProcess: "
74  "Low energy collision failed");
75  return false;
76  }
77  return true;
78  }
79 
80  // Routine to allow user access to low-energy cross sections.
81  double getLowEnergySigma( int idA, int idB, double eCM, double mA,
82  double mB, int type = 0) {
83  return lowEnergySigma.sigmaPartial( idA, idB, eCM, mA, mB, type); }
84 
85  // Give access to b slope in elastic and diffractive interactions.
86  double getLowEnergySlope( int idA, int idB, double eCM, double mA,
87  double mB, int type = 2) {
88  return lowEnergyProcess.bSlope( idA, idB, eCM, mA, mB, type); }
89 
90 protected:
91 
92  virtual void onInitInfoPtr() override{
93  registerSubObject(flavSel);
94  registerSubObject(pTSel);
95  registerSubObject(zSel);
96  registerSubObject(stringFrag);
97  registerSubObject(ministringFrag);
98  registerSubObject(decays);
99  registerSubObject(lowEnergyProcess);
100  registerSubObject(lowEnergySigma);
101  registerSubObject(nucleonExcitations);
102  registerSubObject(boseEinstein);
103  registerSubObject(hiddenvalleyFrag);
104  registerSubObject(junctionSplitting);
105  registerSubObject(deuteronProd);
106  }
107 
108 private:
109 
110  // Constants: could only be changed in the code itself.
111  static const double MTINY;
112 
113  // Initialization data, read from Settings.
114  bool doHadronize{}, doDecay{}, doPartonVertex{}, doBoseEinstein{},
115  doDeuteronProd{}, allowRH{}, closePacking{}, doNonPertAll{};
116  double mStringMin{}, eNormJunction{}, widthSepBE{}, widthSepRescatter{};
117  vector<int> nonPertProc{};
118 
119  // Configuration of colour-singlet systems.
120  ColConfig colConfig;
121 
122  // Colour and mass information.
123  vector<int> iParton{}, iJunLegA{}, iJunLegB{}, iJunLegC{},
124  iAntiLegA{}, iAntiLegB{}, iAntiLegC{}, iGluLeg{};
125  vector<double> m2Pair{};
126 
127  // The generator class for normal string fragmentation.
128  StringFragmentation stringFrag;
129 
130  // The generator class for special low-mass string fragmentation.
131  MiniStringFragmentation ministringFrag;
132 
133  // The generator class for normal decays.
134  ParticleDecays decays;
135 
136  // The generator class for Bose-Einstein effects.
137  BoseEinstein boseEinstein;
138 
139  // The generator class for deuteron production.
140  DeuteronProduction deuteronProd;
141 
142  // Classes for flavour, pT and z generation.
143  StringFlav flavSel;
144  StringPT pTSel;
145  StringZ zSel;
146 
147  // Class for colour tracing.
148  ColourTracing colTrace;
149 
150  // Junction splitting class.
151  JunctionSplitting junctionSplitting;
152 
153  // The RHadrons class is used to fragment off and decay R-hadrons.
154  RHadrons* rHadronsPtr;
155 
156  // Special class for Hidden-Valley hadronization. Not always used.
157  HiddenValleyFragmentation hiddenvalleyFrag;
158  bool useHiddenValley{};
159 
160  // Special case: colour-octet onium decays, to be done initially.
161  bool decayOctetOnia(Event& event);
162 
163  // Trace colour flow in the event to form colour singlet subsystems.
164  // Option to keep junctions, needed for rope hadronization.
165  bool findSinglets(Event& event, bool keepJunctions = false);
166 
167  // Class to displace hadron vertices from parton impact-parameter picture.
168  PartonVertexPtr partonVertexPtr;
169 
170  // Hadronic rescattering.
171  class PriorityNode;
172  bool doRescatter{}, scatterManyTimes{}, scatterQuickCheck{},
173  scatterNeighbours{}, delayRegeneration{};
174  double b2Max, tauRegeneration{};
175  void queueDecResc(Event& event, int iStart,
176  priority_queue<HadronLevel::PriorityNode>& queue);
177  int boostDir;
178  double boost;
179  bool doBoost;
180  bool useVelocityFrame;
181 
182  // The generator class for low-energy hadron-hadron collisions.
183  LowEnergyProcess lowEnergyProcess;
184  int impactModel{};
185  double impactOpacity{};
186 
187  // Cross sections for low-energy processes.
188  LowEnergySigma lowEnergySigma;
189 
190  // Nucleon excitations data.
191  NucleonExcitations nucleonExcitations = {};
192 
193  // Class for event geometry for Rope Hadronization. Production vertices.
194  StringRepPtr stringRepulsionPtr;
195  FragModPtr fragmentationModifierPtr;
196 
197  // Extract rapidity pairs.
198  vector< vector< pair<double,double> > > rapidityPairs(Event& event);
199 
200  // Calculate the rapidity for string ends, protected against too large y.
201  double yMax(Particle pIn, double mTiny) {
202  double temp = log( ( pIn.e() + abs(pIn.pz()) ) / max( mTiny, pIn.mT()) );
203  return (pIn.pz() > 0) ? temp : -temp; }
204 
205 };
206 
207 //==========================================================================
208 
209 } // end namespace Pythia8
210 
211 #endif // Pythia8_HadronLevel_H
Definition: AgUStep.h:26