]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGJE/EMCALJetTasks/AliAnalysisTaskLocalRho.h
From Redmer: Task to be used to calculate rho relative to RP.
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliAnalysisTaskLocalRho.h
CommitLineData
5ce0dc6d 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
15class TF1;
16class THF1;
17class THF2;
18class TProfile;
19class AliLocalRhoParameter;
20
21class 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) {if (fFitModulation) delete fFitModulation;
70 fFitModulation = fit; }
71 void SetModulationFitMinMaxP(Float_t m, Float_t n) {fMinPvalue = m; fMaxPvalue = n; }
72 void SetModulationFitType(fitModulationType type) {fFitModulationType = type; }
73 void SetQCnRecoveryType(qcRecovery type) {fQCRecovery = type; }
74 void SetModulationFitOptions(TString opt) {fFitModulationOptions = opt; }
75 void SetReferenceDetector(detectorType type) {fDetectorType = type; }
76 void SetUsePtWeight(Bool_t w) {fUsePtWeight = w; }
77 void SetRunModeType(runModeType type) {fRunModeType = type; }
78 void SetForceAbsVnHarmonics(Bool_t f) {fAbsVnHarmonics = f; }
79 void SetExcludeLeadingJetsFromFit(Float_t n) {fExcludeLeadingJetsFromFit = n; }
80 void SetRebinSwapHistoOnTheFly(Bool_t r) {fRebinSwapHistoOnTheFly = r; }
81 void SetSaveThisPercentageOfFits(Float_t p) {fPercentageOfFits = p; }
82 void SetUseV0EventPlaneFromHeader(Bool_t h) {fUseV0EventPlaneFromHeader = h;}
83 void SetSoftTrackMinMaxPt(Float_t min, Float_t max) {fSoftTrackMinPt = min; fSoftTrackMaxPt = max;}
84 // getters
85 TString GetLocalRhoName() const {return fLocalRhoName; }
86 // numerical evaluations
87 void CalculateEventPlaneVZERO(Double_t vzero[2][2]) const;
88 void CalculateEventPlaneTPC(Double_t* tpc);
89 void CalculateEventPlaneCombinedVZERO(Double_t* comb) const;
90 Double_t CalculateQC2(Int_t harm);
91 Double_t CalculateQC4(Int_t harm);
92 // helper calculations for the q-cumulant analysis, also used by AliAnalyisTaskJetFlow
93 void QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ);
94 Double_t QCnS(Int_t i, Int_t j);
95 Double_t QCnM();
96 Double_t QCnM11();
97 Double_t QCnM1111();
98 Bool_t QCnRecovery(Double_t psi2, Double_t psi3);
99 // analysis details
100 Bool_t CorrectRho(Double_t psi2, Double_t psi3);
101 void FillEventPlaneHistograms(Double_t psi2, Double_t psi3) const;
102 void FillAnalysisSummaryHistogram() const;
103 // track selection
104 /* inline */ Bool_t PassesCuts(const AliVTrack* track) const {
105 if(!track) return kFALSE;
106 return (track->Pt() < fTrackPtCut || track->Eta() < fTrackMinEta || track->Eta() > fTrackMaxEta || track->Phi() < fTrackMinPhi || track->Phi() > fTrackMaxPhi) ? kFALSE : kTRUE; }
107 /* inline */ Bool_t PassesCuts(AliEmcalJet* jet) const {
108 if(!jet || fJetRadius <= 0) return kFALSE;
109 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; }
110 // filling histograms
111 virtual void Terminate(Option_t* option);
112
113 private:
114 Int_t fDebug; // debug level (0 none, 1 fcn calls, 2 verbose)
115 Bool_t fInitialized; //! is the analysis initialized?
116 Bool_t fAttachToEvent; // attach local rho to the event
117 Bool_t fFillHistograms; // fill qa histograms
118 Bool_t fNoEventWeightsForQC; // don't store event weights for qc analysis
119 TString fLocalRhoName; // name for local rho
120 TArrayI* fCentralityClasses; //-> centrality classes (maximum 10) used for QA
121 TH1F* fUserSuppliedV2; // histo with integrated v2
122 TH1F* fUserSuppliedV3; // histo with integrated v3
123 TH1F* fUserSuppliedR2; // correct the extracted v2 with this r
124 TH1F* fUserSuppliedR3; // correct the extracted v3 with this r
125 Int_t fNAcceptedTracks; //! number of accepted tracks
126 Int_t fNAcceptedTracksQCn; //! accepted tracks for QCn
127 Int_t fInCentralitySelection; //! centrality bin, only for QA plots
128 fitModulationType fFitModulationType; // fit modulation type
129 qcRecovery fQCRecovery; // recovery type for e-by-e qc method
130 Bool_t fUsePtWeight; // use dptdphi instead of dndphi
131 detectorType fDetectorType; // type of detector used for modulation fit
132 TString fFitModulationOptions; // fit options for modulation fit
133 runModeType fRunModeType; // run mode type
134 TF1* fFitModulation; //-> modulation fit for rho
135 Float_t fMinPvalue; // minimum value of p
136 Float_t fMaxPvalue; // maximum value of p
137 AliLocalRhoParameter* fLocalRho; //! local rho
138 // additional jet cuts (most are inherited)
139 Float_t fLocalJetMinEta; // local eta cut for jets
140 Float_t fLocalJetMaxEta; // local eta cut for jets
141 Float_t fLocalJetMinPhi; // local phi cut for jets
142 Float_t fLocalJetMaxPhi; // local phi cut for jets
143 Float_t fSoftTrackMinPt; // min pt for soft tracks
144 Float_t fSoftTrackMaxPt; // max pt for soft tracks
145 // general qa histograms
146 TH1F* fHistPvalueCDF; //! cdf value of chisquare p
147 // general settings
148 Bool_t fAbsVnHarmonics; // force postive local rho
149 Float_t fExcludeLeadingJetsFromFit; // exclude n leading jets from fit
150 Bool_t fRebinSwapHistoOnTheFly; // rebin swap histo on the fly
151 Float_t fPercentageOfFits; // save this percentage of fits
152 Bool_t fUseV0EventPlaneFromHeader; // use the vzero event plane from the header
153 // transient object pointers
154 TList* fOutputList; //! output list
155 TList* fOutputListGood; //! output list for local analysis
156 TList* fOutputListBad; //! output list for local analysis
157 TH1F* fHistSwap; //! swap histogram
158 TH1F* fHistAnalysisSummary; //! flags
159 TProfile* fProfV2; //! extracted v2
160 TProfile* fProfV2Cumulant; //! v2 cumulant
161 TProfile* fProfV3; //! extracted v3
162 TProfile* fProfV3Cumulant; //! v3 cumulant
163 TH1F* fHistPsi2[10]; //! psi 2
164 TH1F* fHistPsi3[10]; //! psi 3
165
166 AliAnalysisTaskLocalRho(const AliAnalysisTaskLocalRho&); // not implemented
167 AliAnalysisTaskLocalRho& operator=(const AliAnalysisTaskLocalRho&); // not implemented
168
169 ClassDef(AliAnalysisTaskLocalRho, 1);
170};
171
172#endif