]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenAfterBurnerFlow.cxx
Coverity 19921, 19913, 19912
[u/mrichter/AliRoot.git] / EVGEN / AliGenAfterBurnerFlow.cxx
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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 // 
20 // AliGenAfterBurnerFlow is a After Burner event generator applying flow.
21 // The generator changes Phi coordinate of the particle momentum.
22 // Flow (directed and elliptical) can be defined on particle type level
23 //
24 // Author:
25 // Sylwester Radomski, 2002
26 // Martin Poghosyan, 2008
27 // Constantin Loizides, 2010
28 //////////////////////////////////////////////////////////////////////////////
29
30 #include <Riostream.h>
31 #include <TParticle.h>
32 #include <TLorentzVector.h>
33 #include <TList.h>
34 #include <TRandom.h>
35 #include "AliStack.h"
36 #include "AliGenAfterBurnerFlow.h"
37 #include "AliGenCocktailAfterBurner.h"
38 #include "AliMC.h"
39 #include "AliRun.h"
40 #include "AliCollisionGeometry.h"
41 #include "AliGenCocktailEntry.h"
42
43
44 ClassImp(AliGenAfterBurnerFlow)
45
46
47 AliGenAfterBurnerFlow::AliGenAfterBurnerFlow():AliGenerator(),
48   fReactionPlane(0),
49   fHow(0),
50   fCounter(0),
51   fStack(0)
52 {
53   //
54   // Default Construction
55   InitPrimaries();
56   SetNpParams();
57 }
58
59 AliGenAfterBurnerFlow::AliGenAfterBurnerFlow(Float_t reactionPlane):AliGenerator(),
60   fReactionPlane(TMath::Pi()*reactionPlane/180.),
61   fHow(1),
62   fCounter(0),
63   fStack(0)
64 {
65   // reactionPlane - Reaction Plane Angle given in Deg [0-360]
66   // but stored and applied in radiants (standard for TParticle & AliCollisionGeometry) 
67
68   InitPrimaries();
69   SetNpParams();
70 }
71
72 ////////////////////////////////////////////////////////////////////////////////////////////////////
73
74 AliGenAfterBurnerFlow::~AliGenAfterBurnerFlow() 
75 {
76   // def. dest.
77 }
78
79 void AliGenAfterBurnerFlow::SetDirectedSimple(Int_t pdg, Float_t v1) 
80 {
81   //
82   // Set Directed Flow 
83   // The same directed flow is applied to all specified particles 
84   // independently on transverse momentum or rapidity
85   //
86   // PDG - particle type to apply directed flow
87   //       if (PDG == 0) use as default  
88   //
89
90   SetFlowParameters(pdg, 1, 0, v1, 0, 0, 0);
91 }
92
93 ////////////////////////////////////////////////////////////////////////////////////////////////////
94
95 void AliGenAfterBurnerFlow::SetDirectedParam(Int_t pdg, Float_t v11, Float_t v12, 
96                                              Float_t v13, Float_t v14) 
97 {
98   //
99   // Set Directed Flow 
100   // Directed flow is parameterised as follows
101   //
102   // V1(Pt,Y) = (V11 + V12*Pt) * sign(Y) * (V13 + V14 * abs(Y)^3)
103   //
104   // where sign = 1 for Y > 0 and -1 for Y < 0
105   // 
106   // Defaults values
107   // v12 = v14 = 0
108   // v13 = 1
109   //  
110   // PDG - particle type to apply directed flow
111   //       if (PDG == 0) use as default  
112   //
113   
114   SetFlowParameters(pdg, 1, 1, v11, v12, v13, v14);
115 }
116
117 ////////////////////////////////////////////////////////////////////////////////////////////////////
118
119 void AliGenAfterBurnerFlow::SetEllipticSimple(Int_t pdg, Float_t v2) 
120 {
121   //
122   // Set Elliptic Flow
123   // The same Elliptic flow is applied to all specified particles
124   // independently on transverse momentum or rapidity
125   //
126   // PDG - particle type to apply directed flow
127   //       if (PDG == 0) use as default  
128   //
129   // V2 - flow coefficient
130   //      
131   // NOTE: for starting playing with FLOW 
132   //       start with this function and values 0.05 - 0.1
133   //
134
135   SetFlowParameters(pdg, 2, 0, v2, 0, 0, 0);
136 }
137
138 ////////////////////////////////////////////////////////////////////////////////////////////////////
139
140 void AliGenAfterBurnerFlow::SetEllipticParam(Int_t pdg, 
141                                              Float_t v00, Float_t v10, Float_t v11,
142                                              Float_t v22) 
143 {
144   //
145   // Set Elliptic Flow
146   //
147   // Elliptic flow is parametrised to reproduce 
148   // V2 of Pions at RHIC energies and is given by:
149   // 
150   // V2 = (v00 + v10*pt + v11*pt^2) * exp (-v22 * y^2) and zero if V2<0.
151   //
152
153   SetFlowParameters(pdg, 2, 3, v00, v10, v11, v22);
154 }
155
156 void AliGenAfterBurnerFlow::SetEllipticParamPion(Int_t pdg, Float_t v21, 
157                                                  Float_t pTmax, Float_t v22) 
158 {
159   //
160   // Set Elliptic Flow
161   //
162   // Elliptic flow is parametrised to reproduce 
163   // V2 of Pions at RHIC energies and is given by:
164   // 
165   // V2 = v21 * (pT/pTMax ) * exp (-v22 * y^2)    where pT <= pTmax  
166   //      v21 * exp (-v22 * y^2)                  where pT > pTmax  
167   //
168   // v21   - value at saturation
169   // pTmax - saturation transverse momentum
170   // v22   - rapidity decreasing
171   //
172
173   SetFlowParameters(pdg, 2, 1, v21, pTmax, v22, 0);
174 }
175
176 ////////////////////////////////////////////////////////////////////////////////////////////////////
177
178 void AliGenAfterBurnerFlow::SetEllipticParamOld(Int_t pdg, Float_t v21, Float_t v22, Float_t v23) 
179 {
180   //
181   // Set Elliptic Flow
182   //
183   // Elliptic flow is parameterised using 
184   // old MevSim parameterisation
185   // 
186   // V2 = (V21 + V22 pT^2) * exp (-v22 * y^2)
187   //
188
189   SetFlowParameters(pdg, 2, 2, v21, v22, v23, 0);
190 }
191
192 ////////////////////////////////////////////////////////////////////////////////////////////////////
193
194 void AliGenAfterBurnerFlow::SetNpParams(Int_t order, Float_t p0, Float_t p1, Float_t p2, Float_t p3)
195 {
196   //
197   // Set npart parameterization.
198   //
199
200   fNpParams[0] = order;
201   fNpParams[1] = p0;
202   fNpParams[2] = p1;
203   fNpParams[3] = p2;
204   fNpParams[4] = p3;
205 }
206
207 ////////////////////////////////////////////////////////////////////////////////////////////////////
208
209 void AliGenAfterBurnerFlow::SetFlowParameters(Int_t pdg, Int_t order, Int_t type, 
210                                               Float_t v1, Float_t v2,Float_t v3,Float_t v4) 
211 {
212   // 
213   // private function
214   // 
215
216   if(TMath::Abs(pdg)>=fgkPDG){
217     Error("AliAfterBurnerFlow","Overflow");
218     return;
219   }
220   fIsPrim[TMath::Abs(pdg)]=kTRUE;
221   
222   Int_t index = 0;
223   Bool_t newEntry = kTRUE;
224
225   // Defaults
226   if (pdg == 0) {
227     index = fgkN - order;
228     newEntry = kFALSE;
229   }
230
231   // try to find existing entry
232   for (Int_t i=0; i<fCounter; i++) {
233     if (pdg == (Int_t)fParams[i][0] && 
234         order == (Int_t)fParams[i][1]) {
235       
236       index = i;
237       newEntry = kFALSE;
238     }
239   }
240   
241   // check fCounter
242
243   if (newEntry && (fCounter > fgkN-3)) {
244     Error("AliAfterBurnerFlow","Overflow");
245     return;
246   }
247   
248   if (newEntry) {
249     index = fCounter;
250     fCounter++;
251   }
252   
253   // Set new particle type
254
255   fParams[index][0] = pdg;
256   fParams[index][1] = order;
257   fParams[index][2] = type;
258   fParams[index][3] = v1;
259   fParams[index][4] = v2;
260   fParams[index][5] = v3;
261   fParams[index][6] = v4;  
262 }
263
264 ////////////////////////////////////////////////////////////////////////////////////////////////////
265
266 void AliGenAfterBurnerFlow::Init() 
267 {
268   // 
269   // Standard AliGenerator Initializer
270   //
271
272   if(fHow == 0)     { Info("AliGenAfterBurnerFlow", "Using the Hijing R.P. Angle event by event "); }
273   else if(fHow == 1){ Info("AliGenAfterBurnerFlow", "Using a fixed R.P. Angle for every event ") ; }
274   else { Info("AliGenAfterBurnerFlow", 
275               "Using a random R.P. Angle event by event ( ! not the same used by Hijing ! ) "); }
276 }
277
278 ////////////////////////////////////////////////////////////////////////////////////////////////////
279
280 Float_t AliGenAfterBurnerFlow::GetCoefficient(Int_t pdg, Int_t n, Float_t Pt, Float_t Y) const
281 {
282   //
283   // private function
284   // Return Flow Coefficient for a given particle type flow order
285   // and particle momentum (Pt, Y)
286   //
287
288   Int_t index = fgkN - n;  // default index (for all pdg)
289   Float_t v = 0;
290
291   // try to find specific parametrs
292
293   for (Int_t i=0; i<fCounter; i++) {
294     if ((Int_t)fParams[i][0] == pdg &&
295         (Int_t)fParams[i][1] == n) {
296       
297       index = i;
298       break;
299     }
300   } 
301   
302   // calculate v
303   
304   Int_t type = (Int_t)fParams[index][2];
305
306   if ((Int_t)fParams[index][1] == 1) { // Directed
307
308     if (type == 0 )
309       v = fParams[index][3];
310     else 
311       v = (fParams[index][3] + fParams[index][4] * Pt) * TMath::Sign((Float_t)1.,Y) *
312         (fParams[index][5] + fParams[index][6] * TMath::Abs(Y*Y*Y) );
313
314   } else {  // Elliptic
315
316     if (type == 0) v = fParams[index][3];
317
318     // Pion parameterisation 
319     if (type == 1) { 
320       if (Pt < fParams[index][4]) 
321         v = fParams[index][3] * (Pt / fParams[index][4]) ;
322       else 
323         v = fParams[index][3];
324       
325       v *= TMath::Exp( - fParams[index][5] * Y * Y);
326     }
327
328     // Old parameterisation
329     if (type == 2) 
330       v = (fParams[index][3] + fParams[index][4] * Pt * Pt) *
331         TMath::Exp( - fParams[index][5] * Y * Y);
332
333     // New v2 parameterisation
334     if (type == 3) {
335       v = (fParams[index][3] + fParams[index][4] *Pt + fParams[index][5] *Pt*Pt) *
336         TMath::Exp( - fParams[index][6] * Y*Y);
337       if (v<0)
338         v = 0;
339     }
340   }
341   
342   return v;
343 }
344
345 ////////////////////////////////////////////////////////////////////////////////////////////////////
346
347 Float_t AliGenAfterBurnerFlow::GetNpNorm(Int_t npart) const
348 {
349   //
350   // Calculate npart norm.
351   //
352
353   if (npart<0)
354     return 1;
355
356   Int_t order = (Int_t)fNpParams[0];
357   if (order<0)
358     return 1;
359
360   Float_t ret = 0;
361   Int_t npp = 1;
362   for (Int_t i=0; i<=order; i++) {
363     ret += npp*fNpParams[i+1];
364     npp *= npart;
365   } 
366   return ret;
367 }
368
369 ////////////////////////////////////////////////////////////////////////////////////////////////////
370
371 Bool_t AliGenAfterBurnerFlow::IsPrimary(Int_t pdg) const
372 {
373   if(pdg>=fgkPDG) return kFALSE;
374   return fIsPrim[pdg];
375 }
376
377 ////////////////////////////////////////////////////////////////////////////////////////////////////
378
379 Double_t CalcAngle(Double_t phi, Double_t phi0, Double_t phiRP, Double_t v2, Double_t v1=0.)
380 {
381     // Calculate relative angle
382   Double_t phi1 = phi-(phi+2*v1*TMath::Sin(phi-phiRP)+v2*TMath::Sin(2*(phi-phiRP))-phi0)/
383     (1.+2*v1*TMath::Cos(phi-phiRP)+ 2*v2*TMath::Cos(2*(phi-phiRP)));
384   if(TMath::Abs(phi/phi1-1.)<0.00001) return phi1;
385   return CalcAngle(phi1, phi0, phiRP, v2, v1);
386 }
387
388 ////////////////////////////////////////////////////////////////////////////////////////////////////
389
390 void AliGenAfterBurnerFlow::InitPrimaries()
391 {
392   // Init the primary particle list
393   for(Int_t i=0; i<fgkPDG; i++) fIsPrim[i]=kFALSE;
394
395   //mesons
396   fIsPrim[211]=kTRUE;
397   fIsPrim[311]=kTRUE;
398   fIsPrim[321]=kTRUE;
399   fIsPrim[411]=kTRUE;
400   fIsPrim[421]=kTRUE;
401   fIsPrim[431]=kTRUE;
402   fIsPrim[511]=kTRUE;
403   fIsPrim[521]=kTRUE;
404   fIsPrim[531]=kTRUE;
405   fIsPrim[541]=kTRUE;
406   fIsPrim[111]=kTRUE;
407   fIsPrim[221]=kTRUE;
408   fIsPrim[331]=kTRUE;
409   fIsPrim[441]=kTRUE;
410   fIsPrim[551]=kTRUE;
411   fIsPrim[130]=kTRUE;
412   fIsPrim[310]=kTRUE;
413   fIsPrim[213]=kTRUE;
414   fIsPrim[313]=kTRUE;
415   fIsPrim[323]=kTRUE;
416   fIsPrim[413]=kTRUE;
417   fIsPrim[423]=kTRUE;
418   fIsPrim[433]=kTRUE;
419   fIsPrim[513]=kTRUE;
420   fIsPrim[523]=kTRUE;
421   fIsPrim[533]=kTRUE;
422   fIsPrim[543]=kTRUE;
423   fIsPrim[113]=kTRUE;
424   fIsPrim[223]=kTRUE;
425   fIsPrim[333]=kTRUE;
426   fIsPrim[443]=kTRUE;
427   fIsPrim[553]=kTRUE;
428
429   //baryons
430   fIsPrim[2112]=kTRUE;
431   fIsPrim[2212]=kTRUE;
432   fIsPrim[3112]=kTRUE;
433   fIsPrim[3122]=kTRUE;
434   fIsPrim[3212]=kTRUE;
435   fIsPrim[3222]=kTRUE;
436   fIsPrim[3312]=kTRUE;
437   fIsPrim[3322]=kTRUE;
438   fIsPrim[4112]=kTRUE;
439   fIsPrim[4122]=kTRUE;
440   fIsPrim[4212]=kTRUE;
441   fIsPrim[4222]=kTRUE;
442   fIsPrim[4132]=kTRUE;
443   fIsPrim[4312]=kTRUE;
444   fIsPrim[4232]=kTRUE;
445   fIsPrim[4322]=kTRUE;
446   fIsPrim[4332]=kTRUE;
447   fIsPrim[5112]=kTRUE;
448   fIsPrim[5122]=kTRUE;
449   fIsPrim[5212]=kTRUE;
450   fIsPrim[5222]=kTRUE;
451   fIsPrim[1114]=kTRUE;
452   fIsPrim[2114]=kTRUE;
453   fIsPrim[2214]=kTRUE;
454   fIsPrim[2224]=kTRUE;
455   fIsPrim[3114]=kTRUE;
456   fIsPrim[3214]=kTRUE;
457   fIsPrim[3224]=kTRUE;
458   fIsPrim[3314]=kTRUE;
459   fIsPrim[3324]=kTRUE;
460   fIsPrim[3334]=kTRUE;
461   fIsPrim[4114]=kTRUE;
462   fIsPrim[4214]=kTRUE;
463   fIsPrim[4224]=kTRUE;
464   fIsPrim[4314]=kTRUE;
465   fIsPrim[4324]=kTRUE;
466   fIsPrim[4334]=kTRUE;
467   fIsPrim[5114]=kTRUE;
468   fIsPrim[5214]=kTRUE;
469   fIsPrim[5224]=kTRUE;
470 }
471
472 ////////////////////////////////////////////////////////////////////////////////////////////////////
473
474 void AliGenAfterBurnerFlow::Generate() 
475 {
476   // 
477   // AliGenerator generate function doing actual job.
478   // Algorythm:
479   //
480   // 1. loop over particles on the stack and choose primaries
481   // 2. calculate delta phi
482   // 3. change phi of primary particle and if it is non-stable 
483   //    then its daughters' phi and vertex also
484   // 
485   // For more details see :
486   // M.G. Poghosyan 
487   // PWG2 meeting on 06.05.2008 and 03.06.2008
488
489
490   if (0) 
491     for(Int_t ii=0; ii<fCounter;ii++)
492     {
493       printf("%d  %f  %f  %f  %f\n",ii,fParams[ii][0],fParams[ii][1],fParams[ii][2],fParams[ii][3]);
494     } 
495
496   AliGenCocktailAfterBurner *gen;
497   
498   TParticle *particle;
499   TParticle *particleM;
500   TLorentzVector momentum;
501   TLorentzVector vertex;
502
503   Int_t pdg;
504   Float_t phi;
505   Float_t pt, y;
506
507   // Get Stack of the first Generator
508   //  gen = (AliGenCocktailAfterBurner *)gAlice->Generator();
509   gen = (AliGenCocktailAfterBurner *)gAlice->GetMCApp()->Generator();
510
511
512   AliGenerator* genHijing = 0 ;
513   AliCollisionGeometry* geom = 0 ;
514   AliGenCocktailEntry* entry = 0 ;
515   TList* fEntries = 0 ;
516
517   TRandom* rand = new TRandom(0) ;
518   for(Int_t ns=0;ns<gen->GetNumberOfEvents();ns++) 
519   {
520     gen->SetActiveEventNumber(ns) ;
521    
522     fStack = gen->GetStack(ns);
523     fEntries = gen->Entries() ;
524
525     TIter next(fEntries) ;
526     Int_t npart = -1;
527
528     if(fHow == 0) // hijing R.P.
529     {
530       while((entry = (AliGenCocktailEntry*)next())) 
531       {
532         Info("Generate (e)","Using R.P. from HIJING ... ");       
533         genHijing = entry->Generator() ;        
534         if(genHijing->ProvidesCollisionGeometry()) 
535         { 
536           geom = gen->GetCollisionGeometry(ns) ;
537           fReactionPlane = geom->ReactionPlaneAngle() ; 
538           npart =  geom->ProjectileParticipants() + geom->TargetParticipants();
539           break;
540         }
541         else 
542         {
543           Error("Generate (e)", "NO CollisionGeometry !!!  -  using fixed R.P. angle = 0. ") ; 
544           fReactionPlane = 0. ; 
545         }
546       }
547     }
548     else if(fHow ==1 ) //  fixed R.P. 
549     {
550       Info("Generate (e)","Using fixed R.P. ...");       
551     }
552     else
553     {
554       Info("Generate (e)","Using random R.P.s ... ");       
555       fReactionPlane = 2 * TMath::Pi() * rand->Rndm() ;
556     }    
557     
558     cout << " * Reaction Plane Angle (event " << ns << ") = " << fReactionPlane << 
559       " rad. ( = " << (360*fReactionPlane/(2*TMath::Pi())) << " deg.) Npart = " << npart << "* " << endl ;
560
561     Int_t nParticles = fStack->GetNprimary();
562     for (Int_t i=0; i<nParticles; i++) 
563     {
564       particle = fStack->Particle(i);
565  
566       Int_t iM=particle->GetMother(0);
567       pdg = particle->GetPdgCode();
568
569       //exclude incoming protons in PYTHIA     
570       if(particle->GetPdgCode()==21) continue;
571
572       if(TMath::Abs(pdg)>=fgkPDG) continue; 
573       // is particle primary?
574       if(!fIsPrim[TMath::Abs(pdg)]) continue; 
575
576       if(iM>0)
577       {
578         particleM = fStack->Particle(iM);
579         Int_t pdgM = TMath::Abs(particleM->GetPdgCode());
580         // is mother primary?
581         if((TMath::Abs(pdgM)<fgkPDG)&&fIsPrim[TMath::Abs(pdgM)]) continue;
582       }
583    
584       particle->Momentum(momentum);
585       phi = particle->Phi();
586
587       // get Pt, y    
588       pt = momentum.Pt() ; 
589       y = 10000.;
590
591       if(TMath::Abs(momentum.Z()) != TMath::Abs(momentum.T())) 
592         y = momentum.Rapidity() ; 
593       
594       Double_t v1 = GetCoefficient(pdg, 1, pt, y);
595       Double_t v2 = GetCoefficient(pdg, 2, pt, y);
596       Double_t npartnorm = GetNpNorm(npart);
597       v2 *= npartnorm;
598
599       //printf("ntup %d %f  %f %f %f %f\n ",npart, v1, v2, pt, y, npartnorm);
600      
601       Double_t phi1 = CalcAngle(phi, phi, fReactionPlane,v2,v1);
602      
603       Rotate(i, phi1-phi);
604     }
605   }
606
607   Info("Generate","Flow After Burner: DONE");  
608 }
609
610 ////////////////////////////////////////////////////////////////////////////////////////////////////
611
612 void AliGenAfterBurnerFlow::Rotate(Int_t i, Double_t phi, Bool_t IsPrim)
613 {
614   // Rotation
615   TParticle*  particle = fStack->Particle(i);
616   
617   TLorentzVector momentum;
618   particle->Momentum(momentum);
619   momentum.RotateZ(phi);
620   particle->SetMomentum(momentum);
621  
622   if(!IsPrim)
623   {
624     TLorentzVector vertex;
625     particle->ProductionVertex(vertex);
626     vertex.RotateZ(phi);
627     particle->SetProductionVertex(vertex);
628   }
629   
630   if(particle->GetFirstDaughter()<0) return;     
631   for(Int_t iD=particle->GetFirstDaughter(); iD<=particle->GetLastDaughter(); iD++) Rotate(iD, phi, kFALSE);
632
633   return;
634 }
635
636