]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskRhoVnModulation.h
Up from Redmer
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskRhoVnModulation.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 #ifndef ALIANALYSISTASKRHOVNMODULATION_H
6 #define ALIANALYSISTASKRHOVNMODULATION_H
7
8 #include <AliAnalysisTaskEmcalJet.h>
9 #include <AliEmcalJet.h>
10 #include <AliVEvent.h>
11 #include <AliVTrack.h>
12 #include <AliVCluster.h>
13 #include <TClonesArray.h>
14 #include <TMath.h>
15
16 class TF1;
17 class THF1;
18 class THF2;
19 class TProfile;
20 class TRandom3;
21
22 class AliAnalysisTaskRhoVnModulation : public AliAnalysisTaskEmcalJet
23 {
24     public:
25          // enumerators
26         enum fitModulationType  { kNoFit, kV2, kV3, kCombined, kUser, kFourierSeries }; // fit type
27         enum runModeType        { kLocal, kGrid };                      // run mode type
28         enum dataType           { kESD, kAOD, kESDMC, kAODMC };         // data type
29         enum detectorType       { kTPC, kVZEROA, kVZEROC };             // detector that was used
30         // constructors, destructor
31                                 AliAnalysisTaskRhoVnModulation();
32                                 AliAnalysisTaskRhoVnModulation(const char *name, runModeType type);
33         virtual                 ~AliAnalysisTaskRhoVnModulation();
34        
35         // setting up the task and technical aspects
36         Bool_t                  InitializeAnalysis();
37         virtual void            UserCreateOutputObjects();
38         virtual Bool_t          Run();
39         TH1F*                   BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c = -1);
40         TH2F*                   BookTH2F(const char* name, const char* x, const char* y, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t c = -1);
41         /* inline */    Double_t PhaseShift(Double_t x) const {  
42             while (x>=TMath::TwoPi())x-=TMath::TwoPi();
43             while (x<0.)x+=TMath::TwoPi();
44             return x; }
45         /* inline */    Double_t ChiSquare(Int_t ndf, Double_t x) const {
46             Double_t n(ndf/2.), denom(TMath::Power(2, n)*TMath::Gamma(n));
47             if (denom!=0)  return ((1./denom)*TMath::Power(x, n-1)*TMath::Exp(-x/2.)); 
48             return -999; }
49         /* inline */    Double_t RhoVal() const { return (fRho) ? fRho->GetVal(): -999.;}                 
50         /* inline */    Double_t RhoVal(Double_t phi, Double_t r, Double_t n) const {
51             if(!fFitModulation) return -999; // coverity
52             switch (fFitModulationType) {
53                 case kNoFit : return RhoVal();
54                 default : {
55                     Double_t denom(2*r*fFitModulation->GetParameter(0));
56                     return  (denom <= 0.) ? RhoVal() : n*(fFitModulation->Integral(phi-r, phi+r)/denom); 
57                 }
58             }
59         }
60         // setters - analysis setup
61         void                    SetDebugMode(Int_t d)                           {fDebug = d;}
62         void                    SetFillQAHistograms(Bool_t qa)                  {fFillQAHistograms = qa;}
63         void                    SetCentralityClasses(TArrayI* c)                {fCentralityClasses = c;}
64         void                    SetNameJetClones(const char* name)              {fNameJetClones = name; }
65         void                    SetNamePicoTrackClones(const char* name)        {fNamePicoTrackClones = name; }
66         void                    SetNameRho(const char* name)                    {fNameRho = name; }
67         void                    SetRandomSeed(TRandom3* r)                      {if (fRandom) delete fRandom; fRandom = r; }
68         void                    SetModulationFit(TF1* fit)                      {if (fFitModulation) delete fFitModulation;
69                                                                                  fFitModulation = fit; }
70         void                    SetModulationFitMinimumPvalue(Float_t p)        {fMinPvalue = p;}
71         void                    SetModulationFitType(fitModulationType type)    {fFitModulationType = type; }
72         void                    SetModulationFitOptions(TString opt)            {fFitModulationOptions = opt; }
73         void                    SetModulationFitDetector(detectorType type)     {fDetectorType = type; }
74         void                    SetRunModeType(runModeType type)                {fRunModeType = type; }
75         void                    SetAbsVertexZ(Float_t v)                        {fAbsVertexZ = v; }
76         void                    SetMinDistanceRctoLJ(Float_t m)                 {fMinDisanceRCtoLJ = m; }
77         void                    SetRandomConeRadius(Float_t r)                  {fRandomConeRadius = r; }
78         // getters
79         /* FIXME implement getters */
80         // 'trivial' helper calculations
81         void                    CalculateEventPlaneVZERO(Double_t vzero[2][2]) const;
82         void                    CalculateEventPlaneTPC(Double_t* tpc) const;
83         void                    CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi, AliEmcalJet* jet = 0x0, Bool_t randomize = 0) const;
84         // analysis details
85         void                    CorrectRho(Double_t* params, Double_t psi2, Double_t psi3) const;
86         // event and track selection
87         /* inline */    Bool_t PassesCuts(const AliVTrack* track) const {
88             if(!track) return kFALSE;
89             return kTRUE; }
90         Bool_t                  PassesCuts(AliVEvent* event);
91         Bool_t                  PassesCuts(const AliVCluster* track) const;
92         Bool_t                  PassesCuts(const AliEmcalJet* jet) const;
93         // filling histograms
94         void                    FillHistogramsAfterSubtraction(Double_t vzero[2][2], Double_t* tpc) const;
95         void                    FillTrackHistograms() const;
96         void                    FillClusterHistograms() const;
97         void                    FillCorrectedClusterHistograms() const;
98         void                    FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* tpc) const;
99         void                    FillRhoHistograms() const;
100         void                    FillDeltaPtHistograms(Double_t* tpc) const; 
101         void                    FillJetHistograms(Double_t vzero[2][2], Double_t* tpc) const;
102         void                    FillDeltaPhiHistograms(Double_t vzero[2][2], Double_t* tpc) const;
103         void                    FillQAHistograms(AliVTrack* vtrack) const;
104         void                    FillQAHistograms(AliVEvent* vevent);
105         virtual void            Terminate(Option_t* option);
106     private:
107         // analysis flags and settings
108         Int_t                   fDebug;                 // debug level (0 none, 1 fcn calls, 2 verbose)
109         Bool_t                  fInitialized;           //! is the analysis initialized?
110         Bool_t                  fFillQAHistograms;      // fill qa histograms
111         TArrayI*                fCentralityClasses;     //-> centrality classes (maximum 10)
112         // members
113         fitModulationType       fFitModulationType;     // fit modulation type
114         detectorType            fDetectorType;          // type of detector used for modulation fit
115         TString                 fFitModulationOptions;  // fit options for modulation fit
116         runModeType             fRunModeType;           // run mode type 
117         dataType                fDataType;              // datatype 
118         TRandom3*               fRandom;                //-> dont use gRandom to not interfere with other tasks
119         Int_t                   fMappedRunNumber;       //! mapped runnumer (for QA)
120         Int_t                   fInCentralitySelection; //! centrality bin
121         TF1*                    fFitModulation;         //-> modulation fit for rho
122         Float_t                 fMinPvalue;             // minimum value of p
123         const char*             fNameJetClones;         //! collection of tclones array with jets
124         const char*             fNamePicoTrackClones;   //! collection of tclones with pico tracks
125         const char*             fNameRho;               //! name of rho
126         // event cuts
127         Float_t                 fAbsVertexZ;            // cut on zvertex
128         // general qa histograms
129         TH1F*                   fHistCentrality;        //! accepted centrality
130         TH1F*                   fHistVertexz;           //! accepted verte
131         TH2F*                   fHistRunnumbersPhi;     //! run numbers averaged phi
132         TH2F*                   fHistRunnumbersEta;     //! run numbers averaged eta
133         // general settings
134         Float_t                 fMinDisanceRCtoLJ;      //! min distance between rc and leading jet
135         Float_t                 fRandomConeRadius;      //! radius of random cone
136         // transient object pointers
137         TList*                  fOutputList;            //! output list
138         TList*                  fOutputListGood;        //! output list for local analysis
139         TList*                  fOutputListBad;         //! output list for local analysis
140         TH1F*                   fHistAnalysisSummary;   //! analysis summary
141         TH1F*                   fHistSwap;              //! swap histogram
142         TProfile*               fProfVn;                //! extracted vn
143         // qa histograms for accepted pico tracks
144         TH1F*                   fHistPicoTrackPt[10];    //! pt of all charged tracks
145         TH2F*                   fHistPicoCat1[10];       //! pico tracks spd hit and refit
146         TH2F*                   fHistPicoCat2[10];       //! pico tracks wo spd hit w refit, constrained
147         TH2F*                   fHistPicoCat3[10];       //! pico tracks wo spd hit wo refit, constrained
148         // qa histograms for accepted emcal clusters
149         /* TH1F*                   fHistClusterPt[10];      //! pt uncorrected emcal clusters */
150         /* TH1F*                   fHistClusterPhi[10];     //! phi uncorrected emcal clusters */
151         /* TH1F*                   fHistClusterEta[10];     //! eta uncorrected emcal clusters */
152         //// qa histograms for accepted emcal clusters aftehadronic correction
153         /* TH1F*                   fHistClusterCorrPt[10];  //! pt corrected emcal clusters */
154         /* TH1F*                   fHistClusterCorrPhi[10]; //! phi corrected emcal clusters */
155         /* TH1F*                   fHistClusterCorrEta[10]; //! eta corrected emcal clusters */
156         // qa event planes
157         TProfile*               fHistPsi2;               //! Psi_2 estimates
158         TProfile*               fHistPsi2Spread;         //! spread between 
159         TH1F*                   fHistPsiVZEROA;          //! psi 2 from vzero a
160         TH1F*                   fHistPsiVZEROC;          //! psi 2 from vzero c
161         TH1F*                   fHistPsiTPC;             //! psi 2 from tpc
162         // background
163         TH1F*                   fHistRhoPackage[10];     //! rho as estimated by emcal jet package
164         TH1F*                   fHistRho[10];            //! background
165         TH2F*                   fHistRhoVsMult;          //! rho versus multiplicity
166         TH2F*                   fHistRhoVsCent;          //! rho veruss centrality
167         TH2F*                   fHistRhoAVsMult;         //! rho * A vs multiplicity for all jets
168         TH2F*                   fHistRhoAVsCent;         //! rho * A vs centrality for all jets
169         // delta pt distributions
170         TH2F*                   fHistRCPhiEta[10];              //! random cone eta and phi
171         TH2F*                   fHistRhoVsRCPt[10];             //! rho * A vs rcpt
172         TH1F*                   fHistRCPt[10];                  //! rcpt
173         TH2F*                   fHistDeltaPtDeltaPhi[10];       //! dpt vs dphi
174         TH2F*                   fHistRCPhiEtaExLJ[10];          //! random cone eta and phi, excl leading jet
175         TH2F*                   fHistRhoVsRCPtExLJ[10];         //! rho * A vs rcpt, excl leading jet
176         TH1F*                   fHistRCPtExLJ[10];              //! rcpt, excl leading jet
177         TH2F*                   fHistDeltaPtDeltaPhiExLJ[10];   //! dpt vs dphi, excl leading jet
178         TH2F*                   fHistRCPhiEtaRand[10];          //! random cone eta and phi, randomized
179         TH2F*                   fHistRhoVsRCPtRand[10];         //! rho * A vs rcpt, randomized
180         TH1F*                   fHistRCPtRand[10];              //! rcpt, randomized
181         TH2F*                   fHistDeltaPtDeltaPhiRand[10];   //! dpt vs dphi, randomized
182         // jet histograms (after kinematic cuts)
183         TH1F*                   fHistJetPtRaw[10];              //! jet pt - no background subtraction
184         TH1F*                   fHistJetPt[10];                 //! pt of found jets (background subtracted)
185         TH2F*                   fHistJetEtaPhi[10];             //! eta and phi correlation
186         TH2F*                   fHistJetPtArea[10];             //! jet pt versus area
187         TH2F*                   fHistJetPtConstituents[10];     //! jet pt versus number of constituents
188         // in plane, out of plane jet spectra
189         TH2F*                   fHistJetPsiTPCPt[10];            //! psi tpc versus pt
190         TH2F*                   fHistJetPsiVZEROAPt[10];         //! psi vzeroa versus pt
191         TH2F*                   fHistJetPsiVZEROCPt[10];         //! psi vzeroc versus pt
192         // phi minus psi 
193         TH1F*                   fHistDeltaPhiVZEROA[10]; //! phi minus psi_A
194         TH1F*                   fHistDeltaPhiVZEROC[10]; //! phi minus psi_C
195         TH1F*                   fHistDeltaPhiTPC[10];    //! phi minus psi_TPC
196
197         AliAnalysisTaskRhoVnModulation(const AliAnalysisTaskRhoVnModulation&);                  // not implemented
198         AliAnalysisTaskRhoVnModulation& operator=(const AliAnalysisTaskRhoVnModulation&);       // not implemented
199
200         ClassDef(AliAnalysisTaskRhoVnModulation, 2);
201 };
202
203 #endif