StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
ParticleDecays.h
1 // ParticleDecays.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 classes to perform a particle decay.
7 // DecayHandler: base class for external handling of decays.
8 // ParticleDecays: decay a particle.
9 
10 #ifndef Pythia8_ParticleDecays_H
11 #define Pythia8_ParticleDecays_H
12 
13 #include "Pythia8/Basics.h"
14 #include "Pythia8/Event.h"
15 #include "Pythia8/FragmentationFlavZpT.h"
16 #include "Pythia8/Info.h"
17 #include "Pythia8/ParticleData.h"
18 #include "Pythia8/PythiaStdlib.h"
19 #include "Pythia8/Settings.h"
20 #include "Pythia8/TimeShower.h"
21 #include "Pythia8/TauDecays.h"
22 
23 namespace Pythia8 {
24 
25 //==========================================================================
26 
27 // DecayHandler is base class for the external handling of decays.
28 // There is only one pure virtual method, that should do the decay.
29 
30 class DecayHandler {
31 
32 public:
33 
34  // Destructor.
35  virtual ~DecayHandler() {}
36 
37  // A virtual method, wherein the derived class method does a decay.
38  // Usage: decay( idProd, mProd, pProd, iDec, event).
39  virtual bool decay(vector<int>& , vector<double>& , vector<Vec4>& ,
40  int , const Event& ) {return false;}
41 
42  // A virtual method, to do sequential decay chains.
43  // Usage: decay( idProd, motherProd, mProd, pProd, iDec, event).
44  virtual bool chainDecay(vector<int>& , vector<int>& , vector<double>& ,
45  vector<Vec4>& , int , const Event& ) {return false;}
46 
47 };
48 
49 //==========================================================================
50 
51 // The ParticleDecays class contains the routines to decay a particle.
52 
53 class ParticleDecays {
54 
55 public:
56 
57  // Constructor.
58  ParticleDecays() {}
59 
60  // Initialize: store pointers and find settings
61  void init(Info* infoPtrIn, Settings& settings,
62  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
63  Couplings* couplingsPtrIn, TimeShower* timesDecPtrIn,
64  StringFlav* flavSelPtrIn, DecayHandler* decayHandlePtrIn,
65  vector<int> handledParticles);
66 
67  // Perform a decay of a single particle.
68  bool decay(int iDec, Event& event);
69 
70  // Did decay result in new partons to hadronize?
71  bool moreToDo() const {return hasPartons && keepPartons;}
72 
73 private:
74 
75  // Constants: could only be changed in the code itself.
76  static const int NTRYDECAY, NTRYPICK, NTRYMEWT, NTRYDALITZ;
77  static const double MSAFEDALITZ, WTCORRECTION[11];
78 
79  // Pointer to various information on the generation.
80  Info* infoPtr;
81 
82  // Pointer to the particle data table.
83  ParticleData* particleDataPtr;
84 
85  // Pointer to the random number generator.
86  Rndm* rndmPtr;
87 
88  // Pointers to Standard Model couplings.
89  Couplings* couplingsPtr;
90 
91  // Pointers to timelike showers, for decays to partons (e.g. Upsilon).
92  TimeShower* timesDecPtr;
93 
94  // Pointer to class for flavour generation; needed when to pick hadrons.
95  StringFlav* flavSelPtr;
96 
97  // Pointer to a handler of external decays.
98  DecayHandler* decayHandlePtr;
99 
100  // Initialization data, read from Settings.
101  bool limitTau0, limitTau, limitRadius, limitCylinder, limitDecay,
102  mixB, doFSRinDecays, doGammaRad;
103  int tauMode;
104  double mSafety, tau0Max, tauMax, rMax, xyMax, zMax, xBdMix, xBsMix,
105  sigmaSoft, multIncrease, multIncreaseWeak, multRefMass, multGoffset,
106  colRearrange, stopMass, sRhoDal, wRhoDal;
107 
108  // Multiplicity. Decay products positions and masses.
109  bool hasPartons, keepPartons;
110  int idDec, meMode, mult;
111  double scale;
112  vector<int> iProd, idProd, motherProd, cols, acols, idPartons;
113  vector<double> mProd, mInv, rndmOrd;
114  vector<Vec4> pInv, pProd;
115  vector<FlavContainer> flavEnds;
116 
117  // Pointer to particle data for currently decaying particle
118  ParticleDataEntry* decDataPtr;
119 
120  // Tau particle decayer.
121  TauDecays tauDecayer;
122 
123  // Check whether a decay is allowed, given the upcoming decay vertex.
124  bool checkVertex(Particle& decayer);
125 
126  // Check for oscillations B0 <-> B0bar or B_s0 <-> B_s0bar.
127  bool oscillateB(Particle& decayer);
128 
129  // Do a one-body decay.
130  bool oneBody(Event& event);
131 
132  // Do a two-body decay;
133  bool twoBody(Event& event);
134 
135  // Do a three-body decay;
136  bool threeBody(Event& event);
137 
138  // Do a multibody decay using the M-generator algorithm.
139  bool mGenerator(Event& event);
140 
141  // Select mass of lepton pair in a Dalitz decay.
142  bool dalitzMass();
143 
144  // Do kinematics of gamma* -> l- l+ in Dalitz decay.
145  bool dalitzKinematics(Event& event);
146 
147  // Translate a partonic content into a set of actual hadrons.
148  bool pickHadrons();
149 
150  // Set colour flow and scale in a decay explicitly to partons.
151  bool setColours(Event& event);
152 
153 };
154 
155 //==========================================================================
156 
157 } // end namespace Pythia8
158 
159 #endif // Pythia8_ParticleDecays_H
Definition: AgUStep.h:26