]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PYTHIA8/AliTPythia8.cxx
Added plugin->SetProofParameter(const char *name, const char *value) to set special...
[u/mrichter/AliRoot.git] / PYTHIA8 / AliTPythia8.cxx
CommitLineData
b584e2f5 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16// Author: Andreas Morsch 27/10/2007
17////////////////////////////////////////////////////////////////////////////////
18// //
19// TPythia8 //
20// //
21// TPythia is an interface class to C++ version of Pythia 8.1 //
22// event generators, written by T.Sjostrand. //
23// //
24// The user is assumed to be familiar with the Pythia package. //
25// This class includes only a basic interface to Pythia8. Because Pythia8 is //
26// also written in C++, its functions/classes can be called directly from a //
27// compiled C++ script. //
28// To call Pythia functions not available in this interface a dictionary must //
29// be generated. //
30// see $ROOTSYS/tutorials/pythia/pythia8.C for an example of use from CINT. //
31// //
32////////////////////////////////////////////////////////////////////////////////
33/*
34*------------------------------------------------------------------------------------*
35 | |
36 | *------------------------------------------------------------------------------* |
37 | | | |
38 | | | |
39 | | PPP Y Y TTTTT H H III A Welcome to the Lund Monte Carlo! | |
40 | | P P Y Y T H H I A A This is PYTHIA version 8.100 | |
41 | | PPP Y T HHHHH I AAAAA Last date of change: 20 Oct 2007 | |
42 | | P Y T H H I A A | |
43 | | P Y T H H III A A Now is 27 Oct 2007 at 18:26:53 | |
44 | | | |
45 | | Main author: Torbjorn Sjostrand; CERN/PH, CH-1211 Geneva, Switzerland, | |
46 | | and Department of Theoretical Physics, Lund University, Lund, Sweden; | |
47 | | phone: + 41 - 22 - 767 82 27; e-mail: torbjorn@thep.lu.se | |
48 | | Author: Stephen Mrenna; Computing Division, Simulations Group, | |
49 | | Fermi National Accelerator Laboratory, MS 234, Batavia, IL 60510, USA; | |
50 | | phone: + 1 - 630 - 840 - 2556; e-mail: mrenna@fnal.gov | |
51 | | Author: Peter Skands; CERN/PH, CH-1211 Geneva, Switzerland, | |
52 | | and Theoretical Physics Department, | |
53 | | Fermi National Accelerator Laboratory, MS 106, Batavia, IL 60510, USA; | |
54 | | phone: + 41 - 22 - 767 24 59; e-mail: skands@fnal.gov | |
55 | | | |
56 | | The main program reference is the 'Brief Introduction to PYTHIA 8.1', | |
57 | | T. Sjostrand, S. Mrenna and P. Skands, arXiv:0710.3820 | |
58 | | | |
59 | | The main physics reference is the 'PYTHIA 6.4 Physics and Manual', | |
60 | | T. Sjostrand, S. Mrenna and P. Skands, JHEP05 (2006) 026 [hep-ph/0603175]. | |
61 | | | |
62 | | An archive of program versions and documentation is found on the web: | |
63 | | http://www.thep.lu.se/~torbjorn/Pythia.html | |
64 | | | |
65 | | This program is released under the GNU General Public Licence version 2. | |
66 | | Please respect the MCnet Guidelines for Event Generator Authors and Users. | |
67 | | | |
68 | | Disclaimer: this program comes without any guarantees. | |
69 | | Beware of errors and use common sense when interpreting results. | |
70 | | | |
71 | | Copyright (C) 2007 Torbjorn Sjostrand | |
72 | | | |
73 | | | |
74 | *------------------------------------------------------------------------------* |
75 | |
76 *------------------------------------------------------------------------------------*
77*/
78
79#include "AliTPythia8.h"
80
81#include "TClonesArray.h"
82#include "TParticle.h"
83#include "TDatabasePDG.h"
84#include "TLorentzVector.h"
85
86ClassImp(AliTPythia8)
87
88AliTPythia8* AliTPythia8::fgInstance = 0;
89
90//___________________________________________________________________________
91AliTPythia8::AliTPythia8():
92 TGenerator("AliTPythia8", "AliTPythia8"),
93 fPythia(0),
94 fNumberOfParticles(0)
95{
96 // Constructor
97 if (fgInstance)
98 Fatal("AliTPythia8", "There's already an instance of AliTPythia8");
99
100 delete fParticles; // was allocated as TObjArray in TGenerator
101
102 fParticles = new TClonesArray("TParticle",50);
103 fPythia = new Pythia8::Pythia();
104}
105
106//___________________________________________________________________________
107AliTPythia8::AliTPythia8(const char *xmlDir):
108 TGenerator("AliTPythia8", "AliTPythia8"),
109 fPythia(0),
110 fNumberOfParticles(0)
111{
112 // Constructor with an xmlDir (eg "../xmldoc"
113 if (fgInstance)
114 Fatal("AliTPythia8", "There's already an instance of AliTPythia8");
115
116 delete fParticles; // was allocated as TObjArray in TGenerator
117
118 fParticles = new TClonesArray("TParticle",50);
119 fPythia = new Pythia8::Pythia(xmlDir);
120}
121
122//___________________________________________________________________________
123AliTPythia8::~AliTPythia8()
124{
125 // Destructor
126 if (fParticles) {
127 fParticles->Delete();
128 delete fParticles;
129 fParticles = 0;
130 }
131 delete fPythia;
132}
133
134//___________________________________________________________________________
135AliTPythia8* AliTPythia8::Instance()
136{
137 // Return an instance of AliTPythia8
138 return fgInstance ? fgInstance : (fgInstance = new AliTPythia8()) ;
139}
140
141//___________________________________________________________________________
142Bool_t AliTPythia8::Initialize(Int_t idAin, Int_t idBin, Double_t ecms)
143{
144 // Initialization
145 AddParticlesToPdgDataBase();
146 return fPythia->init(idAin, idBin, ecms);
147}
148
149//___________________________________________________________________________
150void AliTPythia8::GenerateEvent()
151{
152 // Generate the next event
153 fPythia->next();
154 fNumberOfParticles = fPythia->event.size() - 1;
155 ImportParticles();
156}
157
158//___________________________________________________________________________
159Int_t AliTPythia8::ImportParticles(TClonesArray *particles, Option_t *option)
160{
161 // Import particles from Pythia stack
162 if (particles == 0) return 0;
163 TClonesArray &clonesParticles = *particles;
164 clonesParticles.Clear();
165 Int_t nparts=0;
166 Int_t i;
167 fNumberOfParticles = fPythia->event.size() - 1;
168
169 if (!strcmp(option,"") || !strcmp(option,"Final")) {
170 for (i = 0; i <= fNumberOfParticles; i++) {
171 if (fPythia->event[i].id() == 90) continue;
172 if (fPythia->event[i].isFinal()) {
173 new(clonesParticles[nparts]) TParticle(
174 fPythia->event[i].id(),
175 fPythia->event[i].isFinal(),
176 fPythia->event[i].mother1() - 1,
177 fPythia->event[i].mother2() - 1,
178 fPythia->event[i].daughter1() - 1,
179 fPythia->event[i].daughter2() - 1,
180 fPythia->event[i].px(), // [GeV/c]
181 fPythia->event[i].py(), // [GeV/c]
182 fPythia->event[i].pz(), // [GeV/c]
183 fPythia->event[i].e(), // [GeV]
184 fPythia->event[i].xProd(), // [mm]
185 fPythia->event[i].yProd(), // [mm]
186 fPythia->event[i].zProd(), // [mm]
187 fPythia->event[i].tProd()); // [mm/c]
188 nparts++;
189 } // final state partice
190 } // particle loop
191 } else if (!strcmp(option,"All")) {
192 for (i = 0; i <= fNumberOfParticles; i++) {
193 if (fPythia->event[i].id() == 90) continue;
194 new(clonesParticles[nparts]) TParticle(
195 fPythia->event[i].id(),
196 fPythia->event[i].isFinal(),
197 fPythia->event[i].mother1() - 1,
198 fPythia->event[i].mother2() - 1,
199 fPythia->event[i].daughter1() - 1,
200 fPythia->event[i].daughter2() - 1,
201 fPythia->event[i].px(), // [GeV/c]
202 fPythia->event[i].py(), // [GeV/c]
203 fPythia->event[i].pz(), // [GeV/c]
204 fPythia->event[i].e(), // [GeV]
205 fPythia->event[i].xProd(), // [mm]
206 fPythia->event[i].yProd(), // [mm]
207 fPythia->event[i].zProd(), // [mm]
208 fPythia->event[i].tProd()); // [mm/c]
209 nparts++;
210 } // particle loop
211 }
212 return nparts;
213}
214
215//___________________________________________________________________________
216TObjArray* AliTPythia8::ImportParticles(Option_t* /* option */)
217{
218 // Import particles from Pythia stack
219 fParticles->Clear();
220 Int_t numpart = fPythia->event.size() - 1;
221 TClonesArray &a = *((TClonesArray*)fParticles);
222 for (Int_t i = 1; i <= numpart; i++) {
223 new(a[i]) TParticle(
224 fPythia->event[i].id(),
225 fPythia->event[i].isFinal(),
226 fPythia->event[i].mother1() - 1,
227 fPythia->event[i].mother2() - 1,
228 fPythia->event[i].daughter1() - 1,
229 fPythia->event[i].daughter2() - 1,
230 fPythia->event[i].px(), // [GeV/c]
231 fPythia->event[i].py(), // [GeV/c]
232 fPythia->event[i].pz(), // [GeV/c]
233 fPythia->event[i].e(), // [GeV]
234 fPythia->event[i].xProd(), // [mm]
235 fPythia->event[i].yProd(), // [mm]
236 fPythia->event[i].zProd(), // [mm]
237 fPythia->event[i].tProd()); // [mm/c]
238 }
239 return fParticles;
240}
241
242//___________________________________________________________________________
243Int_t AliTPythia8::GetN() const
244{
245 // Initialization
246 return (fPythia->event.size() - 1);
247}
248
249//___________________________________________________________________________
250void AliTPythia8::ReadString(const char* string) const
251{
252 // Configuration
253 fPythia->readString(string);
254}
255
256//___________________________________________________________________________
257void AliTPythia8::ReadConfigFile(const char* string) const
258{
259 // Configuration
260 fPythia->readFile(string);
261}
262
263//___________________________________________________________________________
264void AliTPythia8::PrintStatistics() const
265{
266 // Print end of run statistics
267 fPythia->statistics();
268}
269
270//___________________________________________________________________________
271void AliTPythia8::EventListing() const
272{
273 // Event listing
274 fPythia->event.list();
275}
276
277//___________________________________________________________________________
c014d45a 278void AliTPythia8::AddParticlesToPdgDataBase() const
b584e2f5 279{
280 // Add some pythia specific particle code to the data base
281
282 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
283 pdgDB->AddParticle("string","string", 0, kTRUE,
284 0, 0, "QCD string", 90);
285 pdgDB->AddParticle("rho_diff0", "rho_diff0", 0, kTRUE,
286 0, 0, "QCD diffr. state", 9900110);
287 pdgDB->AddParticle("pi_diffr+", "pi_diffr+", 0, kTRUE,
288 0, 1, "QCD diffr. state", 9900210);
289 pdgDB->AddParticle("omega_di", "omega_di", 0, kTRUE,
290 0, 0, "QCD diffr. state", 9900220);
291 pdgDB->AddParticle("phi_diff","phi_diff", 0, kTRUE,
292 0, 0, "QCD diffr. state", 9900330);
293 pdgDB->AddParticle("J/psi_di", "J/psi_di", 0, kTRUE,
294 0, 0, "QCD diffr. state", 9900440);
295 pdgDB->AddParticle("n_diffr0","n_diffr0",0,kTRUE,
296 0, 0, "QCD diffr. state", 9902110);
297 pdgDB->AddParticle("p_diffr+","p_diffr+", 0, kTRUE,
298 0, 1, "QCD diffr. state", 9902210);
299}
300