1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
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
25 // Sylwester Radomski, 2002
26 // Martin Poghosyan, 2008
27 // Constantin Loizides, 2010
28 //////////////////////////////////////////////////////////////////////////////
30 #include <Riostream.h>
31 #include <TParticle.h>
32 #include <TLorentzVector.h>
36 #include "AliGenAfterBurnerFlow.h"
37 #include "AliGenCocktailAfterBurner.h"
39 #include "AliCollisionGeometry.h"
40 #include "AliGenCocktailEntry.h"
43 ClassImp(AliGenAfterBurnerFlow)
46 AliGenAfterBurnerFlow::AliGenAfterBurnerFlow():AliGenerator(),
53 // Default Construction
58 AliGenAfterBurnerFlow::AliGenAfterBurnerFlow(Float_t reactionPlane):AliGenerator(),
59 fReactionPlane(TMath::Pi()*reactionPlane/180.),
64 // reactionPlane - Reaction Plane Angle given in Deg [0-360]
65 // but stored and applied in radiants (standard for TParticle & AliCollisionGeometry)
71 ////////////////////////////////////////////////////////////////////////////////////////////////////
73 AliGenAfterBurnerFlow::~AliGenAfterBurnerFlow()
78 void AliGenAfterBurnerFlow::SetDirectedSimple(Int_t pdg, Float_t v1)
82 // The same directed flow is applied to all specified particles
83 // independently on transverse momentum or rapidity
85 // PDG - particle type to apply directed flow
86 // if (PDG == 0) use as default
89 SetFlowParameters(pdg, 1, 0, v1, 0, 0, 0);
92 ////////////////////////////////////////////////////////////////////////////////////////////////////
94 void AliGenAfterBurnerFlow::SetDirectedParam(Int_t pdg, Float_t v11, Float_t v12,
95 Float_t v13, Float_t v14)
99 // Directed flow is parameterised as follows
101 // V1(Pt,Y) = (V11 + V12*Pt) * sign(Y) * (V13 + V14 * abs(Y)^3)
103 // where sign = 1 for Y > 0 and -1 for Y < 0
109 // PDG - particle type to apply directed flow
110 // if (PDG == 0) use as default
113 SetFlowParameters(pdg, 1, 1, v11, v12, v13, v14);
116 ////////////////////////////////////////////////////////////////////////////////////////////////////
118 void AliGenAfterBurnerFlow::SetEllipticSimple(Int_t pdg, Float_t v2)
122 // The same Elliptic flow is applied to all specified particles
123 // independently on transverse momentum or rapidity
125 // PDG - particle type to apply directed flow
126 // if (PDG == 0) use as default
128 // V2 - flow coefficient
130 // NOTE: for starting playing with FLOW
131 // start with this function and values 0.05 - 0.1
134 SetFlowParameters(pdg, 2, 0, v2, 0, 0, 0);
137 ////////////////////////////////////////////////////////////////////////////////////////////////////
139 void AliGenAfterBurnerFlow::SetEllipticParam(Int_t pdg,
140 Float_t v00, Float_t v10, Float_t v11,
146 // Elliptic flow is parametrised to reproduce
147 // V2 of Pions at RHIC energies and is given by:
149 // V2 = (v00 + v10*pt + v11*pt^2) * exp (-v22 * y^2) and zero if V2<0.
152 SetFlowParameters(pdg, 2, 3, v00, v10, v11, v22);
155 void AliGenAfterBurnerFlow::SetEllipticParamPion(Int_t pdg, Float_t v21,
156 Float_t pTmax, Float_t v22)
161 // Elliptic flow is parametrised to reproduce
162 // V2 of Pions at RHIC energies and is given by:
164 // V2 = v21 * (pT/pTMax ) * exp (-v22 * y^2) where pT <= pTmax
165 // v21 * exp (-v22 * y^2) where pT > pTmax
167 // v21 - value at saturation
168 // pTmax - saturation transverse momentum
169 // v22 - rapidity decreasing
172 SetFlowParameters(pdg, 2, 1, v21, pTmax, v22, 0);
175 ////////////////////////////////////////////////////////////////////////////////////////////////////
177 void AliGenAfterBurnerFlow::SetEllipticParamOld(Int_t pdg, Float_t v21, Float_t v22, Float_t v23)
182 // Elliptic flow is parameterised using
183 // old MevSim parameterisation
185 // V2 = (V21 + V22 pT^2) * exp (-v22 * y^2)
188 SetFlowParameters(pdg, 2, 2, v21, v22, v23, 0);
191 ////////////////////////////////////////////////////////////////////////////////////////////////////
193 void AliGenAfterBurnerFlow::SetNpParams(Int_t order, Float_t p0, Float_t p1, Float_t p2, Float_t p3)
196 // Set npart parameterization.
199 fNpParams[0] = order;
206 ////////////////////////////////////////////////////////////////////////////////////////////////////
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)
215 if(TMath::Abs(pdg)>fgkPDG){
216 Error("AliAfterBurnerFlow","Overflow");
219 fIsPrim[TMath::Abs(pdg)]=kTRUE;
222 Bool_t newEntry = kTRUE;
226 index = fgkN - order;
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]) {
242 if (newEntry && (fCounter > fgkN-3)) {
243 Error("AliAfterBurnerFlow","Overflow");
252 // Set new particle type
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;
263 ////////////////////////////////////////////////////////////////////////////////////////////////////
265 void AliGenAfterBurnerFlow::Init()
268 // Standard AliGenerator Initializer
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 ! ) "); }
277 ////////////////////////////////////////////////////////////////////////////////////////////////////
279 Float_t AliGenAfterBurnerFlow::GetCoefficient(Int_t pdg, Int_t n, Float_t Pt, Float_t Y) const
283 // Return Flow Coefficient for a given particle type flow order
284 // and particle momentum (Pt, Y)
287 Int_t index = fgkN - n; // default index (for all pdg)
290 // try to find specific parametrs
292 for (Int_t i=0; i<fCounter; i++) {
293 if ((Int_t)fParams[i][0] == pdg &&
294 (Int_t)fParams[i][1] == n) {
303 Int_t type = (Int_t)fParams[index][2];
305 if ((Int_t)fParams[index][1] == 1) { // Directed
308 v = fParams[index][3];
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) );
315 if (type == 0) v = fParams[index][3];
317 // Pion parameterisation
319 if (Pt < fParams[index][4])
320 v = fParams[index][3] * (Pt / fParams[index][4]) ;
322 v = fParams[index][3];
324 v *= TMath::Exp( - fParams[index][5] * Y * Y);
327 // Old parameterisation
329 v = (fParams[index][3] + fParams[index][4] * Pt * Pt) *
330 TMath::Exp( - fParams[index][5] * Y * Y);
332 // New v2 parameterisation
334 v = (fParams[index][3] + fParams[index][4] *Pt + fParams[index][5] *Pt*Pt) *
335 TMath::Exp( - fParams[index][6] * Y*Y);
344 ////////////////////////////////////////////////////////////////////////////////////////////////////
346 Float_t AliGenAfterBurnerFlow::GetNpNorm(Int_t npart)
349 // Calculate npart norm.
355 Int_t order = (Int_t)fNpParams[0];
361 for (Int_t i=0; i<=order; i++) {
362 ret += npp*fNpParams[i+1];
368 ////////////////////////////////////////////////////////////////////////////////////////////////////
370 Bool_t AliGenAfterBurnerFlow::IsPrimary(Int_t pdg)
372 if(pdg>fgkPDG) return kFALSE;
376 ////////////////////////////////////////////////////////////////////////////////////////////////////
378 Double_t CalcAngle(Double_t phi, Double_t phi0, Double_t phiRP, Double_t v2, Double_t v1=0.)
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);
386 ////////////////////////////////////////////////////////////////////////////////////////////////////
388 void AliGenAfterBurnerFlow::InitPrimaries()
390 for(Int_t i=0; i<fgkPDG; i++) fIsPrim[i]=kFALSE;
469 ////////////////////////////////////////////////////////////////////////////////////////////////////
471 void AliGenAfterBurnerFlow::Generate()
474 // AliGenerator generate function doing actual job.
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
482 // For more details see :
484 // PWG2 meeting on 06.05.2008 and 03.06.2008
488 for(Int_t ii=0; ii<fCounter;ii++)
490 printf("%d %f %f %f %f\n",ii,fParams[ii][0],fParams[ii][1],fParams[ii][2],fParams[ii][3]);
493 AliGenCocktailAfterBurner *gen;
496 TParticle *particleM;
497 TLorentzVector momentum;
498 TLorentzVector vertex;
504 // Get Stack of the first Generator
505 // gen = (AliGenCocktailAfterBurner *)gAlice->Generator();
506 gen = (AliGenCocktailAfterBurner *)gAlice->GetMCApp()->Generator();
509 AliGenerator* genHijing = 0 ;
510 AliCollisionGeometry* geom = 0 ;
511 AliGenCocktailEntry* entry = 0 ;
512 TList* fEntries = 0 ;
514 TRandom* rand = new TRandom(0) ;
515 for(Int_t ns=0;ns<gen->GetNumberOfEvents();ns++)
517 gen->SetActiveEventNumber(ns) ;
519 fStack = gen->GetStack(ns);
520 fEntries = gen->Entries() ;
522 TIter next(fEntries) ;
525 if(fHow == 0) // hijing R.P.
527 while((entry = (AliGenCocktailEntry*)next()))
529 Info("Generate (e)","Using R.P. from HIJING ... ");
530 genHijing = entry->Generator() ;
531 if(genHijing->ProvidesCollisionGeometry())
533 geom = gen->GetCollisionGeometry(ns) ;
534 fReactionPlane = geom->ReactionPlaneAngle() ;
535 npart = geom->ProjectileParticipants() + geom->TargetParticipants();
540 Error("Generate (e)", "NO CollisionGeometry !!! - using fixed R.P. angle = 0. ") ;
541 fReactionPlane = 0. ;
545 else if(fHow ==1 ) // fixed R.P.
547 Info("Generate (e)","Using fixed R.P. ...");
551 Info("Generate (e)","Using random R.P.s ... ");
552 fReactionPlane = 2 * TMath::Pi() * rand->Rndm() ;
555 cout << " * Reaction Plane Angle (event " << ns << ") = " << fReactionPlane <<
556 " rad. ( = " << (360*fReactionPlane/(2*TMath::Pi())) << " deg.) Npart = " << npart << "* " << endl ;
558 Int_t nParticles = fStack->GetNprimary();
559 for (Int_t i=0; i<nParticles; i++)
561 particle = fStack->Particle(i);
563 Int_t iM=particle->GetMother(0);
564 pdg = particle->GetPdgCode();
566 //exclude incoming protons in PYTHIA
567 if(particle->GetPdgCode()==21) continue;
569 if(TMath::Abs(pdg)>fgkPDG) continue;
570 // is particle primary?
571 if(!fIsPrim[TMath::Abs(pdg)]) continue;
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;
581 particle->Momentum(momentum);
582 phi = particle->Phi();
588 if(TMath::Abs(momentum.Z()) != TMath::Abs(momentum.T()))
589 y = momentum.Rapidity() ;
591 Double_t v1 = GetCoefficient(pdg, 1, pt, y);
592 Double_t v2 = GetCoefficient(pdg, 2, pt, y);
593 Double_t npartnorm = GetNpNorm(npart);
596 //printf("ntup %d %f %f %f %f %f\n ",npart, v1, v2, pt, y, npartnorm);
598 Double_t phi1 = CalcAngle(phi, phi, fReactionPlane,v2,v1);
604 Info("Generate","Flow After Burner: DONE");
607 ////////////////////////////////////////////////////////////////////////////////////////////////////
609 void AliGenAfterBurnerFlow::Rotate(Int_t i, Double_t phi, Bool_t IsPrim)
611 TParticle* particle = fStack->Particle(i);
613 TLorentzVector momentum;
614 particle->Momentum(momentum);
615 momentum.RotateZ(phi);
616 particle->SetMomentum(momentum);
620 TLorentzVector vertex;
621 particle->ProductionVertex(vertex);
623 particle->SetProductionVertex(vertex);
626 if(particle->GetFirstDaughter()<0) return;
627 for(Int_t iD=particle->GetFirstDaughter(); iD<=particle->GetLastDaughter(); iD++) Rotate(iD, phi, kFALSE);