]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskRhoVnModulation.h
small update to jet cuts (in eta) to accommodate running on semi-good
[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 #include <TRandom3.h>
16
17 class TF1;
18 class THF1;
19 class THF2;
20 class TProfile;
21
22 class AliAnalysisTaskRhoVnModulation : public AliAnalysisTaskEmcalJet
23 {
24     public:
25          // enumerators
26         enum fitModulationType  { kNoFit, kV2, kV3, kCombined, kUser, kFourierSeries, kIntegratedFlow }; // 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, Bool_t append = kTRUE);
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, Bool_t append = kTRUE);
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 PhaseShift(Double_t x, Double_t n) const {
46             x = PhaseShift(x);
47             if(TMath::Nint(n)==2) while (x>TMath::Pi()) x-=TMath::Pi();
48             if(TMath::Nint(n)==3) {
49                 if(x>2.*TMath::TwoPi()/n) x = TMath::TwoPi() - x;
50                 if(x>TMath::TwoPi()/n) x = TMath::TwoPi()-(x+TMath::TwoPi()/n);
51             }
52             return x; }
53         /* inline */    Double_t ChiSquarePDF(Int_t ndf, Double_t x) const {
54             Double_t n(ndf/2.), denom(TMath::Power(2, n)*TMath::Gamma(n));
55             if (denom!=0)  return ((1./denom)*TMath::Power(x, n-1)*TMath::Exp(-x/2.)); 
56             return -999; }
57         // note that the cdf of the chisquare distribution is the normalized lower incomplete gamma function
58         /* inline */    Double_t ChiSquareCDF(Int_t ndf, Double_t x) const { return TMath::Gamma(ndf/2., x/2.); }
59         /* inline */    Double_t RhoVal() const { return (fRho) ? fRho->GetVal(): -999.;}                 
60         /* inline */    Double_t RhoVal(Double_t phi, Double_t r, Double_t n) const {
61             if(!fFitModulation) return RhoVal(); // coverity
62             switch (fFitModulationType) {
63                 case kNoFit : return RhoVal();
64                 default : {
65                     Double_t denom(2*r*fFitModulation->GetParameter(0));
66                     return  (denom <= 0.) ? RhoVal() : n*(fFitModulation->Integral(phi-r, phi+r)/denom); 
67                 }
68             }
69         }
70         // setters - analysis setup
71         void                    SetDebugMode(Int_t d)                           {fDebug = d;}
72         void                    SetFillQAHistograms(Bool_t qa)                  {fFillQAHistograms = qa;}
73         void                    SetCentralityClasses(TArrayI* c)                {fCentralityClasses = c;}
74         void                    SetIntegratedFlow(TH1F* i, TH1F* j)             {fUserSuppliedV2 = i;
75                                                                                  fUserSuppliedV3 = j; }
76         void                    SetNameJetClones(const char* name)              {fNameJetClones = name; }
77         void                    SetNamePicoTrackClones(const char* name)        {fNamePicoTrackClones = name; }
78         void                    SetNameRho(const char* name)                    {fNameRho = name; }
79         void                    SetRandomSeed(TRandom3* r)                      {if (fRandom) delete fRandom; fRandom = r; }
80         void                    SetModulationFit(TF1* fit)                      {if (fFitModulation) delete fFitModulation;
81                                                                                  fFitModulation = fit; }
82         void                    SetModulationFitMinMaxP(Float_t m, Float_t n)   {fMinPvalue = m; fMaxPvalue = n; }
83         void                    SetModulationFitType(fitModulationType type)    {fFitModulationType = type; }
84         void                    SetModulationFitOptions(TString opt)            {fFitModulationOptions = opt; }
85         void                    SetReferenceDetector(detectorType type)         {fDetectorType = type; }
86         void                    SetUsePtWeight(Bool_t w)                        {fUsePtWeight = w; }
87         void                    SetRunModeType(runModeType type)                {fRunModeType = type; }
88         void                    SetAbsVertexZ(Float_t v)                        {fAbsVertexZ = v; }
89         void                    SetMinDistanceRctoLJ(Float_t m)                 {fMinDisanceRCtoLJ = m; }
90         void                    SetRandomConeRadius(Float_t r)                  {fRandomConeRadius = r; }
91         void                    SetMinLeadingHadronPt(Double_t m)               {fMinLeadingHadronPt = m; }
92         void                    SetForceAbsVnHarmonics(Bool_t f)                {fAbsVnHarmonics = f; }
93         void                    SetExcludeLeadingJetsFromFit(Float_t n)         {fExcludeLeadingJetsFromFit = n; }
94         void                    SetRebinSwapHistoOnTheFly(Bool_t r)             {fRebinSwapHistoOnTheFly = r; }
95         void                    SetSaveThisPercentageOfFits(Float_t p)          {fPercentageOfFits = p; }
96         void                    SetUseV0EventPlaneFromHeader(Bool_t h)          {fUseV0EventPlaneFromHeader = h;}
97         void                    SetSetPtSub(Bool_t s)                           {fSetPtSub = s; }
98         void                    SetExplicitOutlierCutForYear(Int_t y)           {fExplicitOutlierCut = y;}
99         // local cuts
100         void                    SetLocalJetMinMaxEta(Float_t min, Float_t max)  {fLocalJetMinEta = min; fLocalJetMaxEta = max;}
101         void                    SetLocalJetMinMaxEta(Float_t R)                 {fLocalJetMinEta = - 0.9 + R; fLocalJetMaxEta = 0.9 - R; }
102         void                    SetLocalJetMinMaxPhi(Float_t min, Float_t max)  {fLocalJetMinPhi = min; fLocalJetMaxEta = max;}
103         // 'trivial' helper calculations
104         void                    CalculateEventPlaneVZERO(Double_t vzero[2][2]) const;
105         void                    CalculateEventPlaneTPC(Double_t* tpc);
106         void                    CalculateEventPlaneResolution(Double_t vzero[2][2], Double_t* tpc) const;
107         void                    CalculateRandomCone(Float_t &pt, Float_t &eta, Float_t &phi, AliEmcalJet* jet = 0x0, Bool_t randomize = 0) const;
108         // analysis details
109         Bool_t                  CorrectRho(Double_t psi2, Double_t psi3);
110         // event and track selection
111         /* inline */    Bool_t PassesCuts(const AliVTrack* track) const {
112             if(!track) return kFALSE;
113             return (track->Pt() < fTrackPtCut || track->Eta() < fTrackMinEta || track->Eta() > fTrackMaxEta || track->Phi() < fTrackMinPhi || track->Phi() > fTrackMaxPhi) ? kFALSE : kTRUE; }
114         /* inline */    Bool_t PassesCuts(AliEmcalJet* jet) const {
115             if(!jet || fJetRadius <= 0) return kFALSE;
116             return (GetLeadingHadronPt(jet) < fMinLeadingHadronPt || jet->Pt() < fJetPtCut || jet->Area()/(fJetRadius*fJetRadius*TMath::Pi()) < fPercAreaCut || jet->Eta() < fJetMinEta || jet->Eta() > fJetMaxEta || jet->Phi() < fJetMinPhi || jet->Phi() > fJetMaxPhi) ? kFALSE : kTRUE; }
117         Bool_t                  PassesCuts(AliVEvent* event);
118         Bool_t                  PassesCuts(Int_t year);
119         Bool_t                  PassesCuts(const AliVCluster* track) const;
120         // filling histograms
121         void                    FillHistogramsAfterSubtraction(Double_t vzero[2][2], Double_t* tpc) const;
122         void                    FillTrackHistograms() const;
123         void                    FillClusterHistograms() const;
124         void                    FillCorrectedClusterHistograms() const;
125         void                    FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* tpc) const;
126         void                    FillRhoHistograms() const;
127         void                    FillDeltaPtHistograms(Double_t psi2, Double_t psi3) const; 
128         void                    FillJetHistograms(Double_t vzero[2][2], Double_t* psi) const;
129         void                    FillDeltaPhiHistograms(Double_t vzero[2][2], Double_t* tpc) const;
130         void                    FillQAHistograms(AliVTrack* vtrack) const;
131         void                    FillQAHistograms(AliVEvent* vevent);
132         void                    FillAnalysisSummaryHistogram() const;
133         virtual void            Terminate(Option_t* option);
134         // interface methods for the output file
135         void                    SetOutputList(TList* l) {fOutputList = l;}
136         TH1F*                   GetResolutionFromOuptutFile(detectorType detector, Int_t h = 2, TArrayD* c = 0x0);
137         TH1F*                   CorrectForResolutionDiff(TH1F* v, detectorType detector, TArrayD* cen, Int_t c, Int_t h = 2);
138         TH1F*                   CorrectForResolutionInt(TH1F* v, detectorType detector, TArrayD* cen, Int_t h = 2);
139     private:
140         // analysis flags and settings
141         Int_t                   fDebug;                 // debug level (0 none, 1 fcn calls, 2 verbose)
142         Bool_t                  fInitialized;           //! is the analysis initialized?
143         Bool_t                  fFillQAHistograms;      // fill qa histograms
144         TArrayI*                fCentralityClasses;     //-> centrality classes (maximum 10)
145         TH1F*                   fUserSuppliedV2;        // histo with integrated v2
146         TH1F*                   fUserSuppliedV3;        // histo with integrated v3
147         // members
148         Int_t                   fNAcceptedTracks;       //! number of accepted tracks
149         fitModulationType       fFitModulationType;     // fit modulation type
150         Bool_t                  fUsePtWeight;           // use dptdphi instead of dndphi
151         detectorType            fDetectorType;          // type of detector used for modulation fit
152         TString                 fFitModulationOptions;  // fit options for modulation fit
153         runModeType             fRunModeType;           // run mode type 
154         dataType                fDataType;              // datatype 
155         TRandom3*               fRandom;                //-> dont use gRandom to not interfere with other tasks
156         Int_t                   fMappedRunNumber;       //! mapped runnumer (for QA)
157         Int_t                   fInCentralitySelection; //! centrality bin
158         TF1*                    fFitModulation;         //-> modulation fit for rho
159         Float_t                 fMinPvalue;             // minimum value of p
160         Float_t                 fMaxPvalue;             // maximum value of p
161         const char*             fNameJetClones;         //! collection of tclones array with jets
162         const char*             fNamePicoTrackClones;   //! collection of tclones with pico tracks
163         const char*             fNameRho;               //! name of rho
164         // additional jet cuts (most are inherited)
165         Float_t                 fLocalJetMinEta;        // local eta cut for jets
166         Float_t                 fLocalJetMaxEta;        // local eta cut for jets
167         Float_t                 fLocalJetMinPhi;        // local phi cut for jets
168         Float_t                 fLocalJetMaxPhi;        // local phi cut for jets
169         // event cuts
170         Float_t                 fAbsVertexZ;            // cut on zvertex
171         // general qa histograms
172         TH1F*                   fHistCentrality;        //! accepted centrality
173         TH1F*                   fHistVertexz;           //! accepted verte
174         TH2F*                   fHistRunnumbersPhi;     //! run numbers averaged phi
175         TH2F*                   fHistRunnumbersEta;     //! run numbers averaged eta
176         TH1F*                   fHistPvaluePDF;         //! pdf value of chisquare p
177         TH1F*                   fHistPvalueCDF;         //! cdf value of chisquare p
178         // general settings
179         Float_t                 fMinDisanceRCtoLJ;      // min distance between rc and leading jet
180         Float_t                 fRandomConeRadius;      // radius of random cone
181         Bool_t                  fAbsVnHarmonics;        // force postive local rho
182         Float_t                 fExcludeLeadingJetsFromFit;    // exclude n leading jets from fit
183         Bool_t                  fRebinSwapHistoOnTheFly;       // rebin swap histo on the fly
184         Float_t                 fPercentageOfFits;      // save this percentage of fits
185         Bool_t                  fUseV0EventPlaneFromHeader;    // use the vzero event plane from the header
186         Bool_t                  fSetPtSub;              // store the subtracted pt in the jet
187         Int_t                   fExplicitOutlierCut;    // cut on correlation of tpc and global multiplicity
188         Double_t                fMinLeadingHadronPt;    // minimum pt for leading hadron
189         // transient object pointers
190         TList*                  fOutputList;            //! output list
191         TList*                  fOutputListGood;        //! output list for local analysis
192         TList*                  fOutputListBad;         //! output list for local analysis
193         TH1F*                   fHistAnalysisSummary;   //! analysis summary
194         TH1F*                   fHistSwap;              //! swap histogram
195         TProfile*               fProfV2;                //! extracted v2
196         TProfile*               fProfV2Resolution[10];  //! resolution parameters for v2
197         TProfile*               fProfV3;                //! extracted v3
198         TProfile*               fProfV3Resolution[10];  //! resolution parameters for v3
199         // qa histograms for accepted pico tracks
200         TH1F*                   fHistPicoTrackPt[10];    //! pt of all charged tracks
201         TH2F*                   fHistPicoCat1[10];       //! pico tracks spd hit and refit
202         TH2F*                   fHistPicoCat2[10];       //! pico tracks wo spd hit w refit, constrained
203         TH2F*                   fHistPicoCat3[10];       //! pico tracks wo spd hit wo refit, constrained
204         // qa histograms for accepted emcal clusters
205         /* TH1F*                   fHistClusterPt[10];      //! pt uncorrected emcal clusters */
206         /* TH1F*                   fHistClusterPhi[10];     //! phi uncorrected emcal clusters */
207         /* TH1F*                   fHistClusterEta[10];     //! eta uncorrected emcal clusters */
208         // qa histograms for accepted emcal clusters aftehadronic correction
209         /* TH1F*                   fHistClusterCorrPt[10];  //! pt corrected emcal clusters */
210         /* TH1F*                   fHistClusterCorrPhi[10]; //! phi corrected emcal clusters */
211         /* TH1F*                   fHistClusterCorrEta[10]; //! eta corrected emcal clusters */
212         // qa event planes
213         TProfile*               fHistPsiControl;         //! event plane control histogram
214         TProfile*               fHistPsiSpread;          //! event plane spread histogram
215         TH1F*                   fHistPsiVZEROA;          //! psi 2 from vzero a
216         TH1F*                   fHistPsiVZEROC;          //! psi 2 from vzero c
217         TH1F*                   fHistPsiTPC;             //! psi 2 from tpc
218         // background
219         TH1F*                   fHistRhoPackage[10];     //! rho as estimated by emcal jet package
220         TH1F*                   fHistRho[10];            //! background
221         TH2F*                   fHistRhoVsMult;          //! rho versus multiplicity
222         TH2F*                   fHistRhoVsCent;          //! rho veruss centrality
223         TH2F*                   fHistRhoAVsMult;         //! rho * A vs multiplicity for all jets
224         TH2F*                   fHistRhoAVsCent;         //! rho * A vs centrality for all jets
225         // delta pt distributions
226         TH2F*                   fHistRCPhiEta[10];              //! random cone eta and phi
227         TH2F*                   fHistRhoVsRCPt[10];             //! rho * A vs rcpt
228         TH1F*                   fHistRCPt[10];                  //! rcpt
229         TH2F*                   fHistDeltaPtDeltaPhi2[10];      //! dpt vs dphi
230         TH2F*                   fHistDeltaPtDeltaPhi3[10];      //! dpt vs dphi
231         TH2F*                   fHistRCPhiEtaExLJ[10];          //! random cone eta and phi, excl leading jet
232         TH2F*                   fHistRhoVsRCPtExLJ[10];         //! rho * A vs rcpt, excl leading jet
233         TH1F*                   fHistRCPtExLJ[10];              //! rcpt, excl leading jet
234         TH2F*                   fHistDeltaPtDeltaPhi2ExLJ[10];  //! dpt vs dphi, excl leading jet
235         TH2F*                   fHistDeltaPtDeltaPhi3ExLJ[10];  //! dpt vs dphi, excl leading jet
236         /* TH2F*                   fHistRCPhiEtaRand[10];          //! random cone eta and phi, randomized */
237         /* TH2F*                   fHistRhoVsRCPtRand[10];         //! rho * A vs rcpt, randomized */
238         /* TH1F*                   fHistRCPtRand[10];              //! rcpt, randomized */ 
239         /* TH2F*                   fHistDeltaPtDeltaPhi2Rand[10];  //! dpt vs dphi, randomized */
240         /* TH2F*                   fHistDeltaPtDeltaPhi3Rand[10];  //! dpt vs dphi, randomized */
241         // jet histograms (after kinematic cuts)
242         TH1F*                   fHistJetPtRaw[10];              //! jet pt - no background subtraction
243         TH1F*                   fHistJetPt[10];                 //! pt of found jets (background subtracted)
244         TH2F*                   fHistJetEtaPhi[10];             //! eta and phi correlation
245         TH2F*                   fHistJetPtArea[10];             //! jet pt versus area
246         TH2F*                   fHistJetPtConstituents[10];     //! jet pt versus number of constituents
247         TH2F*                   fHistJetEtaRho[10];             //! jet eta versus jet rho
248         // in plane, out of plane jet spectra
249         TH2F*                   fHistJetPsiTPCPt[10];            //! psi tpc versus pt
250         TH2F*                   fHistJetPsiVZEROAPt[10];         //! psi vzeroa versus pt
251         TH2F*                   fHistJetPsiVZEROCPt[10];         //! psi vzeroc versus pt
252         // phi minus psi 
253         TH1F*                   fHistDeltaPhi2VZEROA[10];       //! phi minus psi_A
254         TH1F*                   fHistDeltaPhi2VZEROC[10];       //! phi minus psi_C
255         TH1F*                   fHistDeltaPhi2TPC[10];          //! phi minus psi_TPC
256         TH1F*                   fHistDeltaPhi3VZEROA[10];       //! phi minus psi_A
257         TH1F*                   fHistDeltaPhi3VZEROC[10];       //! phi minus psi_C
258         TH1F*                   fHistDeltaPhi3TPC[10];          //! phi minus psi_TPC
259
260         AliAnalysisTaskRhoVnModulation(const AliAnalysisTaskRhoVnModulation&);                  // not implemented
261         AliAnalysisTaskRhoVnModulation& operator=(const AliAnalysisTaskRhoVnModulation&);       // not implemented
262
263         ClassDef(AliAnalysisTaskRhoVnModulation, 10);
264 };
265
266 #endif