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