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