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