10 #ifndef Pythia8_ParticleDecays_H 11 #define Pythia8_ParticleDecays_H 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" 40 virtual bool decay(vector<int>& , vector<double>& , vector<Vec4>& ,
41 int ,
const Event& ) {
return false;}
45 virtual bool chainDecay(vector<int>& , vector<int>& , vector<double>& ,
46 vector<Vec4>& ,
int ,
const Event& ) {
return false;}
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() {}
69 void init(TimeShowerPtr timesDecPtrIn,
StringFlav* flavSelPtrIn,
70 DecayHandlerPtr decayHandlePtrIn, vector<int> handledParticles);
76 bool decayAll(
Event& event,
double minWidth = 0.);
79 bool moreToDo()
const {
return hasPartons && keepPartons;}
84 registerSubObject(tauDecayer); }
89 static const int NTRYDECAY, NTRYPICK, NTRYMEWT, NTRYDALITZ;
90 static const double MSAFEDALITZ, WTCORRECTION[11];
93 TimeShowerPtr timesDecPtr;
99 DecayHandlerPtr decayHandlePtr;
102 bool limitTau0, limitTau, limitRadius, limitCylinder, limitDecay,
103 mixB, doFSRinDecays, doGammaRad;
105 double mSafety, tau0Max, tauMax, rMax, xyMax, zMax, xBdMix, xBsMix,
106 sigmaSoft, multIncrease, multIncreaseWeak, multRefMass, multGoffset,
107 colRearrange, stopMass, sRhoDal, wRhoDal;
110 bool hasPartons, keepPartons;
111 int idDec, meMode, mult;
113 vector<int> iProd, idProd, motherProd, cols, acols, idPartons;
114 vector<double> mProd, mInv, rndmOrd;
115 vector<Vec4> pInv, pProd;
116 vector<FlavContainer> flavEnds;
125 bool checkVertex(
Particle& decayer);
131 bool oneBody(
Event& event);
134 bool twoBody(
Event& event);
137 bool threeBody(
Event& event);
140 bool mGenerator(
Event& event);
146 bool dalitzKinematics(
Event& event);
152 bool setColours(
Event& event);
virtual void onInitInfoPtr() override
Definition: ParticleDecays.h:83
virtual ~DecayHandler()
Destructor.
Definition: ParticleDecays.h:36
This class holds info on a single particle species.
Definition: ParticleData.h:125
Definition: PhysicsBase.h:27
The ParticleDecays class contains the routines to decay a particle.
Definition: ParticleDecays.h:54
The Event class holds all info on the generated event.
Definition: Event.h:379
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:79
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
ParticleDecays()
Constructor.
Definition: ParticleDecays.h:59
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:84
Definition: ParticleDecays.h:31