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