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 **************************************************************************/
19 // AliGenGeVSim is a class implementing GeVSim event generator.
21 // GeVSim is a simple Monte-Carlo event generator for testing detector and
22 // algorythm performance especialy concerning flow and event-by-event studies
24 // In this event generator particles are generated from thermal distributions
25 // without any dynamics and addicional constrains. Distribution parameters like
26 // multiplicity, particle type yields, inverse slope parameters, flow coeficients
27 // and expansion velocities are expleicite defined by the user.
29 // GeVSim contains four thermal distributions the same as
30 // MevSim event generator developed for STAR experiment.
32 // In addition custom distributions can be used be the mean
33 // either two dimensional formula (TF2), a two dimensional histogram or
34 // two one dimensional histograms.
36 // Azimuthal distribution is deconvoluted from (Pt,Y) distribution
37 // and is described by two Fourier coefficients representing
38 // Directed and Elliptic flow.
40 ////////////////////////////////////////////////////////////////////////////////
42 // To apply flow to event ganerated by an arbitraly event generator
43 // refer to AliGenAfterBurnerFlow class.
45 ////////////////////////////////////////////////////////////////////////////////
47 // For examples, parameters and testing macros refer to:
48 // http:/home.cern.ch/radomski
50 // for more detailed description refer to ALICE NOTE
51 // "GeVSim Monte-Carlo Event Generator"
52 // S.Radosmki, P. Foka.
55 // Sylwester Radomski,
60 ////////////////////////////////////////////////////////////////////////////////
62 // Updated and revised: September 2002, S. Radomski, GSI
64 ////////////////////////////////////////////////////////////////////////////////
67 #include <Riostream.h>
73 #include <TObjArray.h>
75 #include <TParticle.h>
79 #include "AliGeVSimParticle.h"
80 #include "AliGenGeVSim.h"
81 #include "AliGenGeVSimEventHeader.h"
82 #include "AliGenerator.h"
86 ClassImp(AliGenGeVSim)
88 //////////////////////////////////////////////////////////////////////////////////
90 AliGenGeVSim::AliGenGeVSim() :
103 // Default constructor
106 for (Int_t i=0; i<4; i++)
108 for (Int_t i=0; i<2; i++)
112 //////////////////////////////////////////////////////////////////////////////////
114 AliGenGeVSim::AliGenGeVSim(Float_t psi, Bool_t isMultTotal)
118 fIsMultTotal(isMultTotal),
127 // Standard Constructor.
129 // models - thermal model to be used:
130 // 1 - deconvoluted pt and Y source
131 // 2,3 - thermalized sphericaly symetric sources
132 // 4 - thermalized source with expansion
133 // 5 - custom model defined in TF2 object named "gevsimPtY"
134 // 6 - custom model defined by two 1D histograms
135 // 7 - custom model defined by 2D histogram
137 // psi - reaction plane in degrees
138 // isMultTotal - multiplicity mode
139 // kTRUE - total multiplicity (default)
140 // kFALSE - dN/dY at midrapidity
143 // checking consistancy
145 if (psi < 0 || psi > 360 )
146 Error ("AliGenGeVSim", "Reaction plane angle ( %d )out of range [0..360]", psi);
148 fPsi = psi * TMath::Pi() / 180. ;
149 fIsMultTotal = isMultTotal;
153 fPartTypes = new TObjArray();
157 //////////////////////////////////////////////////////////////////////////////////
159 AliGenGeVSim::~AliGenGeVSim() {
161 // Default Destructor
163 // Removes TObjArray keeping list of registered particle types
166 if (fPartTypes != NULL) delete fPartTypes;
170 //////////////////////////////////////////////////////////////////////////////////
172 Bool_t AliGenGeVSim::CheckPtYPhi(Float_t pt, Float_t y, Float_t phi) const {
174 // private function used by Generate()
176 // Check bounds of Pt, Rapidity and Azimuthal Angle of a track
177 // Used only when generating particles from a histogram
180 if ( TestBit(kPtRange) && ( pt < fPtMin || pt > fPtMax )) return kFALSE;
181 if ( TestBit(kPhiRange) && ( phi < fPhiMin || phi > fPhiMax )) return kFALSE;
182 if ( TestBit(kYRange) && ( y < fYMin || y > fYMax )) return kFALSE;
187 //////////////////////////////////////////////////////////////////////////////////
189 Bool_t AliGenGeVSim::CheckAcceptance(Float_t p[3]) {
191 // private function used by Generate()
193 // Check bounds of a total momentum and theta of a track
196 if ( TestBit(kThetaRange) ) {
198 Double_t theta = TMath::ATan2( TMath::Sqrt(p[0]*p[0]+p[1]*p[1]), p[2]);
199 if ( theta < fThetaMin || theta > fThetaMax ) return kFALSE;
203 if ( TestBit(kMomentumRange) ) {
205 Double_t momentum = TMath::Sqrt(p[0]*p[0] + p[1]*p[1] + p[2]*p[2]);
206 if ( momentum < fPMin || momentum > fPMax) return kFALSE;
212 //////////////////////////////////////////////////////////////////////////////////
214 // Deconvoluted Pt Y formula
216 static Double_t aPtForm(Double_t * x, Double_t * par) {
217 // ptForm: pt -> x[0] , mass -> [0] , temperature -> [1]
218 // Description as string: " x * exp( -sqrt([0]*[0] + x*x) / [1] )"
220 return x[0] * TMath::Exp( -sqrt(par[0]*par[0] + x[0]*x[0]) / par[1]);
223 static Double_t aYForm(Double_t * x, Double_t * par) {
224 // y Form: y -> x[0] , sigmaY -> [0]
225 // Description as string: " exp ( - x*x / (2 * [0]*[0] ) )"
227 return TMath::Exp ( - x[0]*x[0] / (2 * par[0]*par[0] ) );
231 // Description as strings:
233 // const char *kFormE = " ( sqrt([0]*[0] + x*x) * cosh(y) ) ";
234 // const char *kFormG = " ( 1 / sqrt( 1 - [2]*[2] ) ) ";
235 // const char *kFormYp = "( [2]*sqrt(([0]*[0]+x*x)*cosh(y)*cosh(y)-[0]*[0])/([1]*sqrt(1-[2]*[2]))) ";
237 // const char* kFormula[3] = {
238 // " x * %s * exp( -%s / [1]) ",
239 // " (x * %s) / ( exp( %s / [1]) - 1 ) ",
240 // " x*%s*exp(-%s*%s/[1])*((sinh(%s)/%s)+([1]/(%s*%s))*(sinh(%s)/%s-cosh(%s)))"
242 // printf(kFormula[0], kFormE, kFormE);
243 // printf(kFormula[1], kFormE, kFormE);
244 // printf(kFormula[2], kFormE, kFormG, kFormE, kFormYp, kFormYp, kFormG, kFormE, kFormYp, kFormYp, kFormYp);
247 static Double_t aPtYFormula0(Double_t *x, Double_t * par) {
249 // mass -> [0] , temperature -> [1] , expansion velocity -> [2]
251 Double_t aFormE = TMath::Sqrt(par[0]*par[0] + x[0]*x[0]) * TMath::CosH(x[1]);
252 return x[0] * aFormE * TMath::Exp(-aFormE/par[1]);
255 static Double_t aPtYFormula1(Double_t *x, Double_t * par) {
257 // mass -> [0] , temperature -> [1] , expansion velocity -> [2]
259 Double_t aFormE = TMath::Sqrt(par[0]*par[0] + x[0]*x[0]) * TMath::CosH(x[1]);
260 return x[0] * aFormE / ( TMath::Exp( aFormE / par[1]) - 1 );
263 static Double_t aPtYFormula2(Double_t *x, Double_t * par) {
265 // mass -> [0] , temperature -> [1] , expansion velocity -> [2]
267 Double_t aFormE = TMath::Sqrt(par[0]*par[0] + x[0]*x[0]) * TMath::CosH(x[1]);
268 Double_t aFormG = 1 / TMath::Sqrt((1.-par[2])*(1.+par[2]));
269 Double_t aFormYp = par[2]*TMath::Sqrt( (par[0]*par[0] + x[0]*x[0])
270 * (TMath::CosH(x[1])-par[0])*(TMath::CosH(x[1])+par[0]))
271 /( par[1]*TMath::Sqrt((1.-par[2])*(1.+par[2])));
273 return x[0] * aFormE * TMath::Exp( - aFormG * aFormE / par[1])
274 *( TMath::SinH(aFormYp)/aFormYp
275 + par[1]/(aFormG*aFormE)
276 * ( TMath::SinH(aFormYp)/aFormYp-TMath::CosH(aFormYp) ) );
281 static Double_t aPhiForm(Double_t * x, Double_t * par) {
283 // Psi -> [0] , Direct Flow -> [1] , Elliptical Flow -> [2]
284 // Description as string: " 1 + 2*[1]*cos(x-[0]) + 2*[2]*cos(2*(x-[0])) "
286 return 1 + 2*par[1]*TMath::Cos(x[0]-par[0])
287 + 2*par[2]*TMath::Cos(2*(x[0]-par[0]));
290 void AliGenGeVSim::InitFormula() {
294 // Initalizes formulas used in GeVSim.
296 // Deconvoluted Pt Y formula
298 fPtFormula = new TF1("gevsimPt", &aPtForm, 0, 3, 2);
299 fYFormula = new TF1("gevsimRapidity", &aYForm, -3, 3,1);
301 fPtFormula->SetParNames("mass", "temperature");
302 fPtFormula->SetParameters(1., 1.);
304 fYFormula->SetParName(0, "sigmaY");
305 fYFormula->SetParameter(0, 1.);
308 fPtFormula->SetNpx(100);
309 fYFormula->SetNpx(100);
314 fPtYFormula[0] = new TF2("gevsimPtY_2", &aPtYFormula0, 0, 3, -2, 2, 2);
316 fPtYFormula[1] = new TF2("gevsimPtY_3", &aPtYFormula1, 0, 3, -2, 2, 2);
318 fPtYFormula[2] = new TF2("gevsimPtY_4", &aPtYFormula2, 0, 3, -2, 2, 3);
323 // setting names & initialisation
325 const char* kParamNames[3] = {"mass", "temperature", "expVel"};
328 for (i=0; i<3; i++) {
330 fPtYFormula[i]->SetNpx(100); // step 30 MeV
331 fPtYFormula[i]->SetNpy(100); //
333 for (j=0; j<3; j++) {
335 if ( i != 2 && j == 2 ) continue; // ExpVel
336 fPtYFormula[i]->SetParName(j, kParamNames[j]);
337 fPtYFormula[i]->SetParameter(j, 0.5);
343 fPhiFormula = new TF1("gevsimPhi", &aPhiForm, 0, 2*TMath::Pi(), 3);
345 fPhiFormula->SetParNames("psi", "directed", "elliptic");
346 fPhiFormula->SetParameters(0., 0., 0.);
348 fPhiFormula->SetNpx(180);
352 //////////////////////////////////////////////////////////////////////////////////
354 void AliGenGeVSim::AddParticleType(AliGeVSimParticle *part) {
356 // Adds new type of particles.
358 // Parameters are defeined in AliGeVSimParticle object
359 // This method has to be called for every particle type
362 if (fPartTypes == NULL)
363 fPartTypes = new TObjArray();
365 fPartTypes->Add(part);
368 //////////////////////////////////////////////////////////////////////////////////
370 void AliGenGeVSim::SetMultTotal(Bool_t isTotal) {
375 fIsMultTotal = isTotal;
378 //////////////////////////////////////////////////////////////////////////////////
380 Float_t AliGenGeVSim::FindScaler(Int_t paramId, Int_t pdg) {
383 // Finds Scallar for a given parameter.
384 // Function used in event-by-event mode.
386 // There are two types of scallars: deterministic and random
387 // and they can work on either global or particle type level.
388 // For every variable there are four scallars defined.
390 // Scallars are named as folowa
391 // deterministic global level : "gevsimParam" (eg. "gevsimTemp")
392 // deterinistig type level : "gevsimPdgParam" (eg. "gevsim211Mult")
393 // random global level : "gevsimParamRndm" (eg. "gevsimMultRndm")
394 // random type level : "gevsimPdgParamRndm" (eg. "gevsim-211V2Rndm");
396 // Pdg - code of a particle type in PDG standard (see: http://pdg.lbl.gov)
397 // Param - parameter name. Allowed parameters:
399 // "Temp" - inverse slope parameter
400 // "SigmaY" - rapidity width - for model (1) only
401 // "ExpVel" - expansion velocity - for model (4) only
402 // "V1" - directed flow
403 // "V2" - elliptic flow
404 // "Mult" - multiplicity
408 static const char* params[] = {"Temp", "SigmaY", "ExpVel", "V1", "V2", "Mult"};
409 static const char* ending[] = {"", "Rndm"};
411 static const char* patt1 = "gevsim%s%s";
412 static const char* patt2 = "gevsim%d%s%s";
419 // Scaler evoluation: i - global/local, j - determ/random
423 for (i=0; i<2; i++) {
424 for (j=0; j<2; j++) {
428 if (i == 0) sprintf(buffer, patt1, params[paramId], ending[j]);
429 else sprintf(buffer, patt2, pdg, params[paramId], ending[j]);
431 form = (TF1 *)gROOT->GetFunction(buffer);
434 if (j == 0) scaler *= form->Eval(gAlice->GetEvNumber());
436 form->SetParameter(0, gAlice->GetEvNumber());
437 scaler *= form->GetRandom();
446 //////////////////////////////////////////////////////////////////////////////////
448 void AliGenGeVSim::DetermineReactionPlane() {
450 // private function used by Generate()
452 // Retermines Reaction Plane angle and set this value
453 // as a parameter [0] in fPhiFormula
455 // Note: if "gevsimPsiRndm" function is found it override both
456 // "gevsimPhi" function and initial fPsi value
462 form = (TF1 *)gROOT->GetFunction("gevsimPsi");
463 if (form) fPsi = form->Eval(gAlice->GetEvNumber()) * TMath::Pi() / 180;
466 form = (TF1 *)gROOT->GetFunction("gevsimPsiRndm");
467 if (form) fPsi = form->GetRandom() * TMath::Pi() / 180;
470 cout << "Psi = " << fPsi << "\t" << (Int_t)(fPsi*180./TMath::Pi()) << endl;
472 fPhiFormula->SetParameter(0, fPsi);
475 //////////////////////////////////////////////////////////////////////////////////
477 Float_t AliGenGeVSim::GetdNdYToTotal() {
479 // Private, helper function used by Generate()
481 // Returns total multiplicity to dN/dY ration using current distribution.
482 // The function have to be called after setting distribution and its
483 // parameters (like temperature).
487 const Double_t kMaxPt = 3.0, kMaxY = 2.;
491 integ = fYFormula->Integral(-kMaxY, kMaxY);
492 mag = fYFormula->Eval(0);
496 // 2D formula standard or custom
498 if (fModel > 1 && fModel < 6) {
500 integ = ((TF2*)fCurrentForm)->Integral(0,kMaxPt, -kMaxY, kMaxY);
501 mag = ((TF2*)fCurrentForm)->Integral(0, kMaxPt, -0.1, 0.1) / 0.2;
509 integ = fHist[1]->Integral();
510 mag = fHist[0]->GetBinContent(fHist[0]->FindBin(0.));
511 mag /= fHist[0]->GetBinWidth(fHist[0]->FindBin(0.));
520 Int_t yBins = fPtYHist->GetNbinsY();
521 Int_t ptBins = fPtYHist->GetNbinsX();
523 integ = fPtYHist->Integral(0, ptBins, 0, yBins);
524 mag = fPtYHist->Integral(0, ptBins, (yBins/2)-1, (yBins/2)+1 ) / 2;
531 //////////////////////////////////////////////////////////////////////////////////
533 void AliGenGeVSim::SetFormula(Int_t pdg) {
535 // Private function used by Generate()
537 // Configure a formula for a given particle type and model Id (in fModel).
538 // If custom formula or histogram was selected the function tries
541 // The function implements naming conventions for custom distributions names
545 const char* msg[4] = {
546 "Custom Formula for Pt Y distribution not found [pdg = %d]",
547 "Histogram for Pt distribution not found [pdg = %d]",
548 "Histogram for Y distribution not found [pdg = %d]",
549 "HIstogram for Pt Y dostribution not found [pdg = %d]"
552 const char* pattern[8] = {
553 "gevsimDistPtY", "gevsimDist%dPtY",
554 "gevsimHistPt", "gevsimHist%dPt",
555 "gevsimHistY", "gevsimHist%dY",
556 "gevsimHistPtY", "gevsimHist%dPtY"
559 const char *where = "SetFormula";
562 if (fModel < 1 || fModel > 7)
563 Error("SetFormula", "Model Id (%d) out of range [1-7]", fModel);
568 if (fModel == 1) fCurrentForm = fPtFormula;
569 if (fModel > 1 && fModel < 5) fCurrentForm = fPtYFormula[fModel-2];
572 // custom model defined by a formula
577 fCurrentForm = (TF2*)gROOT->GetFunction(pattern[0]);
581 sprintf(buff, pattern[1], pdg);
582 fCurrentForm = (TF2*)gROOT->GetFunction(buff);
584 if (!fCurrentForm) Error(where, msg[0], pdg);
592 for (Int_t i=0; i<2; i++) {
595 fHist[i] = (TH1D*)gROOT->FindObject(pattern[2+2*i]);
599 sprintf(buff, pattern[3+2*i], pdg);
600 fHist[i] = (TH1D*)gROOT->FindObject(buff);
602 if (!fHist[i]) Error(where, msg[1+i], pdg);
612 fPtYHist = (TH2D*)gROOT->FindObject(pattern[6]);
616 sprintf(buff, pattern[7], pdg);
617 fPtYHist = (TH2D*)gROOT->FindObject(buff);
620 if (!fPtYHist) Error(where, msg[3], pdg);
625 //////////////////////////////////////////////////////////////////////////////////
627 void AliGenGeVSim:: AdjustFormula() {
630 // Adjust fomula bounds according to acceptance cuts.
632 // Since GeVSim is producing "thermal" particles Pt
633 // is cut at 3 GeV even when acceptance extends to grater momenta.
636 // If custom formula was provided function preserves
640 const Double_t kMaxPt = 3.0;
641 const Double_t kMaxY = 2.0;
642 Double_t minPt, maxPt, minY, maxY;
645 if (fModel > 4) return;
648 if (TestBit(kPtRange) && fPtMax < kMaxPt ) maxPt = fPtMax;
652 if (TestBit(kPtRange)) minPt = fPtMin;
655 if (TestBit(kPtRange) && fPtMin > kMaxPt )
656 Warning("Acceptance", "Minimum Pt (%3.2f GeV) greater that 3.0 GeV ", fPtMin);
659 if (TestBit(kMomentumRange) && fPtMax < maxPt) maxPt = fPtMax;
661 // max and min rapidity
662 if (TestBit(kYRange)) {
673 fPtFormula->SetRange(fPtMin, maxPt);
674 fYFormula->SetRange(fYMin, fYMax);
678 ((TF2*)fCurrentForm)->SetRange(minPt, minY, maxPt, maxY);
682 if (TestBit(kPhiRange))
683 fPhiFormula->SetRange(fPhiMin, fPhiMax);
687 //////////////////////////////////////////////////////////////////////////////////
689 void AliGenGeVSim::GetRandomPtY(Double_t &pt, Double_t &y) {
691 // Private function used by Generate()
693 // Returns random values of Pt and Y corresponding to selected
698 pt = fPtFormula->GetRandom();
699 y = fYFormula->GetRandom();
703 if (fModel > 1 && fModel < 6) {
704 ((TF2*)fCurrentForm)->GetRandom2(pt, y);
709 pt = fHist[0]->GetRandom();
710 y = fHist[1]->GetRandom();
714 fPtYHist->GetRandom2(pt, y);
719 //////////////////////////////////////////////////////////////////////////////////
721 void AliGenGeVSim::Init() {
723 // Standard AliGenerator initializer.
728 //////////////////////////////////////////////////////////////////////////////////
730 void AliGenGeVSim::Generate() {
732 // Standard AliGenerator function
733 // This function do actual job and puts particles on stack.
736 PDG_t pdg; // particle type
737 Float_t mass; // particle mass
738 Float_t orgin[3] = {0,0,0}; // particle orgin [cm]
739 Float_t polar[3] = {0,0,0}; // polarisation
740 Float_t time = 0; // time of creation
742 Float_t multiplicity = 0;
743 Bool_t isMultTotal = kTRUE;
746 Float_t directedScaller = 1., ellipticScaller = 1.;
748 TLorentzVector *v = new TLorentzVector(0,0,0,0);
750 const Int_t kParent = -1;
755 orgin[0] = fVertex[0];
756 orgin[1] = fVertex[1];
757 orgin[2] = fVertex[2];
760 // Particle params database
762 TDatabasePDG *db = TDatabasePDG::Instance();
764 AliGeVSimParticle *partType;
766 Int_t nType, nParticle, nParam;
767 const Int_t kNParams = 6;
769 // reaction plane determination and model
770 DetermineReactionPlane();
772 // loop over particle types
774 for (nType = 0; nType < fPartTypes->GetEntries(); nType++) {
776 partType = (AliGeVSimParticle *)fPartTypes->At(nType);
778 pdg = (PDG_t)partType->GetPdgCode();
779 type = db->GetParticle(pdg);
782 fModel = partType->GetModel();
784 fCurrentForm->SetParameter("mass", mass);
787 // Evaluation of parameters - loop over parameters
789 for (nParam = 0; nParam < kNParams; nParam++) {
791 paramScaler = FindScaler(nParam, pdg);
794 fCurrentForm->SetParameter("temperature", paramScaler * partType->GetTemperature());
796 if (nParam == 1 && fModel == 1)
797 fYFormula->SetParameter("sigmaY", paramScaler * partType->GetSigmaY());
799 if (nParam == 2 && fModel == 4) {
801 Double_t totalExpVal = paramScaler * partType->GetExpansionVelocity();
803 if (totalExpVal == 0.0) totalExpVal = 0.0001;
804 if (totalExpVal == 1.0) totalExpVal = 9.9999;
806 fCurrentForm->SetParameter("expVel", totalExpVal);
811 if (nParam == 3) directedScaller = paramScaler;
812 if (nParam == 4) ellipticScaller = paramScaler;
818 if (partType->IsMultForced()) isMultTotal = partType->IsMultTotal();
819 else isMultTotal = fIsMultTotal;
821 multiplicity = paramScaler * partType->GetMultiplicity();
822 multiplicity *= (isMultTotal)? 1 : GetdNdYToTotal();
826 // Flow defined on the particle type level (not parameterised)
827 if (partType->IsFlowSimple()) {
828 fPhiFormula->SetParameter(1, partType->GetDirectedFlow(0,0) * directedScaller);
829 fPhiFormula->SetParameter(2, partType->GetEllipticFlow(0,0) * ellipticScaller);
835 Info("Generate","PDG = %d \t Mult = %d", pdg, (Int_t)multiplicity);
837 // loop over particles
840 while (nParticle < multiplicity) {
842 Double_t pt, y, phi; // momentum in [pt,y,phi]
843 Float_t p[3] = {0,0,0}; // particle momentum
847 // phi distribution configuration when differential flow defined
848 // to be optimised in future release
850 if (!partType->IsFlowSimple()) {
851 fPhiFormula->SetParameter(1, partType->GetDirectedFlow(pt,y) * directedScaller);
852 fPhiFormula->SetParameter(2, partType->GetEllipticFlow(pt,y) * ellipticScaller);
855 phi = fPhiFormula->GetRandom();
857 if (!isMultTotal) nParticle++;
858 if (fModel > 4 && !CheckPtYPhi(pt,y,phi) ) continue;
860 // coordinate transformation
861 v->SetPtEtaPhiM(pt, y, phi, mass);
867 // momentum range test
868 if ( !CheckAcceptance(p) ) continue;
870 // putting particle on the stack
872 PushTrack(fTrackIt, kParent, pdg, p, orgin, polar, time, kPPrimary, id, fTrackIt);
873 if (isMultTotal) nParticle++;
877 // prepare and store header
879 AliGenGeVSimEventHeader *header = new AliGenGeVSimEventHeader("GeVSim header");
880 TArrayF eventVertex(3,orgin);
882 header->SetPrimaryVertex(eventVertex);
883 header->SetEventPlane(fPsi);
884 header->SetEllipticFlow(fPhiFormula->GetParameter(2));
886 gAlice->SetGenEventHeader(header);
891 //////////////////////////////////////////////////////////////////////////////////