ROOT  6.06/09
Reference Guide
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) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
6  * *
7  * Author: The ALICE Off-line Project. *
8  * Contributors are mentioned in the code where appropriate. *
9  * *
10  * Permission to use, copy, modify and distribute this software and its *
11  * documentation strictly for non-commercial purposes is hereby granted *
12  * without fee, provided that the above copyright notice appears in all *
13  * copies and that both the copyright notice and this permission notice *
14  * appear in the supporting documentation. The authors make no claims *
15  * about the suitability of this software for any purpose. It is *
16  * provided "as is" without express or implied warranty. *
17  **************************************************************************/
18 
19 #include "TLorentzVector.h"
20 #include "TPythia8.h"
21 #include "TPythia8Decayer.h"
22 
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();
49  AppendParticle(pdg, p);
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 
60 {
61  return (fPythia8->ImportParticles(particles, "All"));
62 }
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 
81 {
82  return 0.0;
83 }
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 
96 {
97 }
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 
virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="")
Import particles from Pythia stack.
Definition: TPythia8.cxx:191
virtual void ReadDecayTable()
to read a decay table (not yet implemented)
float Float_t
Definition: RtypesCore.h:53
int Int_t
Definition: RtypesCore.h:41
virtual void Decay(Int_t pdg, TLorentzVector *p)
Decay a single particle.
void ClearEvent()
Clear the event stack.
Double_t M() const
virtual Float_t GetPartialBranchingRatio(Int_t ipart)
virtual void Init()
Initialize the decayer.
virtual void SetForceDecay(Int_t type)
Set forced decay mode.
virtual Float_t GetLifetime(Int_t kf)
return lifetime in seconds of teh particle with PDG number pdg
virtual Int_t ImportParticles(TClonesArray *particles)
import the decay products into particles array
void AppendParticle(Int_t pdg, TLorentzVector *p)
Append a particle to the stack.
virtual void ForceDecay()
ForceDecay not yet implemented.
TPythia8 * fPythia8
ClassImp(TPythia8Decayer) TPythia8Decayer
constructor
void EventListing() const
Event listing.
Definition: TPythia8.cxx:361
Double_t Px() const
Double_t Py() const
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Double_t E() const
Pythia8::Pythia * Pythia8()
Definition: TPythia8.h:89
An array of clone (identical) objects.
Definition: TClonesArray.h:32
Double_t Pz() const
const Bool_t kTRUE
Definition: Rtypes.h:91