]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/EBYE/BalanceFunctions/AliAnalysisTaskBF.h
adding optional checks on 1st event in chunk and pileup events (MW)
[u/mrichter/AliRoot.git] / PWGCF / EBYE / BalanceFunctions / AliAnalysisTaskBF.h
1 #ifndef ALIANALYSISTASKBF_H
2 #define ALIANALYSISTASKBF_H
3
4 // Analysis task for the BF code
5 // Authors: Panos Cristakoglou@cern.ch
6
7 class TList;
8 class TH1F;
9 class TH2F;
10 class TH2D;
11 class TF1;
12
13 class AliBalance;
14 class AliESDtrackCuts;
15
16 #include "AliAnalysisTaskSE.h"
17 #include "AliBalance.h"
18
19 #include "AliPID.h"  
20 #include "AliPIDResponse.h"
21 #include "AliPIDCombined.h"
22  
23
24 class AliAnalysisTaskBF : public AliAnalysisTaskSE {
25  public:
26   AliAnalysisTaskBF(const char *name = "AliAnalysisTaskBF");
27   virtual ~AliAnalysisTaskBF(); 
28   
29   
30   virtual void   UserCreateOutputObjects();
31   virtual void   UserExec(Option_t *option);
32   virtual void   FinishTaskOutput();
33   virtual void   Terminate(Option_t *);
34
35   void SetAnalysisObject(AliBalance *const analysis) {
36     fBalance         = analysis;
37     }
38   void SetShufflingObject(AliBalance *const analysisShuffled) {
39     fRunShuffling = kTRUE;
40     fShuffledBalance = analysisShuffled;
41   }
42   void SetAnalysisCutObject(AliESDtrackCuts *const trackCuts) {
43     fESDtrackCuts = trackCuts;}
44   void SetVertexDiamond(Double_t vx, Double_t vy, Double_t vz) {
45     fVxMax = vx;
46     fVyMax = vy;
47     fVzMax = vz;
48   }
49
50   //==============AOD analysis==============//
51   void SetAODtrackCutBit(Int_t bit){
52     fAODtrackCutBit = bit;
53   }
54
55   void SetKinematicsCutsAOD(Double_t ptmin, Double_t ptmax, Double_t etamin, Double_t etamax){
56     fPtMin  = ptmin;    fPtMax  = ptmax;
57     fEtaMin = etamin;   fEtaMax = etamax;
58   }
59
60   void SetExtraDCACutsAOD(Double_t DCAxy, Double_t DCAz){
61     fDCAxyCut  = DCAxy;
62     fDCAzCut = DCAz;
63   }
64
65    void SetExtraTPCCutsAOD(Double_t maxTPCchi2, Int_t minNClustersTPC){
66     fTPCchi2Cut      = maxTPCchi2;
67     fNClustersTPCCut = minNClustersTPC;
68   }
69
70   //==============MC analysis==============//
71    void SetKinematicsCutsMC(Double_t ptmin, Double_t ptmax,Double_t etamin, Double_t etamax){
72     fPtMin  = ptmin;   fPtMax  = ptmax;
73     fEtaMin = etamin;  fEtaMax = etamax;
74   }
75   void UseFlowAfterBurner(TF1 *gDifferentialV2) {
76     fDifferentialV2 = gDifferentialV2;
77     fUseFlowAfterBurner = kTRUE;
78   }
79   void ExcludeResonancesInMC() {fExcludeResonancesInMC = kTRUE;}
80
81   void SetPDGCode(Int_t gPdgCode) {
82     fUseMCPdgCode = kTRUE;
83     fPDGCodeToBeAnalyzed = gPdgCode;
84   }
85
86   //Centrality
87   void SetCentralityEstimator(const char* centralityEstimator) {fCentralityEstimator = centralityEstimator;}
88   const char* GetCentralityEstimator (void) const {return fCentralityEstimator;}
89   void SetCentralityPercentileRange(Double_t min, Double_t max) { 
90     fUseCentrality = kTRUE;
91     fCentralityPercentileMin=min;
92     fCentralityPercentileMax=max;
93   }
94   void SetImpactParameterRange(Double_t min, Double_t max) { 
95     fUseCentrality = kTRUE;
96     fImpactParameterMin=min;
97     fImpactParameterMax=max;
98   }
99
100   //multiplicity
101   void SetMultiplicityRange(Int_t min, Int_t max) {
102     fUseMultiplicity = kTRUE;
103     fNumberOfAcceptedTracksMin = min;
104     fNumberOfAcceptedTracksMax = max;}
105   
106   void UseOfflineTrigger() {fUseOfflineTrigger = kTRUE;}
107   
108   //Acceptance filter
109   void SetAcceptanceParameterization(TF1 *parameterization) {
110     fAcceptanceParameterization = parameterization;}
111
112   //pid
113   enum kDetectorUsedForPID { kTPCpid, kTOFpid, kTPCTOF }; // default TPC & TOF pid (via GetTPCpid & GetTOFpid)  
114   enum kParticleOfInterest { kElectron, kMuon, kPion, kKaon, kProton, kPhoton, kPi0, kNeutron, kKaon0, kEleCon, kDeuteron, kTriton, kHe3, kAlpha, kUnknown};  
115
116   void SetUseBayesianPID(Double_t gMinProbabilityValue) {
117     fUsePID = kTRUE; fUsePIDnSigma = kFALSE; fUsePIDPropabilities = kTRUE;
118     fMinAcceptedPIDProbability = gMinProbabilityValue; }
119
120   void SetUseNSigmaPID(Double_t gMaxNSigma) {
121     fUsePID = kTRUE; fUsePIDPropabilities = kFALSE; fUsePIDnSigma = kTRUE;
122     fPIDNSigma = gMaxNSigma; }
123
124   void SetParticleOfInterest(kParticleOfInterest poi) {
125     fParticleOfInterest = poi;}
126   void SetDetectorUsedForPID(kDetectorUsedForPID detConfig) {
127     fPidDetectorConfig = detConfig;}
128
129  private:
130   AliBalance *fBalance; //BF object
131   Bool_t fRunShuffling;//run shuffling or not
132   AliBalance *fShuffledBalance; //BF object (shuffled)
133   TList *fList; //fList object
134   TList *fListBF; //fList object
135   TList *fListBFS; //fList object
136   TList *fHistListPIDQA;  //! list of histograms
137
138   TH2D *fHistEventStats; //event stats
139   TH2F *fHistCentStats; //centrality stats
140   TH1F *fHistTriggerStats; //trigger stats
141   TH1F *fHistTrackStats; //Track filter bit stats
142   TH1F *fHistVx; //x coordinate of the primary vertex
143   TH1F *fHistVy; //y coordinate of the primary vertex
144   TH1F *fHistVz; //z coordinate of the primary vertex
145
146   TH2F *fHistClus;//number of clusters (QA histogram)
147   TH2F *fHistDCA;//DCA  (QA histogram)
148   TH1F *fHistChi2;//track chi2 (QA histogram)
149   TH1F *fHistPt;//transverse momentum (QA histogram)
150   TH1F *fHistEta;//pseudorapidity (QA histogram)
151   TH1F *fHistRapidity;//rapidity (QA histogram)
152   TH1F *fHistPhi;//phi (QA histogram)
153   TH1F *fHistPhiBefore;//phi before v2 afterburner (QA histogram)
154   TH1F *fHistPhiAfter;//phi after v2 afterburner (QA histogram)
155   TH1F *fHistPhiPos;//phi for positive particles (QA histogram)
156   TH1F *fHistPhiNeg;//phi for negative particles (QA histogram)
157   TH2F *fHistV0M;//V0 multiplicities (QA histogram)
158   TH2F *fHistRefTracks;//reference track multiplicities (QA histogram)
159
160   //============PID============//
161   TH2D *fHistdEdxVsPTPCbeforePID;//TPC dEdx vs momentum before PID cuts (QA histogram)
162   TH2D *fHistBetavsPTOFbeforePID;//beta vs momentum before PID cuts (QA histogram)
163   TH2D *fHistProbTPCvsPtbeforePID; //TPC probability vs pT before PID cuts (QA histogram)
164   TH2D *fHistProbTOFvsPtbeforePID;//TOF probability vs pT before PID cuts (QA histogram)
165   TH2D *fHistProbTPCTOFvsPtbeforePID;//TOF/TPC probability vs pT before PID cuts (QA histogram)
166   TH2D *fHistNSigmaTPCvsPtbeforePID;//TPC nsigma vs pT before PID cuts (QA histogram)
167   TH2D *fHistNSigmaTOFvsPtbeforePID;//TOF nsigma vs pT before PID cuts (QA histogram)
168   TH2D *fHistdEdxVsPTPCafterPID;//TPC dEdx vs momentum after PID cuts (QA histogram)
169   TH2D *fHistBetavsPTOFafterPID;//beta vs momentum after PID cuts (QA histogram)
170   TH2D *fHistProbTPCvsPtafterPID; //TPC probability vs pT after PID cuts (QA histogram)
171   TH2D *fHistProbTOFvsPtafterPID;//TOF probability vs pT after PID cuts (QA histogram)
172   TH2D *fHistProbTPCTOFvsPtafterPID;//TOF/TPC probability vs pT after PID cuts (QA histogram)
173   TH2D *fHistNSigmaTPCvsPtafterPID;//TPC nsigma vs pT after PID cuts (QA histogram)
174   TH2D *fHistNSigmaTOFvsPtafterPID;//TOF nsigma vs pT after PID cuts (QA histogram)
175
176
177   AliPIDResponse *fPIDResponse;     //! PID response object
178   AliPIDCombined       *fPIDCombined;     //! combined PID object
179   
180   kParticleOfInterest  fParticleOfInterest;//analyzed particle
181   kDetectorUsedForPID   fPidDetectorConfig;//used detector for PID
182
183   Bool_t fUsePID; //flag to use PID 
184   Bool_t fUsePIDnSigma;//flag to use nsigma method for PID
185   Bool_t fUsePIDPropabilities;//flag to use probability method for PID
186   Double_t fPIDNSigma;//nsigma cut for PID
187   Double_t fMinAcceptedPIDProbability;//probability cut for PID
188   //============PID============//
189
190   AliESDtrackCuts *fESDtrackCuts; //ESD track cuts
191
192   TString fCentralityEstimator;      //"V0M","TRK","TKL","ZDC","FMD"
193   Bool_t fUseCentrality;//use the centrality (PbPb) or not (pp)
194   Double_t fCentralityPercentileMin;//centrality percentile min
195   Double_t fCentralityPercentileMax;//centrality percentile max
196   Double_t fImpactParameterMin;//impact parameter min (used for MC)
197   Double_t fImpactParameterMax;//impact parameter max (used for MC)
198
199   Bool_t fUseMultiplicity;//use the multiplicity cuts
200   Int_t fNumberOfAcceptedTracksMin;//min. number of number of accepted tracks (used for the multiplicity dependence study - pp)
201   Int_t fNumberOfAcceptedTracksMax;//max. number of number of accepted tracks (used for the multiplicity dependence study - pp)
202  
203   TH2D *fHistNumberOfAcceptedTracks;//hisot to store the number of accepted tracks
204
205   Bool_t fUseOfflineTrigger;//Usage of the offline trigger selection
206
207   Double_t fVxMax;//vxmax
208   Double_t fVyMax;//vymax
209   Double_t fVzMax;//vzmax
210
211   Int_t fAODtrackCutBit;//track cut bit from track selection (only used for AODs)
212
213   Double_t fPtMin;//only used for AODs
214   Double_t fPtMax;//only used for AODs
215   Double_t fEtaMin;//only used for AODs
216   Double_t fEtaMax;//only used for AODs
217
218   Double_t fDCAxyCut;//only used for AODs
219   Double_t fDCAzCut;//only used for AODs
220
221   Double_t fTPCchi2Cut;//only used for AODs
222   Int_t fNClustersTPCCut;//only used for AODs
223
224   TF1 *fAcceptanceParameterization;//acceptance filter used for MC
225
226   TF1 *fDifferentialV2;//pt-differential v2 (from real data)
227   Bool_t fUseFlowAfterBurner;//Usage of a flow after burner
228
229   Bool_t fExcludeResonancesInMC;//flag to exclude the resonances' decay products from the MC analysis
230   Bool_t fUseMCPdgCode; //Boolean to analyze a set of particles in MC
231   Int_t fPDGCodeToBeAnalyzed; //Analyze a set of particles in MC
232
233   
234
235   AliAnalysisTaskBF(const AliAnalysisTaskBF&); // not implemented
236   AliAnalysisTaskBF& operator=(const AliAnalysisTaskBF&); // not implemented
237   
238   ClassDef(AliAnalysisTaskBF, 5); // example of analysis
239 };
240
241 #endif