]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Base/AliFlowOnTheFlyEventGenerator.h
adds void AliFlowEventSimple::ShuffleTracks() to allow for explicit reshuffling of...
[u/mrichter/AliRoot.git] / PWG / FLOW / Base / AliFlowOnTheFlyEventGenerator.h
1 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
2 * See cxx source for full Copyright notice */
3 /* $Id$ */
4
5 // AliFlowOnTheFlyEventGenerator
6 // origin: Redmer Alexander Bertens (rbertens@cern.ch, rbertens@nikhef.nl, rbertens@uu.nl)
7
8 #ifndef ALIFLOWONTHEFLYEVENTGENERATOR_H
9 #define ALIFLOWONTHEFLYEVENTGENERATOR_H
10
11 #include "TObject.h"
12 #include "TClonesArray.h"
13 #include "TRandom.h"
14 #include "TH1F.h"
15 #include "TH2F.h"
16 #include "TF1.h"
17 #include "TParticle.h"
18
19 class AliFlowEventSimple;
20 class TVirtualMCDecayer;
21
22 class AliFlowOnTheFlyEventGenerator : public TObject {
23     public:
24         // general
25         AliFlowOnTheFlyEventGenerator();                
26         AliFlowOnTheFlyEventGenerator(Bool_t qa, Int_t ff, Int_t mult, TVirtualMCDecayer* decayer, Bool_t a, Bool_t b, Bool_t c, Bool_t d);
27         virtual ~AliFlowOnTheFlyEventGenerator();
28         class NaiveFlowAndSpectrumGenerator : public TObject { // small nested helper class. use methods of AliFlowOnTheFlyEventGenerator to get access to members
29             public: 
30                 NaiveFlowAndSpectrumGenerator(Short_t pdg, Bool_t qa, Int_t ff) : fPdg(pdg), fQA(qa), fFF(ff), fpt(0), fv2(0), fv3(0), fQApt(0), fQAv2(0), fQAv3(0) {
31                     TParticle* t = new TParticle(); t->SetPdgCode(fPdg);
32                     fpt = new TF1(Form("pt_%i", pdg), Form("x/TMath::Power(1+2*TMath::Sqrt(%f*%f+x*x),4)", t->GetMass(), t->GetMass()),0.,20);
33                     fv2 = new TF1(Form("v2_%i", pdg), Form("TMath::Log(x+1)*.2/TMath::Power(%f*%f+x*x,.2)", t->GetMass(), t->GetMass()), 0., 20);
34                     fv3 = new TF1(Form("v3_%i", pdg), Form("TMath::Log(x+1)*.1/TMath::Power(%f*%f+x*x,.2)", t->GetMass(), t->GetMass()), 0., 20);
35                     if(fQA) { fQApt = new TH1F(Form("pt_%i", fPdg),Form("pt_%i", fPdg),400,0,20);
36                               fQAv2 = new TH2F(Form("v2_%i", fPdg),Form("v2_%i", fPdg),400,0,20, 400, -.5, .5);
37                               fQAv3 = new TH2F(Form("v3_%i", fPdg),Form("v3_%i", fPdg),400,0,20, 400, -.5, .5);      }
38                     delete t;
39                 }
40                 virtual ~NaiveFlowAndSpectrumGenerator() {if(fpt) delete fpt; if(fv2) delete fv2; if(fv3) delete fv3; if(fQA) {delete fQApt; delete fQAv2; delete fQAv3;};}
41                 Short_t         GetPDGCode()            const {return fPdg;}
42                 Double_t        GetPt()                 const {double _pt(fpt->GetRandom()); if(fQA) fQApt->Fill(_pt); return _pt;}
43                 Double_t        GetV2(Double_t pt)      const {double _v2(fv2->Eval(pt)); if(fQA&&fFF==0) fQAv2->Fill(pt, _v2); return _v2;} 
44                 Double_t        GetV3(Double_t pt)      const {double _v3(fv3->Eval(pt)); if(fQA&&fFF==0) fQAv3->Fill(pt, _v3); return _v3;}
45                 TF1*            GetPtSpectrum()         const {return fpt;}                     // return pt fSpectrum
46                 TF1*            GetDifferentialV2()     const {return fv2;}                     // return fDifferential v2
47                 TF1*            GetDifferentialV3()     const {return fv3;}                     // return fDifferential v3
48                 TH1*            GetQAType(Int_t t)      const { if(t==0) return (TH1*)fQApt; 
49                                                                 if(t==1) return (TH1*)fQAv2; 
50                                                                 if(t==2) return (TH1*)fQAv3; 
51                                                                 return 0x0; }                   // return base pointer to QA histo class
52                 void            FillV2(Double_t p, Double_t v)  {fQAv2->Fill(p,v);}             // fill QA histo in case of fluctuations
53                 void            SetPtSpectrum(TF1* s)           {fpt = s;}                      // set custom fSpectrum
54                 void            SetDifferentialV2(TF1* v2)      {fv2 = v2; }                    // set custom fDifferential v2
55                 void            SetDifferentialV3(TF1* v3)      {fv3 = v3; }                    // set custom fDifferential v3
56             private:
57                 Short_t                 fPdg;                                                   // pdg value of track
58                 Bool_t                  fQA;                                                    // make fQA plots for all generated species
59                 Int_t                   fFF;                                                    // introduce e-by-e flow fluctuations
60                 TF1*                    fpt;                                                    // !pt fSpectrum
61                 TF1*                    fv2;                                                    // !fDifferential v2
62                 TF1*                    fv3;                                                    // !fDifferential v3
63                 TH1F*                   fQApt;                                                  // !pt fSpectrum for fQA
64                 TH2F*                   fQAv2;                                                  // !v2 for fQA
65                 TH2F*                   fQAv3;                                                  // !v3 for fQA
66                 NaiveFlowAndSpectrumGenerator(const NaiveFlowAndSpectrumGenerator& dummy);              // not implemented
67                 NaiveFlowAndSpectrumGenerator& operator =(const NaiveFlowAndSpectrumGenerator& dummy);  // not implemented
68         };      // end of NaiveFlowAndSpectrumGenerator
69         // access to some members of the hidden nested class
70         NaiveFlowAndSpectrumGenerator*  Find(Short_t pdg, Bool_t make);
71         TObjArray*                      GetGenerators()  {return fGenerators; }
72         void                            SetPtSpectrum(const char* func, Short_t pdg);
73         void                            SetPtDependentV2(const char* func, Short_t pdg);  
74         void                            SetPtDependentV3(const char* func, Short_t pdg);
75         TF1*                            GetPtSpectrum(Short_t pdg)              {return Find(pdg, kTRUE)->GetPtSpectrum();}
76         TF1*                            GetDifferentialV2(Short_t pdg)          {return Find(pdg, kTRUE)->GetDifferentialV2();}
77         TF1*                            GetDifferentialV3(Short_t pdg)          {return Find(pdg, kTRUE)->GetDifferentialV3();}
78         TH1*                            GetQAType(Short_t pdg, Int_t type)      {return Find(pdg, kTRUE)->GetQAType(type);}
79         // event generator 
80         void                    AddV2(TParticle* particle, Double_t v2, Double_t fluc);
81         void                    AddV2(TClonesArray* event);
82         void                    SetAfterBurnerPrecision(Double_t a, Int_t b)    { fPrecisionPhi = a; fMaxNumberOfIterations = b; }
83         void                    GenerateOnTheFlyTracks(Int_t mult, Int_t pid, TClonesArray* event, Double_t fluc);
84         void                    DecayOnTheFlyTracks(TClonesArray* event); 
85         void                    ForceGammaDecay(TClonesArray* arr, TParticle* part); 
86         AliFlowEventSimple*     GenerateOnTheFlyEvent(TClonesArray* event, Int_t nSpecies, Int_t species[], Int_t mult[], Int_t bg, Bool_t fluc);
87         void                    EmbedEvent(TClonesArray* embedMe)               {fEmbedMe = embedMe;}
88         AliFlowEventSimple*     ConvertTClonesToFlowEvent(TClonesArray* event, Int_t totalMultiplicity);
89         void                    AddV2Mothers(Bool_t b)                          { fAddV2Mothers = b; }
90         void                    AddV3Mothers(Bool_t b)                          { fAddV3Mothers = b; }
91         void                    AddV2Daughters(Bool_t b)                        { fAddV2Daughters = b; }
92         void                    AddV3Daughters(Bool_t b)                        { fAddV3Daughters = b; }
93         // 'bookkeeping'
94         void                    InitGenerators();
95         void                    PrintGenerators();
96         void                    DoGeneratorQA(Bool_t v2, Bool_t v3);
97
98     private:
99         // (transient) class members
100         TObjArray*              fGenerators;                    // array with generators
101         TClonesArray*           fEmbedMe;                       // tclones array for embedding
102         AliFlowEventSimple*     fFlowEvent;                     //! flow event simple for output
103         TVirtualMCDecayer*      fDecayer;                       // virtual decayer, needs to be set in macro
104         Bool_t                  fAddV2Mothers;                  // add v2 to mother tracks
105         Bool_t                  fAddV3Mothers;                  // add v3 to mother tracks
106         Bool_t                  fAddV2Daughters;                // add v2 to daughter tracks (not implemented)
107         Bool_t                  fAddV3Daughters;                // add v3 to daughter tracks (not implemented)
108         Double_t                fPsi2;                          // 2nd order symmetry plane
109         Double_t                fPsi3;                          // 3rd order symmetry plane  (not implemented)
110         Double_t                fPrecisionPhi;                  // afterburner convergence precision
111         Int_t                   fMaxNumberOfIterations;         // afterburner convergence precision
112         Bool_t                  fQA;                            // save qa histograms for all generated species
113         Int_t                   fFF;                            // introduce e-by-e flow fluctuations
114         // assignment operator and copy constructor
115         AliFlowOnTheFlyEventGenerator(const AliFlowOnTheFlyEventGenerator& dummy);              // not implemented
116         AliFlowOnTheFlyEventGenerator& operator =(const AliFlowOnTheFlyEventGenerator& dummy);  // not implemented
117
118         ClassDef(AliFlowOnTheFlyEventGenerator, 0)
119 };
120
121 #endif