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