Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TPythia8Decayer.cxx
Go to the documentation of this file.
1// @(#)root/pythia8:$Name$:$Id$
2// Author: Andreas Morsch 04/07/2008
3
4/*************************************************************************
5 * Copyright (C) 1995-2023, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TPythia8Decayer
13 \ingroup pythia8
14
15This class implements the TVirtualMCDecayer interface using TPythia8.
16
17Author: Andreas Morsch 04/07/2008
18*/
19
20#include "TLorentzVector.h"
21#include "TPythia8.h"
22#include "TPythia8Decayer.h"
23
24
25////////////////////////////////////////////////////////////////////////////////
26///constructor
27
29 fPythia8(new TPythia8()),
30 fDebug(0)
31{
32 fPythia8->Pythia8()->readString("SoftQCD:elastic = on");
33 fPythia8->Pythia8()->init();
34}
35
36////////////////////////////////////////////////////////////////////////////////
37/// Initialize the decayer
38
40{
41}
42
43////////////////////////////////////////////////////////////////////////////////
44/// Decay a single particle
45
47{
48 ClearEvent();
50 Int_t idPart = fPythia8->Pythia8()->event[0].id();
51 fPythia8->Pythia8()->particleData.mayDecay(idPart,kTRUE);
52 fPythia8->Pythia8()->moreDecays();
53 if (fDebug > 0) fPythia8->EventListing();
54}
55
56////////////////////////////////////////////////////////////////////////////////
57///import the decay products into particles array
58
63
64////////////////////////////////////////////////////////////////////////////////
65/// Set forced decay mode
66
68{
69 printf("SetForceDecay not yet implemented !\n");
70}
71////////////////////////////////////////////////////////////////////////////////
72/// ForceDecay not yet implemented
73
75{
76 printf("ForceDecay not yet implemented !\n");
77}
78////////////////////////////////////////////////////////////////////////////////
79
84////////////////////////////////////////////////////////////////////////////////
85///return lifetime in seconds of teh particle with PDG number pdg
86
88{
89 return (fPythia8->Pythia8()->particleData.tau0(pdg) * 3.3333e-12) ;
90}
91
92////////////////////////////////////////////////////////////////////////////////
93///to read a decay table (not yet implemented)
94
98
99
100////////////////////////////////////////////////////////////////////////////////
101/// Append a particle to the stack
102
104{
105 fPythia8->Pythia8()->event.append(pdg, 11, 0, 0, p->Px(), p->Py(), p->Pz(), p->E(), p->M());
106}
107
108
109////////////////////////////////////////////////////////////////////////////////
110/// Clear the event stack
111
113{
114 fPythia8->Pythia8()->event.clear();
115}
116
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
An array of clone (identical) objects.
TPythia8 * fPythia8
Float_t GetPartialBranchingRatio(Int_t ipart) override
Get the partial branching ratio for a particle of type IPART (a PDG code).
void ForceDecay() override
ForceDecay not yet implemented.
TPythia8Decayer()
constructor
void ReadDecayTable() override
to read a decay table (not yet implemented)
void Init() override
Initialize the decayer.
Float_t GetLifetime(Int_t kf) override
return lifetime in seconds of teh particle with PDG number pdg
void AppendParticle(Int_t pdg, TLorentzVector *p)
Append a particle to the stack.
void Decay(Int_t pdg, TLorentzVector *p) override
Decay a single particle.
void ClearEvent()
Clear the event stack.
void SetForceDecay(Int_t type) override
Set forced decay mode.
Int_t ImportParticles(TClonesArray *particles) override
import the decay products into particles array
TPythia8 is an interface class to C++ version of Pythia 8.1 event generators, written by T....
Definition TPythia8.h:77
Pythia8::Pythia * Pythia8()
Definition TPythia8.h:89
Int_t ImportParticles(TClonesArray *particles, Option_t *option="") override
Import particles from Pythia stack.
Definition TPythia8.cxx:193
void EventListing() const
Event listing.
Definition TPythia8.cxx:363