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 // AliGenGeVSim is a class implementing GeVSim event generator.
20 // GeVSim is a simple Monte-Carlo event generator for testing detector and
21 // algorythm performance especialy concerning flow and event-by-event studies
23 // In this event generator particles are generated from thermal distributions
24 // without any dynamics and addicional constrains. Distribution parameters like
25 // multiplicity, particle type yields, inverse slope parameters, flow coeficients
26 // and expansion velocities are expleicite defined by the user.
28 // GeVSim contains four thermal distributions the same as
29 // MevSim event generator developed for STAR experiment.
31 // In addition custom distributions can be used be the mean
32 // either two dimensional formula (TF2), a two dimensional histogram or
33 // two one dimensional histograms.
35 // Azimuthal distribution is deconvoluted from (Pt,Y) distribution
36 // and is described by two Fourier coefficients representing
37 // Directed and Elliptic flow.
39 ////////////////////////////////////////////////////////////////////////////////
41 // To apply flow to event ganerated by an arbitraly event generator
42 // refer to AliGenAfterBurnerFlow class.
44 ////////////////////////////////////////////////////////////////////////////////
46 // For examples, parameters and testing macros refer to:
47 // http:/home.cern.ch/radomski
49 // for more detailed description refer to ALICE NOTE
50 // "GeVSim Monte-Carlo Event Generator"
51 // S.Radosmki, P. Foka.
54 // Sylwester Radomski,
59 ////////////////////////////////////////////////////////////////////////////////
61 // Updated and revised: September 2002, S. Radomski, GSI
63 ////////////////////////////////////////////////////////////////////////////////
66 #include <Riostream.h>
72 #include <TObjArray.h>
74 #include <TParticle.h>
78 #include "AliGeVSimParticle.h"
79 #include "AliGenGeVSim.h"
80 #include "AliGenGeVSimEventHeader.h"
81 #include "AliGenerator.h"
85 ClassImp(AliGenGeVSim);
87 //////////////////////////////////////////////////////////////////////////////////
89 AliGenGeVSim::AliGenGeVSim() : AliGenerator(-1) {
91 // Default constructor
98 for (Int_t i=0; i<4; i++)
103 //////////////////////////////////////////////////////////////////////////////////
105 AliGenGeVSim::AliGenGeVSim(Float_t psi, Bool_t isMultTotal) : AliGenerator(-1) {
107 // Standard Constructor.
109 // models - thermal model to be used:
110 // 1 - deconvoluted pt and Y source
111 // 2,3 - thermalized sphericaly symetric sources
112 // 4 - thermalized source with expansion
113 // 5 - custom model defined in TF2 object named "gevsimPtY"
114 // 6 - custom model defined by two 1D histograms
115 // 7 - custom model defined by 2D histogram
117 // psi - reaction plane in degrees
118 // isMultTotal - multiplicity mode
119 // kTRUE - total multiplicity (default)
120 // kFALSE - dN/dY at midrapidity
123 // checking consistancy
125 if (psi < 0 || psi > 360 )
126 Error ("AliGenGeVSim", "Reaction plane angle ( %d )out of range [0..360]", psi);
128 fPsi = psi * TMath::Pi() / 180. ;
129 fIsMultTotal = isMultTotal;
133 fPartTypes = new TObjArray();
137 //////////////////////////////////////////////////////////////////////////////////
139 AliGenGeVSim::~AliGenGeVSim() {
141 // Default Destructor
143 // Removes TObjArray keeping list of registered particle types
146 if (fPartTypes != NULL) delete fPartTypes;
150 //////////////////////////////////////////////////////////////////////////////////
152 Bool_t AliGenGeVSim::CheckPtYPhi(Float_t pt, Float_t y, Float_t phi) const {
154 // private function used by Generate()
156 // Check bounds of Pt, Rapidity and Azimuthal Angle of a track
157 // Used only when generating particles from a histogram
160 if ( TestBit(kPtRange) && ( pt < fPtMin || pt > fPtMax )) return kFALSE;
161 if ( TestBit(kPhiRange) && ( phi < fPhiMin || phi > fPhiMax )) return kFALSE;
162 if ( TestBit(kYRange) && ( y < fYMin || y > fYMax )) return kFALSE;
167 //////////////////////////////////////////////////////////////////////////////////
169 Bool_t AliGenGeVSim::CheckAcceptance(Float_t p[3]) {
171 // private function used by Generate()
173 // Check bounds of a total momentum and theta of a track
176 if ( TestBit(kThetaRange) ) {
178 Double_t theta = TMath::ATan2( TMath::Sqrt(p[0]*p[0]+p[1]*p[1]), p[2]);
179 if ( theta < fThetaMin || theta > fThetaMax ) return kFALSE;
183 if ( TestBit(kMomentumRange) ) {
185 Double_t momentum = TMath::Sqrt(p[0]*p[0] + p[1]*p[1] + p[2]*p[2]);
186 if ( momentum < fPMin || momentum > fPMax) return kFALSE;
192 //////////////////////////////////////////////////////////////////////////////////
194 void AliGenGeVSim::InitFormula() {
198 // Initalizes formulas used in GeVSim.
199 // Manages strings and creates TFormula objects from strings
202 // Deconvoluted Pt Y formula
204 // ptForm: pt -> x , mass -> [0] , temperature -> [1]
205 // y Form: y -> x , sigmaY -> [0]
207 const char* kPtForm = " x * exp( -sqrt([0]*[0] + x*x) / [1] )";
208 const char* kYForm = " exp ( - x*x / (2 * [0]*[0] ) )";
210 fPtFormula = new TF1("gevsimPt", kPtForm, 0, 3);
211 fYFormula = new TF1("gevsimRapidity", kYForm, -3, 3);
213 fPtFormula->SetParNames("mass", "temperature");
214 fPtFormula->SetParameters(1., 1.);
216 fYFormula->SetParName(0, "sigmaY");
217 fYFormula->SetParameter(0, 1.);
220 fPtFormula->SetNpx(100);
221 fYFormula->SetNpx(100);
227 // mass -> [0] , temperature -> [1] , expansion velocity -> [2]
230 const char *kFormE = " ( sqrt([0]*[0] + x*x) * cosh(y) ) ";
231 const char *kFormG = " ( 1 / sqrt( 1 - [2]*[2] ) ) ";
232 const char *kFormYp = "( [2]*sqrt(([0]*[0]+x*x)*cosh(y)*cosh(y)-[0]*[0])/([1]*sqrt(1-[2]*[2]))) ";
234 const char* kFormula[3] = {
235 " x * %s * exp( -%s / [1]) ",
236 " (x * %s) / ( exp( %s / [1]) - 1 ) ",
237 " x*%s*exp(-%s*%s/[1])*((sinh(%s)/%s)+([1]/(%s*%s))*(sinh(%s)/%s-cosh(%s)))"
240 const char* kParamNames[3] = {"mass", "temperature", "expVel"};
244 sprintf(buffer, kFormula[0], kFormE, kFormE);
245 fPtYFormula[0] = new TF2("gevsimPtY_2", buffer, 0, 3, -2, 2);
247 sprintf(buffer, kFormula[1], kFormE, kFormE);
248 fPtYFormula[1] = new TF2("gevsimPtY_3", buffer, 0, 3, -2, 2);
250 sprintf(buffer, kFormula[2], kFormE, kFormG, kFormE, kFormYp, kFormYp, kFormG, kFormE, kFormYp, kFormYp, kFormYp);
251 fPtYFormula[2] = new TF2("gevsimPtY_4", buffer, 0, 3, -2, 2);
256 // setting names & initialisation
259 for (i=0; i<3; i++) {
261 fPtYFormula[i]->SetNpx(100); // step 30 MeV
262 fPtYFormula[i]->SetNpy(100); //
264 for (j=0; j<3; j++) {
266 if ( i != 2 && j == 2 ) continue; // ExpVel
267 fPtYFormula[i]->SetParName(j, kParamNames[j]);
268 fPtYFormula[i]->SetParameter(j, 0.5);
275 // Psi -> [0] , Direct Flow -> [1] , Ellipticla Flow -> [2]
277 const char* kPhiForm = " 1 + 2*[1]*cos(x-[0]) + 2*[2]*cos(2*(x-[0])) ";
278 fPhiFormula = new TF1("gevsimPhi", kPhiForm, 0, 2*TMath::Pi());
280 fPhiFormula->SetParNames("psi", "directed", "elliptic");
281 fPhiFormula->SetParameters(0., 0., 0.);
283 fPhiFormula->SetNpx(180);
287 //////////////////////////////////////////////////////////////////////////////////
289 void AliGenGeVSim::AddParticleType(AliGeVSimParticle *part) {
291 // Adds new type of particles.
293 // Parameters are defeined in AliGeVSimParticle object
294 // This method has to be called for every particle type
297 if (fPartTypes == NULL)
298 fPartTypes = new TObjArray();
300 fPartTypes->Add(part);
303 //////////////////////////////////////////////////////////////////////////////////
305 void AliGenGeVSim::SetMultTotal(Bool_t isTotal) {
310 fIsMultTotal = isTotal;
313 //////////////////////////////////////////////////////////////////////////////////
315 Float_t AliGenGeVSim::FindScaler(Int_t paramId, Int_t pdg) {
318 // Finds Scallar for a given parameter.
319 // Function used in event-by-event mode.
321 // There are two types of scallars: deterministic and random
322 // and they can work on either global or particle type level.
323 // For every variable there are four scallars defined.
325 // Scallars are named as folowa
326 // deterministic global level : "gevsimParam" (eg. "gevsimTemp")
327 // deterinistig type level : "gevsimPdgParam" (eg. "gevsim211Mult")
328 // random global level : "gevsimParamRndm" (eg. "gevsimMultRndm")
329 // random type level : "gevsimPdgParamRndm" (eg. "gevsim-211V2Rndm");
331 // Pdg - code of a particle type in PDG standard (see: http://pdg.lbl.gov)
332 // Param - parameter name. Allowed parameters:
334 // "Temp" - inverse slope parameter
335 // "SigmaY" - rapidity width - for model (1) only
336 // "ExpVel" - expansion velocity - for model (4) only
337 // "V1" - directed flow
338 // "V2" - elliptic flow
339 // "Mult" - multiplicity
343 static const char* params[] = {"Temp", "SigmaY", "ExpVel", "V1", "V2", "Mult"};
344 static const char* ending[] = {"", "Rndm"};
346 static const char* patt1 = "gevsim%s%s";
347 static const char* patt2 = "gevsim%d%s%s";
354 // Scaler evoluation: i - global/local, j - determ/random
358 for (i=0; i<2; i++) {
359 for (j=0; j<2; j++) {
363 if (i == 0) sprintf(buffer, patt1, params[paramId], ending[j]);
364 else sprintf(buffer, patt2, pdg, params[paramId], ending[j]);
366 form = (TF1 *)gROOT->GetFunction(buffer);
369 if (j == 0) scaler *= form->Eval(gAlice->GetEvNumber());
371 form->SetParameter(0, gAlice->GetEvNumber());
372 scaler *= form->GetRandom();
381 //////////////////////////////////////////////////////////////////////////////////
383 void AliGenGeVSim::DetermineReactionPlane() {
385 // private function used by Generate()
387 // Retermines Reaction Plane angle and set this value
388 // as a parameter [0] in fPhiFormula
390 // Note: if "gevsimPsiRndm" function is found it override both
391 // "gevsimPhi" function and initial fPsi value
397 form = (TF1 *)gROOT->GetFunction("gevsimPsi");
398 if (form) fPsi = form->Eval(gAlice->GetEvNumber()) * TMath::Pi() / 180;
401 form = (TF1 *)gROOT->GetFunction("gevsimPsiRndm");
402 if (form) fPsi = form->GetRandom() * TMath::Pi() / 180;
405 cout << "Psi = " << fPsi << "\t" << (Int_t)(fPsi*180./TMath::Pi()) << endl;
407 fPhiFormula->SetParameter(0, fPsi);
410 //////////////////////////////////////////////////////////////////////////////////
412 Float_t AliGenGeVSim::GetdNdYToTotal() {
414 // Private, helper function used by Generate()
416 // Returns total multiplicity to dN/dY ration using current distribution.
417 // The function have to be called after setting distribution and its
418 // parameters (like temperature).
422 const Double_t kMaxPt = 3.0, kMaxY = 2.;
426 integ = fYFormula->Integral(-kMaxY, kMaxY);
427 mag = fYFormula->Eval(0);
431 // 2D formula standard or custom
433 if (fModel > 1 && fModel < 6) {
435 integ = ((TF2*)fCurrentForm)->Integral(0,kMaxPt, -kMaxY, kMaxY);
436 mag = ((TF2*)fCurrentForm)->Integral(0, kMaxPt, -0.1, 0.1) / 0.2;
444 integ = fHist[1]->Integral();
445 mag = fHist[0]->GetBinContent(fHist[0]->FindBin(0.));
446 mag /= fHist[0]->GetBinWidth(fHist[0]->FindBin(0.));
455 Int_t yBins = fPtYHist->GetNbinsY();
456 Int_t ptBins = fPtYHist->GetNbinsX();
458 integ = fPtYHist->Integral(0, ptBins, 0, yBins);
459 mag = fPtYHist->Integral(0, ptBins, (yBins/2)-1, (yBins/2)+1 ) / 2;
466 //////////////////////////////////////////////////////////////////////////////////
468 void AliGenGeVSim::SetFormula(Int_t pdg) {
470 // Private function used by Generate()
472 // Configure a formula for a given particle type and model Id (in fModel).
473 // If custom formula or histogram was selected the function tries
476 // The function implements naming conventions for custom distributions names
480 const char* msg[4] = {
481 "Custom Formula for Pt Y distribution not found [pdg = %d]",
482 "Histogram for Pt distribution not found [pdg = %d]",
483 "Histogram for Y distribution not found [pdg = %d]",
484 "HIstogram for Pt Y dostribution not found [pdg = %d]"
487 const char* pattern[8] = {
488 "gevsimDistPtY", "gevsimDist%dPtY",
489 "gevsimHistPt", "gevsimHist%dPt",
490 "gevsimHistY", "gevsimHist%dY",
491 "gevsimHistPtY", "gevsimHist%dPtY"
494 const char *where = "SetFormula";
497 if (fModel < 1 || fModel > 7)
498 Error("SetFormula", "Model Id (%d) out of range [1-7]", fModel);
503 if (fModel == 1) fCurrentForm = fPtFormula;
504 if (fModel > 1 && fModel < 5) fCurrentForm = fPtYFormula[fModel-2];
507 // custom model defined by a formula
512 fCurrentForm = (TF2*)gROOT->GetFunction(pattern[0]);
516 sprintf(buff, pattern[1], pdg);
517 fCurrentForm = (TF2*)gROOT->GetFunction(buff);
519 if (!fCurrentForm) Error(where, msg[0], pdg);
527 for (Int_t i=0; i<2; i++) {
530 fHist[i] = (TH1D*)gROOT->FindObject(pattern[2+2*i]);
534 sprintf(buff, pattern[3+2*i], pdg);
535 fHist[i] = (TH1D*)gROOT->FindObject(buff);
537 if (!fHist[i]) Error(where, msg[1+i], pdg);
547 fPtYHist = (TH2D*)gROOT->FindObject(pattern[6]);
551 sprintf(buff, pattern[7], pdg);
552 fPtYHist = (TH2D*)gROOT->FindObject(buff);
555 if (!fPtYHist) Error(where, msg[3], pdg);
560 //////////////////////////////////////////////////////////////////////////////////
562 void AliGenGeVSim:: AdjustFormula() {
565 // Adjust fomula bounds according to acceptance cuts.
567 // Since GeVSim is producing "thermal" particles Pt
568 // is cut at 3 GeV even when acceptance extends to grater momenta.
571 // If custom formula was provided function preserves
575 const Double_t kMaxPt = 3.0;
576 const Double_t kMaxY = 2.0;
577 Double_t minPt, maxPt, minY, maxY;
580 if (fModel > 4) return;
583 if (TestBit(kPtRange) && fPtMax < kMaxPt ) maxPt = fPtMax;
587 if (TestBit(kPtRange)) minPt = fPtMin;
590 if (TestBit(kPtRange) && fPtMin > kMaxPt )
591 Warning("Acceptance", "Minimum Pt (%3.2f GeV) greater that 3.0 GeV ", fPtMin);
594 if (TestBit(kMomentumRange) && fPtMax < maxPt) maxPt = fPtMax;
596 // max and min rapidity
597 if (TestBit(kYRange)) {
608 fPtFormula->SetRange(fPtMin, maxPt);
609 fYFormula->SetRange(fYMin, fYMax);
613 ((TF2*)fCurrentForm)->SetRange(minPt, minY, maxPt, maxY);
617 if (TestBit(kPhiRange))
618 fPhiFormula->SetRange(fPhiMin, fPhiMax);
622 //////////////////////////////////////////////////////////////////////////////////
624 void AliGenGeVSim::GetRandomPtY(Double_t &pt, Double_t &y) {
626 // Private function used by Generate()
628 // Returns random values of Pt and Y corresponding to selected
633 pt = fPtFormula->GetRandom();
634 y = fYFormula->GetRandom();
638 if (fModel > 1 && fModel < 6) {
639 ((TF2*)fCurrentForm)->GetRandom2(pt, y);
644 pt = fHist[0]->GetRandom();
645 y = fHist[1]->GetRandom();
649 fPtYHist->GetRandom2(pt, y);
654 //////////////////////////////////////////////////////////////////////////////////
656 void AliGenGeVSim::Init() {
658 // Standard AliGenerator initializer.
663 //////////////////////////////////////////////////////////////////////////////////
665 void AliGenGeVSim::Generate() {
667 // Standard AliGenerator function
668 // This function do actual job and puts particles on stack.
671 PDG_t pdg; // particle type
672 Float_t mass; // particle mass
673 Float_t orgin[3] = {0,0,0}; // particle orgin [cm]
674 Float_t polar[3] = {0,0,0}; // polarisation
675 Float_t time = 0; // time of creation
677 Float_t multiplicity = 0;
678 Bool_t isMultTotal = kTRUE;
681 Float_t directedScaller = 1., ellipticScaller = 1.;
683 TLorentzVector *v = new TLorentzVector(0,0,0,0);
685 const Int_t kParent = -1;
690 orgin[0] = fVertex[0];
691 orgin[1] = fVertex[1];
692 orgin[2] = fVertex[2];
695 // Particle params database
697 TDatabasePDG *db = TDatabasePDG::Instance();
699 AliGeVSimParticle *partType;
701 Int_t nType, nParticle, nParam;
702 const Int_t kNParams = 6;
704 // reaction plane determination and model
705 DetermineReactionPlane();
707 // loop over particle types
709 for (nType = 0; nType < fPartTypes->GetEntries(); nType++) {
711 partType = (AliGeVSimParticle *)fPartTypes->At(nType);
713 pdg = (PDG_t)partType->GetPdgCode();
714 type = db->GetParticle(pdg);
717 fModel = partType->GetModel();
719 fCurrentForm->SetParameter("mass", mass);
722 // Evaluation of parameters - loop over parameters
724 for (nParam = 0; nParam < kNParams; nParam++) {
726 paramScaler = FindScaler(nParam, pdg);
729 fCurrentForm->SetParameter("temperature", paramScaler * partType->GetTemperature());
731 if (nParam == 1 && fModel == 1)
732 fYFormula->SetParameter("sigmaY", paramScaler * partType->GetSigmaY());
734 if (nParam == 2 && fModel == 4) {
736 Double_t totalExpVal = paramScaler * partType->GetExpansionVelocity();
738 if (totalExpVal == 0.0) totalExpVal = 0.0001;
739 if (totalExpVal == 1.0) totalExpVal = 9.9999;
741 fCurrentForm->SetParameter("expVel", totalExpVal);
746 if (nParam == 3) directedScaller = paramScaler;
747 if (nParam == 4) ellipticScaller = paramScaler;
753 if (partType->IsMultForced()) isMultTotal = partType->IsMultTotal();
754 else isMultTotal = fIsMultTotal;
756 multiplicity = paramScaler * partType->GetMultiplicity();
757 multiplicity *= (isMultTotal)? 1 : GetdNdYToTotal();
761 // Flow defined on the particle type level (not parameterised)
762 if (partType->IsFlowSimple()) {
763 fPhiFormula->SetParameter(1, partType->GetDirectedFlow(0,0) * directedScaller);
764 fPhiFormula->SetParameter(2, partType->GetEllipticFlow(0,0) * ellipticScaller);
770 Info("Generate","PDG = %d \t Mult = %d", pdg, (Int_t)multiplicity);
772 // loop over particles
775 while (nParticle < multiplicity) {
777 Double_t pt, y, phi; // momentum in [pt,y,phi]
778 Float_t p[3] = {0,0,0}; // particle momentum
782 // phi distribution configuration when differential flow defined
783 // to be optimised in future release
785 if (!partType->IsFlowSimple()) {
786 fPhiFormula->SetParameter(1, partType->GetDirectedFlow(pt,y) * directedScaller);
787 fPhiFormula->SetParameter(2, partType->GetEllipticFlow(pt,y) * ellipticScaller);
790 phi = fPhiFormula->GetRandom();
792 if (!isMultTotal) nParticle++;
793 if (fModel > 4 && !CheckPtYPhi(pt,y,phi) ) continue;
795 // coordinate transformation
796 v->SetPtEtaPhiM(pt, y, phi, mass);
802 // momentum range test
803 if ( !CheckAcceptance(p) ) continue;
805 // putting particle on the stack
807 PushTrack(fTrackIt, kParent, pdg, p, orgin, polar, time, kPPrimary, id, fTrackIt);
808 if (isMultTotal) nParticle++;
812 // prepare and store header
814 AliGenGeVSimEventHeader *header = new AliGenGeVSimEventHeader("GeVSim header");
815 TArrayF eventVertex(3,orgin);
817 header->SetPrimaryVertex(eventVertex);
818 header->SetEventPlane(fPsi);
819 header->SetEllipticFlow(fPhiFormula->GetParameter(2));
821 gAlice->SetGenEventHeader(header);
826 //////////////////////////////////////////////////////////////////////////////////