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() : AliGenerator(-1) {
92 // Default constructor
99 for (Int_t i=0; i<4; i++)
104 //////////////////////////////////////////////////////////////////////////////////
106 AliGenGeVSim::AliGenGeVSim(Float_t psi, Bool_t isMultTotal) : AliGenerator(-1) {
108 // Standard Constructor.
110 // models - thermal model to be used:
111 // 1 - deconvoluted pt and Y source
112 // 2,3 - thermalized sphericaly symetric sources
113 // 4 - thermalized source with expansion
114 // 5 - custom model defined in TF2 object named "gevsimPtY"
115 // 6 - custom model defined by two 1D histograms
116 // 7 - custom model defined by 2D histogram
118 // psi - reaction plane in degrees
119 // isMultTotal - multiplicity mode
120 // kTRUE - total multiplicity (default)
121 // kFALSE - dN/dY at midrapidity
124 // checking consistancy
126 if (psi < 0 || psi > 360 )
127 Error ("AliGenGeVSim", "Reaction plane angle ( %d )out of range [0..360]", psi);
129 fPsi = psi * TMath::Pi() / 180. ;
130 fIsMultTotal = isMultTotal;
134 fPartTypes = new TObjArray();
138 //////////////////////////////////////////////////////////////////////////////////
140 AliGenGeVSim::~AliGenGeVSim() {
142 // Default Destructor
144 // Removes TObjArray keeping list of registered particle types
147 if (fPartTypes != NULL) delete fPartTypes;
151 //////////////////////////////////////////////////////////////////////////////////
153 Bool_t AliGenGeVSim::CheckPtYPhi(Float_t pt, Float_t y, Float_t phi) const {
155 // private function used by Generate()
157 // Check bounds of Pt, Rapidity and Azimuthal Angle of a track
158 // Used only when generating particles from a histogram
161 if ( TestBit(kPtRange) && ( pt < fPtMin || pt > fPtMax )) return kFALSE;
162 if ( TestBit(kPhiRange) && ( phi < fPhiMin || phi > fPhiMax )) return kFALSE;
163 if ( TestBit(kYRange) && ( y < fYMin || y > fYMax )) return kFALSE;
168 //////////////////////////////////////////////////////////////////////////////////
170 Bool_t AliGenGeVSim::CheckAcceptance(Float_t p[3]) {
172 // private function used by Generate()
174 // Check bounds of a total momentum and theta of a track
177 if ( TestBit(kThetaRange) ) {
179 Double_t theta = TMath::ATan2( TMath::Sqrt(p[0]*p[0]+p[1]*p[1]), p[2]);
180 if ( theta < fThetaMin || theta > fThetaMax ) return kFALSE;
184 if ( TestBit(kMomentumRange) ) {
186 Double_t momentum = TMath::Sqrt(p[0]*p[0] + p[1]*p[1] + p[2]*p[2]);
187 if ( momentum < fPMin || momentum > fPMax) return kFALSE;
193 //////////////////////////////////////////////////////////////////////////////////
195 void AliGenGeVSim::InitFormula() {
199 // Initalizes formulas used in GeVSim.
200 // Manages strings and creates TFormula objects from strings
203 // Deconvoluted Pt Y formula
205 // ptForm: pt -> x , mass -> [0] , temperature -> [1]
206 // y Form: y -> x , sigmaY -> [0]
208 const char* kPtForm = " x * exp( -sqrt([0]*[0] + x*x) / [1] )";
209 const char* kYForm = " exp ( - x*x / (2 * [0]*[0] ) )";
211 fPtFormula = new TF1("gevsimPt", kPtForm, 0, 3);
212 fYFormula = new TF1("gevsimRapidity", kYForm, -3, 3);
214 fPtFormula->SetParNames("mass", "temperature");
215 fPtFormula->SetParameters(1., 1.);
217 fYFormula->SetParName(0, "sigmaY");
218 fYFormula->SetParameter(0, 1.);
221 fPtFormula->SetNpx(100);
222 fYFormula->SetNpx(100);
228 // mass -> [0] , temperature -> [1] , expansion velocity -> [2]
231 const char *kFormE = " ( sqrt([0]*[0] + x*x) * cosh(y) ) ";
232 const char *kFormG = " ( 1 / sqrt( 1 - [2]*[2] ) ) ";
233 const char *kFormYp = "( [2]*sqrt(([0]*[0]+x*x)*cosh(y)*cosh(y)-[0]*[0])/([1]*sqrt(1-[2]*[2]))) ";
235 const char* kFormula[3] = {
236 " x * %s * exp( -%s / [1]) ",
237 " (x * %s) / ( exp( %s / [1]) - 1 ) ",
238 " x*%s*exp(-%s*%s/[1])*((sinh(%s)/%s)+([1]/(%s*%s))*(sinh(%s)/%s-cosh(%s)))"
241 const char* kParamNames[3] = {"mass", "temperature", "expVel"};
245 sprintf(buffer, kFormula[0], kFormE, kFormE);
246 fPtYFormula[0] = new TF2("gevsimPtY_2", buffer, 0, 3, -2, 2);
248 sprintf(buffer, kFormula[1], kFormE, kFormE);
249 fPtYFormula[1] = new TF2("gevsimPtY_3", buffer, 0, 3, -2, 2);
251 sprintf(buffer, kFormula[2], kFormE, kFormG, kFormE, kFormYp, kFormYp, kFormG, kFormE, kFormYp, kFormYp, kFormYp);
252 fPtYFormula[2] = new TF2("gevsimPtY_4", buffer, 0, 3, -2, 2);
257 // setting names & initialisation
260 for (i=0; i<3; i++) {
262 fPtYFormula[i]->SetNpx(100); // step 30 MeV
263 fPtYFormula[i]->SetNpy(100); //
265 for (j=0; j<3; j++) {
267 if ( i != 2 && j == 2 ) continue; // ExpVel
268 fPtYFormula[i]->SetParName(j, kParamNames[j]);
269 fPtYFormula[i]->SetParameter(j, 0.5);
276 // Psi -> [0] , Direct Flow -> [1] , Ellipticla Flow -> [2]
278 const char* kPhiForm = " 1 + 2*[1]*cos(x-[0]) + 2*[2]*cos(2*(x-[0])) ";
279 fPhiFormula = new TF1("gevsimPhi", kPhiForm, 0, 2*TMath::Pi());
281 fPhiFormula->SetParNames("psi", "directed", "elliptic");
282 fPhiFormula->SetParameters(0., 0., 0.);
284 fPhiFormula->SetNpx(180);
288 //////////////////////////////////////////////////////////////////////////////////
290 void AliGenGeVSim::AddParticleType(AliGeVSimParticle *part) {
292 // Adds new type of particles.
294 // Parameters are defeined in AliGeVSimParticle object
295 // This method has to be called for every particle type
298 if (fPartTypes == NULL)
299 fPartTypes = new TObjArray();
301 fPartTypes->Add(part);
304 //////////////////////////////////////////////////////////////////////////////////
306 void AliGenGeVSim::SetMultTotal(Bool_t isTotal) {
311 fIsMultTotal = isTotal;
314 //////////////////////////////////////////////////////////////////////////////////
316 Float_t AliGenGeVSim::FindScaler(Int_t paramId, Int_t pdg) {
319 // Finds Scallar for a given parameter.
320 // Function used in event-by-event mode.
322 // There are two types of scallars: deterministic and random
323 // and they can work on either global or particle type level.
324 // For every variable there are four scallars defined.
326 // Scallars are named as folowa
327 // deterministic global level : "gevsimParam" (eg. "gevsimTemp")
328 // deterinistig type level : "gevsimPdgParam" (eg. "gevsim211Mult")
329 // random global level : "gevsimParamRndm" (eg. "gevsimMultRndm")
330 // random type level : "gevsimPdgParamRndm" (eg. "gevsim-211V2Rndm");
332 // Pdg - code of a particle type in PDG standard (see: http://pdg.lbl.gov)
333 // Param - parameter name. Allowed parameters:
335 // "Temp" - inverse slope parameter
336 // "SigmaY" - rapidity width - for model (1) only
337 // "ExpVel" - expansion velocity - for model (4) only
338 // "V1" - directed flow
339 // "V2" - elliptic flow
340 // "Mult" - multiplicity
344 static const char* params[] = {"Temp", "SigmaY", "ExpVel", "V1", "V2", "Mult"};
345 static const char* ending[] = {"", "Rndm"};
347 static const char* patt1 = "gevsim%s%s";
348 static const char* patt2 = "gevsim%d%s%s";
355 // Scaler evoluation: i - global/local, j - determ/random
359 for (i=0; i<2; i++) {
360 for (j=0; j<2; j++) {
364 if (i == 0) sprintf(buffer, patt1, params[paramId], ending[j]);
365 else sprintf(buffer, patt2, pdg, params[paramId], ending[j]);
367 form = (TF1 *)gROOT->GetFunction(buffer);
370 if (j == 0) scaler *= form->Eval(gAlice->GetEvNumber());
372 form->SetParameter(0, gAlice->GetEvNumber());
373 scaler *= form->GetRandom();
382 //////////////////////////////////////////////////////////////////////////////////
384 void AliGenGeVSim::DetermineReactionPlane() {
386 // private function used by Generate()
388 // Retermines Reaction Plane angle and set this value
389 // as a parameter [0] in fPhiFormula
391 // Note: if "gevsimPsiRndm" function is found it override both
392 // "gevsimPhi" function and initial fPsi value
398 form = (TF1 *)gROOT->GetFunction("gevsimPsi");
399 if (form) fPsi = form->Eval(gAlice->GetEvNumber()) * TMath::Pi() / 180;
402 form = (TF1 *)gROOT->GetFunction("gevsimPsiRndm");
403 if (form) fPsi = form->GetRandom() * TMath::Pi() / 180;
406 cout << "Psi = " << fPsi << "\t" << (Int_t)(fPsi*180./TMath::Pi()) << endl;
408 fPhiFormula->SetParameter(0, fPsi);
411 //////////////////////////////////////////////////////////////////////////////////
413 Float_t AliGenGeVSim::GetdNdYToTotal() {
415 // Private, helper function used by Generate()
417 // Returns total multiplicity to dN/dY ration using current distribution.
418 // The function have to be called after setting distribution and its
419 // parameters (like temperature).
423 const Double_t kMaxPt = 3.0, kMaxY = 2.;
427 integ = fYFormula->Integral(-kMaxY, kMaxY);
428 mag = fYFormula->Eval(0);
432 // 2D formula standard or custom
434 if (fModel > 1 && fModel < 6) {
436 integ = ((TF2*)fCurrentForm)->Integral(0,kMaxPt, -kMaxY, kMaxY);
437 mag = ((TF2*)fCurrentForm)->Integral(0, kMaxPt, -0.1, 0.1) / 0.2;
445 integ = fHist[1]->Integral();
446 mag = fHist[0]->GetBinContent(fHist[0]->FindBin(0.));
447 mag /= fHist[0]->GetBinWidth(fHist[0]->FindBin(0.));
456 Int_t yBins = fPtYHist->GetNbinsY();
457 Int_t ptBins = fPtYHist->GetNbinsX();
459 integ = fPtYHist->Integral(0, ptBins, 0, yBins);
460 mag = fPtYHist->Integral(0, ptBins, (yBins/2)-1, (yBins/2)+1 ) / 2;
467 //////////////////////////////////////////////////////////////////////////////////
469 void AliGenGeVSim::SetFormula(Int_t pdg) {
471 // Private function used by Generate()
473 // Configure a formula for a given particle type and model Id (in fModel).
474 // If custom formula or histogram was selected the function tries
477 // The function implements naming conventions for custom distributions names
481 const char* msg[4] = {
482 "Custom Formula for Pt Y distribution not found [pdg = %d]",
483 "Histogram for Pt distribution not found [pdg = %d]",
484 "Histogram for Y distribution not found [pdg = %d]",
485 "HIstogram for Pt Y dostribution not found [pdg = %d]"
488 const char* pattern[8] = {
489 "gevsimDistPtY", "gevsimDist%dPtY",
490 "gevsimHistPt", "gevsimHist%dPt",
491 "gevsimHistY", "gevsimHist%dY",
492 "gevsimHistPtY", "gevsimHist%dPtY"
495 const char *where = "SetFormula";
498 if (fModel < 1 || fModel > 7)
499 Error("SetFormula", "Model Id (%d) out of range [1-7]", fModel);
504 if (fModel == 1) fCurrentForm = fPtFormula;
505 if (fModel > 1 && fModel < 5) fCurrentForm = fPtYFormula[fModel-2];
508 // custom model defined by a formula
513 fCurrentForm = (TF2*)gROOT->GetFunction(pattern[0]);
517 sprintf(buff, pattern[1], pdg);
518 fCurrentForm = (TF2*)gROOT->GetFunction(buff);
520 if (!fCurrentForm) Error(where, msg[0], pdg);
528 for (Int_t i=0; i<2; i++) {
531 fHist[i] = (TH1D*)gROOT->FindObject(pattern[2+2*i]);
535 sprintf(buff, pattern[3+2*i], pdg);
536 fHist[i] = (TH1D*)gROOT->FindObject(buff);
538 if (!fHist[i]) Error(where, msg[1+i], pdg);
548 fPtYHist = (TH2D*)gROOT->FindObject(pattern[6]);
552 sprintf(buff, pattern[7], pdg);
553 fPtYHist = (TH2D*)gROOT->FindObject(buff);
556 if (!fPtYHist) Error(where, msg[3], pdg);
561 //////////////////////////////////////////////////////////////////////////////////
563 void AliGenGeVSim:: AdjustFormula() {
566 // Adjust fomula bounds according to acceptance cuts.
568 // Since GeVSim is producing "thermal" particles Pt
569 // is cut at 3 GeV even when acceptance extends to grater momenta.
572 // If custom formula was provided function preserves
576 const Double_t kMaxPt = 3.0;
577 const Double_t kMaxY = 2.0;
578 Double_t minPt, maxPt, minY, maxY;
581 if (fModel > 4) return;
584 if (TestBit(kPtRange) && fPtMax < kMaxPt ) maxPt = fPtMax;
588 if (TestBit(kPtRange)) minPt = fPtMin;
591 if (TestBit(kPtRange) && fPtMin > kMaxPt )
592 Warning("Acceptance", "Minimum Pt (%3.2f GeV) greater that 3.0 GeV ", fPtMin);
595 if (TestBit(kMomentumRange) && fPtMax < maxPt) maxPt = fPtMax;
597 // max and min rapidity
598 if (TestBit(kYRange)) {
609 fPtFormula->SetRange(fPtMin, maxPt);
610 fYFormula->SetRange(fYMin, fYMax);
614 ((TF2*)fCurrentForm)->SetRange(minPt, minY, maxPt, maxY);
618 if (TestBit(kPhiRange))
619 fPhiFormula->SetRange(fPhiMin, fPhiMax);
623 //////////////////////////////////////////////////////////////////////////////////
625 void AliGenGeVSim::GetRandomPtY(Double_t &pt, Double_t &y) {
627 // Private function used by Generate()
629 // Returns random values of Pt and Y corresponding to selected
634 pt = fPtFormula->GetRandom();
635 y = fYFormula->GetRandom();
639 if (fModel > 1 && fModel < 6) {
640 ((TF2*)fCurrentForm)->GetRandom2(pt, y);
645 pt = fHist[0]->GetRandom();
646 y = fHist[1]->GetRandom();
650 fPtYHist->GetRandom2(pt, y);
655 //////////////////////////////////////////////////////////////////////////////////
657 void AliGenGeVSim::Init() {
659 // Standard AliGenerator initializer.
664 //////////////////////////////////////////////////////////////////////////////////
666 void AliGenGeVSim::Generate() {
668 // Standard AliGenerator function
669 // This function do actual job and puts particles on stack.
672 PDG_t pdg; // particle type
673 Float_t mass; // particle mass
674 Float_t orgin[3] = {0,0,0}; // particle orgin [cm]
675 Float_t polar[3] = {0,0,0}; // polarisation
676 Float_t time = 0; // time of creation
678 Float_t multiplicity = 0;
679 Bool_t isMultTotal = kTRUE;
682 Float_t directedScaller = 1., ellipticScaller = 1.;
684 TLorentzVector *v = new TLorentzVector(0,0,0,0);
686 const Int_t kParent = -1;
691 orgin[0] = fVertex[0];
692 orgin[1] = fVertex[1];
693 orgin[2] = fVertex[2];
696 // Particle params database
698 TDatabasePDG *db = TDatabasePDG::Instance();
700 AliGeVSimParticle *partType;
702 Int_t nType, nParticle, nParam;
703 const Int_t kNParams = 6;
705 // reaction plane determination and model
706 DetermineReactionPlane();
708 // loop over particle types
710 for (nType = 0; nType < fPartTypes->GetEntries(); nType++) {
712 partType = (AliGeVSimParticle *)fPartTypes->At(nType);
714 pdg = (PDG_t)partType->GetPdgCode();
715 type = db->GetParticle(pdg);
718 fModel = partType->GetModel();
720 fCurrentForm->SetParameter("mass", mass);
723 // Evaluation of parameters - loop over parameters
725 for (nParam = 0; nParam < kNParams; nParam++) {
727 paramScaler = FindScaler(nParam, pdg);
730 fCurrentForm->SetParameter("temperature", paramScaler * partType->GetTemperature());
732 if (nParam == 1 && fModel == 1)
733 fYFormula->SetParameter("sigmaY", paramScaler * partType->GetSigmaY());
735 if (nParam == 2 && fModel == 4) {
737 Double_t totalExpVal = paramScaler * partType->GetExpansionVelocity();
739 if (totalExpVal == 0.0) totalExpVal = 0.0001;
740 if (totalExpVal == 1.0) totalExpVal = 9.9999;
742 fCurrentForm->SetParameter("expVel", totalExpVal);
747 if (nParam == 3) directedScaller = paramScaler;
748 if (nParam == 4) ellipticScaller = paramScaler;
754 if (partType->IsMultForced()) isMultTotal = partType->IsMultTotal();
755 else isMultTotal = fIsMultTotal;
757 multiplicity = paramScaler * partType->GetMultiplicity();
758 multiplicity *= (isMultTotal)? 1 : GetdNdYToTotal();
762 // Flow defined on the particle type level (not parameterised)
763 if (partType->IsFlowSimple()) {
764 fPhiFormula->SetParameter(1, partType->GetDirectedFlow(0,0) * directedScaller);
765 fPhiFormula->SetParameter(2, partType->GetEllipticFlow(0,0) * ellipticScaller);
771 Info("Generate","PDG = %d \t Mult = %d", pdg, (Int_t)multiplicity);
773 // loop over particles
776 while (nParticle < multiplicity) {
778 Double_t pt, y, phi; // momentum in [pt,y,phi]
779 Float_t p[3] = {0,0,0}; // particle momentum
783 // phi distribution configuration when differential flow defined
784 // to be optimised in future release
786 if (!partType->IsFlowSimple()) {
787 fPhiFormula->SetParameter(1, partType->GetDirectedFlow(pt,y) * directedScaller);
788 fPhiFormula->SetParameter(2, partType->GetEllipticFlow(pt,y) * ellipticScaller);
791 phi = fPhiFormula->GetRandom();
793 if (!isMultTotal) nParticle++;
794 if (fModel > 4 && !CheckPtYPhi(pt,y,phi) ) continue;
796 // coordinate transformation
797 v->SetPtEtaPhiM(pt, y, phi, mass);
803 // momentum range test
804 if ( !CheckAcceptance(p) ) continue;
806 // putting particle on the stack
808 PushTrack(fTrackIt, kParent, pdg, p, orgin, polar, time, kPPrimary, id, fTrackIt);
809 if (isMultTotal) nParticle++;
813 // prepare and store header
815 AliGenGeVSimEventHeader *header = new AliGenGeVSimEventHeader("GeVSim header");
816 TArrayF eventVertex(3,orgin);
818 header->SetPrimaryVertex(eventVertex);
819 header->SetEventPlane(fPsi);
820 header->SetEllipticFlow(fPhiFormula->GetParameter(2));
822 gAlice->SetGenEventHeader(header);
827 //////////////////////////////////////////////////////////////////////////////////