Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 /** \class TPythia8Decayer
20  \ingroup pythia8
21 
22 This class implements the TVirtualMCDecayer interface using TPythia8.
23 
24 Author: Andreas Morsch 04/07/2008
25 */
26 
27 #include "TLorentzVector.h"
28 #include "TPythia8.h"
29 #include "TPythia8Decayer.h"
30 
31 ClassImp(TPythia8Decayer);
32 
33 ////////////////////////////////////////////////////////////////////////////////
34 ///constructor
35 
36 TPythia8Decayer::TPythia8Decayer():
37  fPythia8(new TPythia8()),
38  fDebug(0)
39 {
40  fPythia8->Pythia8()->readString("SoftQCD:elastic = on");
41  fPythia8->Pythia8()->init();
42 }
43 
44 ////////////////////////////////////////////////////////////////////////////////
45 /// Initialize the decayer
46 
47 void TPythia8Decayer::Init()
48 {
49 }
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 /// Decay a single particle
53 
54 void TPythia8Decayer::Decay(Int_t pdg, TLorentzVector* p)
55 {
56  ClearEvent();
57  AppendParticle(pdg, p);
58  Int_t idPart = fPythia8->Pythia8()->event[0].id();
59  fPythia8->Pythia8()->particleData.mayDecay(idPart,kTRUE);
60  fPythia8->Pythia8()->moreDecays();
61  if (fDebug > 0) fPythia8->EventListing();
62 }
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 ///import the decay products into particles array
66 
67 Int_t TPythia8Decayer::ImportParticles(TClonesArray *particles)
68 {
69  return (fPythia8->ImportParticles(particles, "All"));
70 }
71 
72 ////////////////////////////////////////////////////////////////////////////////
73 /// Set forced decay mode
74 
75 void TPythia8Decayer::SetForceDecay(Int_t /*type*/)
76 {
77  printf("SetForceDecay not yet implemented !\n");
78 }
79 ////////////////////////////////////////////////////////////////////////////////
80 /// ForceDecay not yet implemented
81 
82 void TPythia8Decayer::ForceDecay()
83 {
84  printf("ForceDecay not yet implemented !\n");
85 }
86 ////////////////////////////////////////////////////////////////////////////////
87 
88 Float_t TPythia8Decayer::GetPartialBranchingRatio(Int_t /*ipart*/)
89 {
90  return 0.0;
91 }
92 ////////////////////////////////////////////////////////////////////////////////
93 ///return lifetime in seconds of teh particle with PDG number pdg
94 
95 Float_t TPythia8Decayer::GetLifetime(Int_t pdg)
96 {
97  return (fPythia8->Pythia8()->particleData.tau0(pdg) * 3.3333e-12) ;
98 }
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 ///to read a decay table (not yet implemented)
102 
103 void TPythia8Decayer::ReadDecayTable()
104 {
105 }
106 
107 
108 ////////////////////////////////////////////////////////////////////////////////
109 /// Append a particle to the stack
110 
111 void TPythia8Decayer::AppendParticle(Int_t pdg, TLorentzVector* p)
112 {
113  fPythia8->Pythia8()->event.append(pdg, 11, 0, 0, p->Px(), p->Py(), p->Pz(), p->E(), p->M());
114 }
115 
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 /// Clear the event stack
119 
120 void TPythia8Decayer::ClearEvent()
121 {
122  fPythia8->Pythia8()->event.clear();
123 }
124