Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TPythia8Decayer.h
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#ifndef TPYTHIA8DECAYER_H
13#define TPYTHIA8DECAYER_H
14
15#include "TVirtualMCDecayer.h"
16
17class TClonesArrray;
18class TLorentzVector;
19class TPythia8;
20
22public:
24 ~TPythia8Decayer() override{;}
25 void Init() override;
26 void Decay(Int_t pdg, TLorentzVector* p) override;
27 Int_t ImportParticles(TClonesArray *particles) override;
28 void SetForceDecay(Int_t type) override;
29 void ForceDecay() override;
31 Float_t GetLifetime(Int_t kf) override;
32 void ReadDecayTable() override;
33
34 virtual void SetDebugLevel(Int_t debug) {fDebug = debug;}
35protected:
37 void ClearEvent();
38private:
39 TPythia8* fPythia8; // Pointer to pythia8
40 Int_t fDebug; // Debug level
41
42 ClassDefOverride(TPythia8Decayer, 1) // Particle Decayer using Pythia8
43
44};
45#endif
46
47
48
49
50
51
52
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
An array of clone (identical) objects.
This class implements the TVirtualMCDecayer interface using TPythia8.
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)
~TPythia8Decayer() override
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.
virtual void SetDebugLevel(Int_t debug)
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
Abstract base class for particle decays.