]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliAnalysisTaskLocalRho.h
aa7c2efd734220fd23540ece1eae027156dcf89a
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliAnalysisTaskLocalRho.h
1 #ifndef ALIANALYSISTASKLOCALRHO_H
2 #define ALIANALYSISTASKLOCALRHO_H
3
4 // $Id$
5
6 #include <AliAnalysisTaskEmcalJetDev.h>
7 #include <AliEmcalJet.h>
8 #include <AliVEvent.h>
9 #include <AliVTrack.h>
10 #include <AliVCluster.h>
11 #include <TClonesArray.h>
12 #include <TMath.h>
13 #include <TRandom3.h>
14 #include <AliLog.h>
15 #include <AliJetContainer.h>
16
17 class TF1;
18 class THF1;
19 class THF2;
20 class TProfile;
21 class AliLocalRhoParameter;
22
23 class AliAnalysisTaskLocalRho : public AliAnalysisTaskEmcalJetDev {
24     public:
25          // enumerators
26         enum fitModulationType  { kNoFit, kV2, kV3, kCombined, kFourierSeries, kIntegratedFlow, kQC2, kQC4 }; // fit type
27         enum detectorType       { kTPC, kVZEROA, kVZEROC, kVZEROComb};  // detector that was used
28         enum qcRecovery         { kFixedRho, kNegativeVn, kTryFit };    // how to deal with negative cn value for qcn value
29         enum runModeType        { kLocal, kGrid };                      // run mode type
30         // constructors, destructor
31                                 AliAnalysisTaskLocalRho();
32                                 AliAnalysisTaskLocalRho(const char *name, runModeType type);
33         virtual                 ~AliAnalysisTaskLocalRho();
34         // setting up the task and technical aspects
35         void                    ExecOnce();
36         Bool_t                  InitializeAnalysis();
37         virtual void            UserCreateOutputObjects();
38         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);
39         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);
40         virtual Bool_t          Run();
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         // 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         // setters - setup how to run
60         void                    SetDebugMode(Int_t d)                           {fDebug = d;}
61         void                    SetCentralityClasses(TArrayI* c)                {fCentralityClasses = c;}
62         void                    SetAttachToEvent(Bool_t a)                      {fAttachToEvent = a;}
63         void                    SetUseScaledRho(Bool_t s)                       {fUseScaledRho = s;}
64         void                    SetFillHistograms(Bool_t b)                     {fFillHistograms = b;}
65         // setters - analysis details
66         void                    SetNoEventWeightsForQC(Bool_t e)                {fNoEventWeightsForQC = e;}
67         void                    SetIntegratedFlow(TH1F* i, TH1F* j)             {fUserSuppliedV2 = i;
68                                                                                  fUserSuppliedV3 = j; }
69         void                    SetOnTheFlyResCorrection(TH1F* r2, TH1F* r3)    {fUserSuppliedR2 = r2;
70                                                                                  fUserSuppliedR3 = r3; }
71         void                    SetModulationFit(TF1* fit);
72         void                    SetModulationFitMinMaxP(Float_t m, Float_t n)   {fMinPvalue = m; fMaxPvalue = n; }
73         void                    SetModulationFitType(fitModulationType type)    {fFitModulationType = type; }
74         void                    SetQCnRecoveryType(qcRecovery type)             {fQCRecovery = type; }
75         void                    SetModulationFitOptions(TString opt)            {fFitModulationOptions = opt; }
76         void                    SetReferenceDetector(detectorType type)         {fDetectorType = type; }
77         void                    SetUsePtWeight(Bool_t w)                        {fUsePtWeight = w; }
78         void                    SetRunModeType(runModeType type)                {fRunModeType = type; }
79         void                    SetForceAbsVnHarmonics(Bool_t f)                {fAbsVnHarmonics = f; }
80         void                    SetExcludeLeadingJetsFromFit(Float_t n)         {fExcludeLeadingJetsFromFit = n; }
81         void                    SetRebinSwapHistoOnTheFly(Bool_t r)             {fRebinSwapHistoOnTheFly = r; }
82         void                    SetSaveThisPercentageOfFits(Float_t p)          {fPercentageOfFits = p; }
83         void                    SetUseV0EventPlaneFromHeader(Bool_t h)          {fUseV0EventPlaneFromHeader = h;}
84         void                    SetSoftTrackMinMaxPt(Float_t min, Float_t max)  {fSoftTrackMinPt = min; fSoftTrackMaxPt = max;}
85         // getters
86         TString                 GetLocalRhoName() const                         {return fLocalRhoName; }
87         // numerical evaluations
88         void                    CalculateEventPlaneVZERO(Double_t vzero[2][2]) const;
89         void                    CalculateEventPlaneTPC(Double_t* tpc);
90         void                    CalculateEventPlaneCombinedVZERO(Double_t* comb) const;
91         Double_t                CalculateQC2(Int_t harm);
92         Double_t                CalculateQC4(Int_t harm);
93         // helper calculations for the q-cumulant analysis, also used by AliAnalyisTaskJetFlow
94         void                    QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ);
95         Double_t                QCnS(Int_t i, Int_t j);
96         Double_t                QCnM();
97         Double_t                QCnM11();
98         Double_t                QCnM1111();
99         Bool_t                  QCnRecovery(Double_t psi2, Double_t psi3);
100         // analysis details
101         Bool_t                  CorrectRho(Double_t psi2, Double_t psi3);
102         void                    FillEventPlaneHistograms(Double_t psi2, Double_t psi3) const;
103         void                    FillAnalysisSummaryHistogram() const;
104         // track selection
105         /* inline */    Bool_t PassesCuts(AliVTrack* track) const { return AcceptTrack(track, 0);}
106         /* inline */    Bool_t PassesCuts(AliEmcalJet* jet) { return AcceptJet(jet, 0);}
107         virtual void            Terminate(Option_t* option);
108
109     private: 
110         Int_t                   fDebug;                 // debug level (0 none, 1 fcn calls, 2 verbose)
111         Bool_t                  fInitialized;           //! is the analysis initialized?
112         Bool_t                  fAttachToEvent;         // attach local rho to the event
113         Bool_t                  fFillHistograms;        // fill qa histograms
114         Bool_t                  fNoEventWeightsForQC;   // don't store event weights for qc analysis
115         Bool_t                  fUseScaledRho;          // use scaled rho
116         TArrayI*                fCentralityClasses;     // centrality classes (maximum 10) used for QA
117         TH1F*                   fUserSuppliedV2;        // histo with integrated v2
118         TH1F*                   fUserSuppliedV3;        // histo with integrated v3
119         TH1F*                   fUserSuppliedR2;        // correct the extracted v2 with this r
120         TH1F*                   fUserSuppliedR3;        // correct the extracted v3 with this r
121         Int_t                   fNAcceptedTracks;       //! number of accepted tracks
122         Int_t                   fNAcceptedTracksQCn;    //! accepted tracks for QCn
123         Int_t                   fInCentralitySelection; //! centrality bin, only for QA plots
124         fitModulationType       fFitModulationType;     // fit modulation type
125         qcRecovery              fQCRecovery;            // recovery type for e-by-e qc method
126         Bool_t                  fUsePtWeight;           // use dptdphi instead of dndphi
127         detectorType            fDetectorType;          // type of detector used for modulation fit
128         TString                 fFitModulationOptions;  // fit options for modulation fit
129         runModeType             fRunModeType;           // run mode type 
130         TF1*                    fFitModulation;         // modulation fit for rho
131         Float_t                 fMinPvalue;             // minimum value of p
132         Float_t                 fMaxPvalue;             // maximum value of p
133         // additional jet cuts (most are inherited)
134         Float_t                 fLocalJetMinEta;        // local eta cut for jets
135         Float_t                 fLocalJetMaxEta;        // local eta cut for jets
136         Float_t                 fLocalJetMinPhi;        // local phi cut for jets
137         Float_t                 fLocalJetMaxPhi;        // local phi cut for jets
138         Float_t                 fSoftTrackMinPt;        // min pt for soft tracks
139         Float_t                 fSoftTrackMaxPt;        // max pt for soft tracks
140         // general qa histograms
141         TH1F*                   fHistPvalueCDF;         //! cdf value of chisquare p
142         // general settings
143         Bool_t                  fAbsVnHarmonics;        // force postive local rho
144         Float_t                 fExcludeLeadingJetsFromFit;    // exclude n leading jets from fit
145         Bool_t                  fRebinSwapHistoOnTheFly;       // rebin swap histo on the fly
146         Float_t                 fPercentageOfFits;      // save this percentage of fits
147         Bool_t                  fUseV0EventPlaneFromHeader;    // use the vzero event plane from the header
148         // transient object pointers
149         TList*                  fOutputList;            //! output list
150         TList*                  fOutputListGood;        //! output list for local analysis
151         TList*                  fOutputListBad;         //! output list for local analysis
152         TH1F*                   fHistSwap;              //! swap histogram
153         TH1F*                   fHistAnalysisSummary;   //! flags
154         TProfile*               fProfV2;                //! extracted v2
155         TProfile*               fProfV2Cumulant;        //! v2 cumulant
156         TProfile*               fProfV3;                //! extracted v3
157         TProfile*               fProfV3Cumulant;        //! v3 cumulant
158         TH1F*                   fHistPsi2[10];          //! psi 2
159         TH1F*                   fHistPsi3[10];          //! psi 3
160
161         AliAnalysisTaskLocalRho(const AliAnalysisTaskLocalRho&);                  // not implemented
162         AliAnalysisTaskLocalRho& operator=(const AliAnalysisTaskLocalRho&);       // not implemented
163
164         ClassDef(AliAnalysisTaskLocalRho, 4);
165 };
166
167 #endif