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