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