]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPHIC/AliGenTPHIC.cxx
Default parameter for backward compatibility (Matthias)
[u/mrichter/AliRoot.git] / TPHIC / AliGenTPHIC.cxx
CommitLineData
0b5dd071 1/**************************************************************************
2 * Copyright(c) 1998-2003, 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
17/* $Id$ */
18
0b5dd071 19// Event generator of two-photon processes
20// in ultra-peripheral ion collisions.
21// 5 two-photon process are implemented, see comments to SetProcess().
22//%
23// The example of the generator initialization for the process
24// gamma gamma -> X in CaCa collisions at 7000 TeV/nucleon is the following:
25//%
26// AliGenTPHIC *gener = new AliGenTPHIC();
27// gener->SetProcess(1);
28// gener->SetBeamEnergy(3500.);
29// gener->SetBeamZ(20);
30// gener->SetBeamA(40);
31// gener->SetMggRange(70.,200.);
32// gener->SetYggRange(-5.,5.);
33// gener->SetLumFunName("lum_ca_70_200.dat");
34// gener->SetLumFunFlag(-1);
35// gener->Init();
36//%
37// The two-photon luminosity function needed for the process cross section
38// calculation is time-consuming, therefore it can be calculated once for a
39// selected two-photon energy and rapidity range and selected ion species,
40// saved into a file and then can be reused in further calculation.
41//%
42// The integral cross section of the process is calculated in each event
43// by the Monte Carlo integration procedure, the differential cross section
44// of each thrown event is assigned as a weight to each track of the event
45// which can be then used during the event analysis by retreiving this weight
46// from any of the event tracks.
47//
48// The manual of the fortran generator is available in
49// $ALICE_ROOT/TPHIC/TPHIC_doc.ps.gz
50// For the two-photon physics in heavy ion collisions see
51// G.Baur et al, Phys.Rep. 364 (2002), 359.
52//%
53// Author: Yuri.Kharlov@cern.ch
54// 15 April 2003
55
56#include <TParticle.h>
57#include <TParticlePDG.h>
58#include <TDatabasePDG.h>
37b09b91 59#include <TClonesArray.h>
60
0b5dd071 61#include "AliPythia.h"
62#include "AliRun.h"
63#include <AliGenTPHIC.h>
64#include <TPHICgen.h>
65#include <TPHICcommon.h>
66
67ClassImp(AliGenTPHIC)
68
69//------------------------------------------------------------
70
71AliGenTPHIC::AliGenTPHIC()
72{
73 // Constructor: create generator instance,
74 // create particle array
75 // Set TPHIC parameters to default values:
76 // eta_b production in Ca-Ca collisions at 7 A*TeV
77
78 SetMC(new TPHICgen());
79 fPythia = AliPythia::Instance();
80 fParticles = new TClonesArray("TParticle",100);
81
82 SetProcess ();
83 SetBeamEnergy();
84 SetBeamZ ();
85 SetBeamA ();
86 SetYggRange ();
87 SetMggRange ();
88 SetNgridY ();
89 SetNgridM ();
90 SetLumFunName();
91 SetLumFunFlag();
92 SetKfOnium ();
93}
94
0b5dd071 95//____________________________________________________________
96AliGenTPHIC::~AliGenTPHIC()
97{
98 // Destroys the object, deletes and disposes all TParticles currently on list.
99 if (fParticles) {
100 fParticles->Delete();
101 delete fParticles;
102 fParticles = 0;
103 }
104}
105
106//____________________________________________________________
107void AliGenTPHIC::Init()
108{
109 // Initialize the generator TPHIC
110
111 fTPHICgen->Initialize();
112 fEvent = 0;
113}
114
115//____________________________________________________________
116void AliGenTPHIC::Generate()
117{
118 // Generate one event of two-photon process.
119 // Gaussian smearing on the vertex is done if selected.
120 // All particles from the TPHIC/PYTHIA event listing
121 // are stored in TreeK, and only final-state particles are tracked.
122 // The event differectial cross section is assigned as a weight
123 // to each track of the event.
124
125 Float_t polar[3]= {0,0,0};
126 Float_t origin0[3],origin[3];
127 Float_t p[3], tof;
128 Double_t weight;
129
130 Int_t ks,kf,iparent,nt, trackIt;
131 Float_t random[6];
132 const Float_t kconv=0.001/2.999792458e8;
133
134 fTPHICgen->GenerateEvent();
135 weight = GetXSectionCurrent();
136 if (gAlice->GetEvNumber()>=fDebugEventFirst &&
137 gAlice->GetEvNumber()<=fDebugEventLast) fPythia->Pylist(1);
138 fPythia->ImportParticles(fParticles,"All");
139
140 if (fDebug == 1)
141 Info("Generate()","one event is produced");
142
143 Int_t j;
144 for (j=0;j<3;j++) origin[j]=fOrigin[j];
145 if(fVertexSmear==kPerEvent) {
146 Rndm(random,6);
147 for (j=0;j<3;j++) {
148 origin0[j]+=fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
149 TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
150 }
151 }
152
153 Int_t ip;
154 Int_t np = fParticles->GetEntriesFast();
155 TParticle *iparticle;
156// Int_t* pParent = new Int_t[np];
157 for (ip=0; ip<np; ip++) {
158 iparticle = (TParticle *) fParticles->At(ip);
159 ks = iparticle->GetStatusCode();
160// // No initial state partons
161// if (ks==21) continue;
162 p[0] = iparticle->Px();
163 p[1] = iparticle->Py();
164 p[2] = iparticle->Pz();
165 origin[0] = origin0[0]+iparticle->Vx()/10.;
166 origin[1] = origin0[1]+iparticle->Vy()/10.;
167 origin[2] = origin0[2]+iparticle->Vz()/10.;
168 kf = CheckPDGCode(iparticle->GetPdgCode());
169 iparent = -1;
170 tof = kconv*iparticle->T();
171 if (ks == 1) trackIt = 1;
172 else trackIt = 0;
642f15cf 173 PushTrack(fTrackIt*trackIt,iparent,kf,p,origin,polar,tof,kPPrimary,nt,weight,ks);
0b5dd071 174 KeepTrack(nt);
175
176 if (fDebug == 2)
177 printf("id=%+4d, parent=%3d, ks=%d, p = (%+11.4e,%+11.4e,%+11.4e) GeV\n",
178 kf,iparent,fTrackIt*trackIt,p[0],p[1],p[2]);
179 }
180 fEvent++;
181 fTPHICgen->SetNEVENT(fEvent);
182 if (fDebug == 1 && fEvent%100 == 0) {
183 Info("Generate","Event %d\n",fEvent);
184 fTPHICgen->Finish();
185 }
186}
187
188//____________________________________________________________
189void AliGenTPHIC::SetEventListRange(Int_t eventFirst, Int_t eventLast)
190{
191 // Set a range of event numbers, for which a table
192 // of generated particle will be printed
193 fDebugEventFirst = eventFirst;
194 fDebugEventLast = eventLast;
195 if (fDebugEventLast==-1) fDebugEventLast=fDebugEventFirst;
196}
197
198//____________________________________________________________
199void AliGenTPHIC::SetProcess (Int_t proc )
200{
201 // Set process number:
202 // proc=1 - gamma gamma -> X
203 // proc=2 - gamma gamma -> quarkonium
204 // proc=3 - gamma gamma -> fermion+ fermion-
205 // proc=4 - gamma gamma -> W+ W-
206 // proc=5 - not implemented
207 // proc=6 - gamma gamma -> V1 V2 (vector meson pair)
208
209 fTPHICgen->SetIPROC(proc);
210}
211//____________________________________________________________
212 void AliGenTPHIC::SetBeamEnergy (Float_t energy)
213{
214 // Set energy of the beam ion per nucleon in GeV
215 fTPHICgen->SetEBMN(energy);
216}
217//____________________________________________________________
218 void AliGenTPHIC::SetBeamZ (Int_t z )
219{
220 // Set beam ion charge
221 fTPHICgen->SetIZ(z);
222}
223//____________________________________________________________
224 void AliGenTPHIC::SetBeamA (Int_t a )
225{
226 // Set beam ion atomic number
227 fTPHICgen->SetIA(a);
228}
229//____________________________________________________________
230 void AliGenTPHIC::SetYggRange (Float_t ymin, Float_t ymax)
231{
232 // Set rapidity range of 2-photon system for the
233 // luminosity function calculation
234 fTPHICgen->SetYMIN(ymin);
235 fTPHICgen->SetYMAX(ymax);
236}
237//____________________________________________________________
238 void AliGenTPHIC::SetMggRange (Float_t mmin, Float_t mmax)
239{
240 // Set invariant mass range of 2-photon system for the
241 // luminosity function calculation
242 fTPHICgen->SetAMIN(mmin);
243 fTPHICgen->SetAMAX(mmax);
244}
245//____________________________________________________________
246 void AliGenTPHIC::SetNgridY (Int_t ny )
247{
248 // Set number of nodes on the grid along the rapidity axis
249 // to calculate the 2-photon luminosity function
250 fTPHICgen->SetNY(ny);
251}
252//____________________________________________________________
253 void AliGenTPHIC::SetNgridM (Int_t nm )
254{
255 // Set number of nodes on the grid along the mass axis
256 // to calculate the 2-photon luminosity function
257 fTPHICgen->SetNMAS(nm);
258}
259//____________________________________________________________
260 void AliGenTPHIC::SetLumFunName (TString name )
261{
262 // Set filename to store the 2-photon luminosity
263 // function calculated on the grid
264 fTPHICgen->SetLUMFIL(name);
265}
266//____________________________________________________________
267 void AliGenTPHIC::SetLumFunFlag (Int_t flag )
268{
269 // Set flag to calculate the 2-photon luminosity function:
270 // flag=-1 if a new lumimosity function to be calculated
271 // and stored in a file
272 // flag=+1 if a previously calculated function to be read
273 // from a file
274 fTPHICgen->SetILUMF(flag);
275}
276//____________________________________________________________
277 void AliGenTPHIC::SetKfFermion (Int_t kf )
278{
279 // Set a PDG flavour code of a fermion for the process 3,
280 // gamma gamma -> fermion+ fermion-
281 fTPHICgen->SetKFERM(kf);
282}
283//____________________________________________________________
284 void AliGenTPHIC::SetKfOnium (Int_t kf )
285{
286 // Set a PDG flavour code of a quarkonium for the process 2,
287 // gamma gamma -> quarkonium
288 fTPHICgen->SetKFONIUM(kf);
289}
290//____________________________________________________________
291 void AliGenTPHIC::SetMassOnium (Float_t mass )
292{
293 // Set a quarkonium mass [GeV] for the process 2 if it
294 // differes from the Pythia's particle table.
295 // For the well-known quarkonia no need to set this mass
296 // because it will be taken from the Pythia table
297 fTPHICgen->SetXMRES(mass);
298}
299//____________________________________________________________
300 void AliGenTPHIC::SetGGwidthOnium(Float_t width )
301{
302 // Set 2-photon partial width [GeV] of the quarkonium for
303 // the process 2 if it differes fromthe Pythia's particle table.
304 // For the well-known quarkonia no need to set this width
305 // because it will be taken from the Pythia table
306 fTPHICgen->SetXGGRES(width);
307}
308//____________________________________________________________
309 void AliGenTPHIC::SetKfVmesons (Int_t kf1, Int_t kf2)
310{
311 // Set PDG flavour codes of the two vector vesons
312 // for the process 2: gamma gamma -> V1 V2
313 // So far this processes is implemented for the following
314 // mesons and their combinations only:
315 // pho0 (113), omega (223), phi (333), J/psi (443)
316 fTPHICgen->SetKV1(kf1);
317 fTPHICgen->SetKV2(kf2);
318}
319//____________________________________________________________
320 Float_t AliGenTPHIC::GetGGmass ()
321{
322 // Get invariant mass of generated 2-photon system
323 return fTPHICgen->GetWSQ();
324}
325//____________________________________________________________
326 Float_t AliGenTPHIC::GetGGrapidity()
327{
328 // Get rapidity of generated 2-photon system
329 return fTPHICgen->GetYGG();
330}
331//____________________________________________________________
332 Float_t AliGenTPHIC::GetG1mass ()
333{
334 // Get a mass of the first virtual photon of
335 // the 2-photon process (-sqrt(q1^2)).
336 return fTPHICgen->GetXMG1();
337}
338//____________________________________________________________
339 Float_t AliGenTPHIC::GetG2mass ()
340{
341 // Get a mass of the second virtual photon of
342 // the 2-photon process (-sqrt(q2^2)).
343 return fTPHICgen->GetXMG2();
344}
345//____________________________________________________________
346 TClonesArray* AliGenTPHIC::GetParticleList ()
347{
348 // Get array of particles of the event listing
349 return fParticles;
350}
351//____________________________________________________________
352 TLorentzVector AliGenTPHIC::MomentumRecNucl1()
353{
354 // Get 4-momentum of the first recoil nucleus after
355 // the 2-photon process.
356 return TLorentzVector(fTPHICgen->GetPTAG1(1),
357 fTPHICgen->GetPTAG1(2),
358 fTPHICgen->GetPTAG1(3),
359 fTPHICgen->GetPTAG1(4));
360}
361//____________________________________________________________
362 TLorentzVector AliGenTPHIC::MomentumRecNucl2()
363{
364 // Get 4-momentum of the first recoil nucleus after
365 // the 2-photon process.
366 return TLorentzVector(fTPHICgen->GetPTAG2(1),
367 fTPHICgen->GetPTAG2(2),
368 fTPHICgen->GetPTAG2(3),
369 fTPHICgen->GetPTAG2(4));
370}
371//____________________________________________________________
372 Float_t AliGenTPHIC::GetXSectionCurrent()
373{
374 // Get the cross section of the produced event for the
375 // Monte Carlo integral cross section calculation
376 return fTPHICgen->GetXSCUR();
377}
378//____________________________________________________________
379 Float_t AliGenTPHIC::GetXSection ()
380{
381 // Get the integral cross section of the process
382 // calculated so far
383 return fTPHICgen->GetXSTOT();
384}
385//____________________________________________________________
386 Float_t AliGenTPHIC::GetXSectionError ()
387{
388 // Get the error of the integral cross section of the process
389 // calculated so far
390 return fTPHICgen->GetXSTOTE();
391}