ROOT
6.06/09
Reference Guide
ROOT Home Page
Main Page
Related Pages
User's Classes
Namespaces
All Classes
Files
Release Notes
File List
File Members
montecarlo
pythia8
src
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
23
ClassImp
(
TPythia8Decayer
)
24
25
////////////////////////////////////////////////////////////////////////////////
26
///constructor
27
28
TPythia8Decayer
::
TPythia8Decayer
():
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
39
void
TPythia8Decayer::Init
()
40
{
41
}
42
43
////////////////////////////////////////////////////////////////////////////////
44
/// Decay a single particle
45
46
void
TPythia8Decayer::Decay
(
Int_t
pdg,
TLorentzVector
* p)
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
59
Int_t
TPythia8Decayer::ImportParticles
(
TClonesArray
*particles)
60
{
61
return
(
fPythia8
->
ImportParticles
(particles,
"All"
));
62
}
63
64
////////////////////////////////////////////////////////////////////////////////
65
/// Set forced decay mode
66
67
void
TPythia8Decayer::SetForceDecay
(
Int_t
/*type*/
)
68
{
69
printf
(
"SetForceDecay not yet implemented !\n"
);
70
}
71
////////////////////////////////////////////////////////////////////////////////
72
/// ForceDecay not yet implemented
73
74
void
TPythia8Decayer::ForceDecay
()
75
{
76
printf
(
"ForceDecay not yet implemented !\n"
);
77
}
78
////////////////////////////////////////////////////////////////////////////////
79
80
Float_t
TPythia8Decayer::GetPartialBranchingRatio
(
Int_t
/*ipart*/
)
81
{
82
return
0.0;
83
}
84
////////////////////////////////////////////////////////////////////////////////
85
///return lifetime in seconds of teh particle with PDG number pdg
86
87
Float_t
TPythia8Decayer::GetLifetime
(
Int_t
pdg)
88
{
89
return
(
fPythia8
->
Pythia8
()->particleData.tau0(pdg) * 3.3333e-12) ;
90
}
91
92
////////////////////////////////////////////////////////////////////////////////
93
///to read a decay table (not yet implemented)
94
95
void
TPythia8Decayer::ReadDecayTable
()
96
{
97
}
98
99
100
////////////////////////////////////////////////////////////////////////////////
101
/// Append a particle to the stack
102
103
void
TPythia8Decayer::AppendParticle
(
Int_t
pdg,
TLorentzVector
* p)
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
112
void
TPythia8Decayer::ClearEvent
()
113
{
114
fPythia8
->
Pythia8
()->event.clear();
115
}
116
TPythia8::ImportParticles
virtual Int_t ImportParticles(TClonesArray *particles, Option_t *option="")
Import particles from Pythia stack.
Definition:
TPythia8.cxx:191
TPythia8Decayer::ReadDecayTable
virtual void ReadDecayTable()
to read a decay table (not yet implemented)
Definition:
TPythia8Decayer.cxx:95
Float_t
float Float_t
Definition:
RtypesCore.h:53
TPythia8Decayer::fDebug
Int_t fDebug
Definition:
TPythia8Decayer.h:35
Int_t
int Int_t
Definition:
RtypesCore.h:41
TPythia8Decayer::Decay
virtual void Decay(Int_t pdg, TLorentzVector *p)
Decay a single particle.
Definition:
TPythia8Decayer.cxx:46
TPythia8Decayer::ClearEvent
void ClearEvent()
Clear the event stack.
Definition:
TPythia8Decayer.cxx:112
TLorentzVector::M
Double_t M() const
Definition:
TLorentzVector.h:498
TPythia8Decayer.h
TPythia8Decayer::GetPartialBranchingRatio
virtual Float_t GetPartialBranchingRatio(Int_t ipart)
Definition:
TPythia8Decayer.cxx:80
TPythia8Decayer::Init
virtual void Init()
Initialize the decayer.
Definition:
TPythia8Decayer.cxx:39
TPythia8Decayer::SetForceDecay
virtual void SetForceDecay(Int_t type)
Set forced decay mode.
Definition:
TPythia8Decayer.cxx:67
TPythia8Decayer::GetLifetime
virtual Float_t GetLifetime(Int_t kf)
return lifetime in seconds of teh particle with PDG number pdg
Definition:
TPythia8Decayer.cxx:87
TPythia8Decayer::ImportParticles
virtual Int_t ImportParticles(TClonesArray *particles)
import the decay products into particles array
Definition:
TPythia8Decayer.cxx:59
TPythia8Decayer::AppendParticle
void AppendParticle(Int_t pdg, TLorentzVector *p)
Append a particle to the stack.
Definition:
TPythia8Decayer.cxx:103
TPythia8Decayer::ForceDecay
virtual void ForceDecay()
ForceDecay not yet implemented.
Definition:
TPythia8Decayer.cxx:74
TLorentzVector
Definition:
TLorentzVector.h:38
TPythia8Decayer::fPythia8
TPythia8 * fPythia8
Definition:
TPythia8Decayer.h:34
ClassImp
ClassImp(TPythia8Decayer) TPythia8Decayer
constructor
Definition:
TPythia8Decayer.cxx:23
TPythia8::EventListing
void EventListing() const
Event listing.
Definition:
TPythia8.cxx:361
TLorentzVector::Px
Double_t Px() const
Definition:
TLorentzVector.h:287
TLorentzVector::Py
Double_t Py() const
Definition:
TLorentzVector.h:288
TPythia8.h
TPythia8
Definition:
TPythia8.h:76
printf
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
TLorentzVector::E
Double_t E() const
Definition:
TLorentzVector.h:291
TPythia8Decayer
Definition:
TPythia8Decayer.h:16
TPythia8::Pythia8
Pythia8::Pythia * Pythia8()
Definition:
TPythia8.h:89
TClonesArray
An array of clone (identical) objects.
Definition:
TClonesArray.h:32
TLorentzVector::Pz
Double_t Pz() const
Definition:
TLorentzVector.h:289
TLorentzVector.h
kTRUE
const Bool_t kTRUE
Definition:
Rtypes.h:91