]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/EBYE/BalanceFunctions/AliAnalysisTaskBFPsi.h
Corrected end-of-line behavior
[u/mrichter/AliRoot.git] / PWGCF / EBYE / BalanceFunctions / AliAnalysisTaskBFPsi.h
1 #ifndef ALIANALYSISTASKBFPSI_H
2 #define ALIANALYSISTASKBFPSI_H
3
4 // Analysis task for the BF vs Psi code
5 // Authors: Panos Cristakoglou@cern.ch
6
7 class TList;
8 class TH1F;
9 class TH2F;
10 class TH3F; 
11 class TF1;
12 class TH3D;
13
14 class AliBalancePsi;
15 class AliESDtrackCuts;
16 class AliEventPoolManager;
17
18
19 #include "AliAnalysisTaskSE.h"
20 #include "AliBalancePsi.h"
21
22 #include "AliPID.h"  
23 #include "AliPIDResponse.h"
24 #include "AliPIDCombined.h"
25  
26 //================================correction
27 #define kCENTRALITY 101  
28 //const Double_t centralityArrayForPbPb[kCENTRALITY+1] = {0.,5.,10.,20.,30.,40.,50.,60.,70.,80.};
29 //const TString centralityArrayForPbPb_string[kCENTRALITY] = {"0-5","5-10","10-20","20-30","30-40","40-50","50-60","60-70","70-80"};
30 //================================correction
31
32 class AliAnalysisTaskBFPsi : public AliAnalysisTaskSE {
33  public:
34   AliAnalysisTaskBFPsi(const char *name = "AliAnalysisTaskBFPsi");
35   virtual ~AliAnalysisTaskBFPsi(); 
36    
37   virtual void   UserCreateOutputObjects();
38   virtual void   UserExec(Option_t *option);
39   virtual void   FinishTaskOutput();
40   virtual void   Terminate(Option_t *);
41
42   //========================correction
43   virtual void   SetInputCorrection(TString filename, 
44                                     Int_t nCentralityBins, 
45                                     Double_t *centralityArrayForCorrections);
46   //========================correction
47   // void SetDebugLevel() {fDebugLevel = kTRUE;} //hides overloaded virtual function
48
49   void SetAnalysisObject(AliBalancePsi *const analysis) {
50     fBalance         = analysis;
51     }
52   void SetShufflingObject(AliBalancePsi *const analysisShuffled) {
53     fRunShuffling = kTRUE;
54     fShuffledBalance = analysisShuffled;
55   }
56   void SetMixingObject(AliBalancePsi *const analysisMixed) {
57     fRunMixing = kTRUE;
58     fMixedBalance = analysisMixed;
59   }
60   void SetMixingWithEventPlane(Bool_t bMixingWithEventPlane = kTRUE) { fRunMixingEventPlane = bMixingWithEventPlane; }
61   void SetMixingTracks(Int_t tracks) { fMixingTracks = tracks; }
62   void SetAnalysisCutObject(AliESDtrackCuts *const trackCuts) {
63     fESDtrackCuts = trackCuts;}
64   void SetVertexDiamond(Double_t vx, Double_t vy, Double_t vz) {
65     fVxMax = vx;
66     fVyMax = vy;
67     fVzMax = vz;
68   }
69
70   //==============AOD analysis==============//
71   void SetAODtrackCutBit(Int_t bit){
72     fnAODtrackCutBit = bit;
73   }
74
75   void SetKinematicsCutsAOD(Double_t ptmin, Double_t ptmax, Double_t etamin, Double_t etamax){
76     fPtMin  = ptmin;  fPtMax  = ptmax;
77     fEtaMin = etamin; fEtaMax = etamax;
78   }
79
80   void SetExtraDCACutsAOD(Double_t DCAxy, Double_t DCAz){
81     fDCAxyCut  = DCAxy;
82     fDCAzCut = DCAz;
83   }
84
85    void SetExtraTPCCutsAOD(Double_t maxTPCchi2, Int_t minNClustersTPC){
86     fTPCchi2Cut      = maxTPCchi2;
87     fNClustersTPCCut = minNClustersTPC;
88   }
89
90   //==============MC analysis==============//
91   void SetKinematicsCutsMC(Double_t ptmin, Double_t ptmax,
92                            Double_t etamin, Double_t etamax){
93     fPtMin  = ptmin; fPtMax  = ptmax;
94     fEtaMin = etamin; fEtaMax = etamax;
95   }
96   void UseFlowAfterBurner(TF1 *gDifferentialV2) {
97     fDifferentialV2 = gDifferentialV2;
98     fUseFlowAfterBurner = kTRUE;
99   }
100   void ExcludeResonancesInMC() {fExcludeResonancesInMC = kTRUE;}
101   void ExcludeElectronsInMC()  {fExcludeElectronsInMC = kTRUE;}
102
103   void SetPDGCode(Int_t gPdgCode) {
104     fUseMCPdgCode = kTRUE;
105     fPDGCodeToBeAnalyzed = gPdgCode;
106   }
107
108   //Centrality
109   void SetCentralityEstimator(const char* centralityEstimator) {fCentralityEstimator = centralityEstimator;}
110   const char* GetCentralityEstimator(void)  const              {return fCentralityEstimator;}
111   void SetCentralityPercentileRange(Double_t min, Double_t max) { 
112     fUseCentrality = kTRUE;
113     fCentralityPercentileMin=min;
114     fCentralityPercentileMax=max;
115   }
116   void SetImpactParameterRange(Double_t min, Double_t max) { 
117     fUseCentrality = kTRUE;
118     fImpactParameterMin=min;
119     fImpactParameterMax=max;
120   }
121
122   //multiplicity
123   void SetMultiplicityEstimator(const char* multiplicityEstimator) {fMultiplicityEstimator = multiplicityEstimator;}
124   const char* GetMultiplicityEstimator(void)  const              {return fMultiplicityEstimator;}
125   void SetMultiplicityRange(Double_t min, Double_t max) {
126     fUseMultiplicity = kTRUE;
127     fNumberOfAcceptedTracksMin = min;
128     fNumberOfAcceptedTracksMax = max;}
129   
130   // additional event cuts (default = kFALSE)
131   void UseOfflineTrigger() {fUseOfflineTrigger = kTRUE;}
132   void CheckFirstEventInChunk() {fCheckFirstEventInChunk = kTRUE;}
133   void CheckPileUp() {fCheckPileUp = kTRUE;}
134   void CheckPrimaryFlagAOD() {fCheckPrimaryFlagAOD = kTRUE;}
135   void UseMCforKinematics() {fUseMCforKinematics = kTRUE;}
136   void SetCentralityWeights(TH1* hist) { fCentralityWeights = hist; }
137   Bool_t AcceptEventCentralityWeight(Double_t centrality);
138
139   
140   //Acceptance filter
141   void SetAcceptanceParameterization(TF1 *parameterization) {
142     fAcceptanceParameterization = parameterization;}
143
144   //pid
145   enum kDetectorUsedForPID { kTPCpid, kTOFpid, kTPCTOF }; // default TPC & TOF pid (via GetTPCpid & GetTOFpid)  
146   enum kParticleOfInterest { kMuon, kElectron, kPion, kKaon, kProton };  
147
148   void SetUseBayesianPID(Double_t gMinProbabilityValue) {
149     fUsePID = kTRUE; fUsePIDnSigma = kFALSE; fUsePIDPropabilities = kTRUE;
150     fMinAcceptedPIDProbability = gMinProbabilityValue; }
151
152   void SetUseNSigmaPID(Double_t gMaxNSigma) {
153     fUsePID = kTRUE; fUsePIDPropabilities = kFALSE; fUsePIDnSigma = kTRUE;
154     fPIDNSigma = gMaxNSigma; }
155
156   void SetParticleOfInterest(kParticleOfInterest poi) {
157     fParticleOfInterest = poi;}
158   void SetDetectorUsedForPID(kDetectorUsedForPID detConfig) {
159     fPidDetectorConfig = detConfig;}
160     void SetEventClass(TString receivedEventClass){
161         fEventClass = receivedEventClass;
162     }
163     
164   void SetCustomBinning(TString receivedCustomBinning) { fCustomBinning = receivedCustomBinning; }
165
166
167     // electron rejection
168     void SetElectronRejection(Double_t gMaxNSigma){
169       fElectronRejection = kTRUE;
170       fElectronRejectionNSigma = gMaxNSigma;
171     }
172
173     void SetElectronOnlyRejection(Double_t gMaxNSigma){
174       fElectronRejection       = kTRUE;
175       fElectronOnlyRejection   = kTRUE;
176       fElectronRejectionNSigma = gMaxNSigma;
177     }
178
179     void SetElectronRejectionPt(Double_t minPt,Double_t maxPt){
180       fElectronRejectionMinPt  = minPt;
181       fElectronRejectionMaxPt  = maxPt;
182     }
183
184     void SetVZEROCalibrationFile(const char* filename, const char* lhcPeriod);
185
186  private:
187   Double_t    IsEventAccepted(AliVEvent* event);
188   Double_t    GetRefMultiOrCentrality(AliVEvent* event);
189   Double_t    GetReferenceMultiplicityFromAOD(AliVEvent* event);
190   Double_t    GetEventPlane(AliVEvent* event);
191   //===============================correction
192   Double_t    GetTrackbyTrackCorrectionMatrix(Double_t vEta, 
193                                               Double_t vPhi, 
194                                               Double_t vPt, 
195                                               Short_t vCharge, 
196                                               Double_t gCentrality);
197   //===============================correction
198   TObjArray* GetAcceptedTracks(AliVEvent* event, Double_t gCentrality, Double_t gReactionPlane);
199   TObjArray* GetShuffledTracks(TObjArray* tracks, Double_t gCentrality);
200
201   Double_t GetChannelEqualizationFactor(Int_t run, Int_t channel);
202   Double_t GetEqualizationFactor(Int_t run, const char *side);
203  
204   Bool_t fDebugLevel; // debug level
205
206   TClonesArray* fArrayMC; //! AOD object  //+++++++++++++++++++++
207   AliBalancePsi *fBalance; //BF object
208   Bool_t fRunShuffling;//run shuffling or not
209   AliBalancePsi *fShuffledBalance; //BF object (shuffled)
210   Bool_t fRunMixing;//run mixing or not
211   Bool_t fRunMixingEventPlane;//run mixing with Event Plane
212   Int_t  fMixingTracks;
213   AliBalancePsi *fMixedBalance; //TriggeredBF object (mixed)
214   AliEventPoolManager*     fPoolMgr;         //! event pool manager
215
216   TList *fList; //fList object
217   TList *fListBF; //fList object
218   TList *fListBFS; //fList object
219   TList *fListBFM; //fList object
220   TList *fHistListPIDQA;  //! list of histograms
221
222   TH2F *fHistEventStats; //event stats
223   TH2F *fHistCentStats; //centrality stats
224   TH2F *fHistCentStatsUsed; //centrality stats USED +++++++++++++++++++++++
225   TH1F *fHistTriggerStats; //trigger stats
226   TH1F *fHistTrackStats; //Track filter bit stats
227   TH1F *fHistVx; //x coordinate of the primary vertex
228   TH1F *fHistVy; //y coordinate of the primary vertex
229   TH2F *fHistVz; //z coordinate of the primary vertex
230
231   TH2F *fHistTPCvsVZEROMultiplicity; //VZERO vs TPC reference multiplicity
232   TH2F *fHistVZEROSignal; //VZERO channel vs signal
233
234   TH2F *fHistEventPlane; //event plane distribution
235
236   TH2F *fHistClus;//number of clusters (QA histogram)
237   TH2F *fHistDCA;//DCA  (QA histogram)
238   TH2F *fHistChi2;//track chi2 (QA histogram)
239   TH2F *fHistPt;//transverse momentum (QA histogram)
240   TH2F *fHistEta;//pseudorapidity (QA histogram)
241   TH2F *fHistRapidity;//rapidity (QA histogram)
242   TH2F *fHistPhi;//phi (QA histogram)
243   TH3F *fHistEtaPhiPos;//eta-phi pos particles (QA histogram)                    
244   TH3F *fHistEtaPhiNeg;//eta-phi neg particles (QA histogram)
245   TH2F *fHistPhiBefore;//phi before v2 afterburner (QA histogram)
246   TH2F *fHistPhiAfter;//phi after v2 afterburner (QA histogram)
247   TH2F *fHistPhiPos;//phi for positive particles (QA histogram)
248   TH2F *fHistPhiNeg;//phi for negative particles (QA histogram)
249   TH2F *fHistV0M;//V0 multiplicities (QA histogram)
250   TH2F *fHistRefTracks;//reference track multiplicities (QA histogram)
251
252   //============PID============//
253   TH2D *fHistdEdxVsPTPCbeforePID;//TPC dEdx vs momentum before PID cuts (QA histogram)
254   TH2D *fHistBetavsPTOFbeforePID;//beta vs momentum before PID cuts (QA histogram)
255   TH2D *fHistProbTPCvsPtbeforePID; //TPC probability vs pT before PID cuts (QA histogram)
256   TH2D *fHistProbTOFvsPtbeforePID;//TOF probability vs pT before PID cuts (QA histogram)
257   TH2D *fHistProbTPCTOFvsPtbeforePID;//TOF/TPC probability vs pT before PID cuts (QA histogram)
258   TH2D *fHistNSigmaTPCvsPtbeforePID;//TPC nsigma vs pT before PID cuts (QA histogram)
259   TH2D *fHistNSigmaTOFvsPtbeforePID;//TOF nsigma vs pT before PID cuts (QA histogram)
260   TH2D *fHistBetaVsdEdXbeforePID;//TPCTOF  before PID cuts (QA histogram)//+++++++++++++++++++++
261   TH2D *fHistNSigmaTPCTOFvsPtbeforePID;//TPCTOF  before PID cuts (QA histogram)//+++++++++++++++++++++
262   TH2D *fHistdEdxVsPTPCafterPID;//TPC dEdx vs momentum after PID cuts (QA histogram)
263   TH2D *fHistBetavsPTOFafterPID;//beta vs momentum after PID cuts (QA histogram)
264   TH2D *fHistProbTPCvsPtafterPID; //TPC probability vs pT after PID cuts (QA histogram)
265   TH2D *fHistProbTOFvsPtafterPID;//TOF probability vs pT after PID cuts (QA histogram)
266   TH2D *fHistProbTPCTOFvsPtafterPID;//TOF/TPC probability vs pT after PID cuts (QA histogram)
267   TH2D *fHistNSigmaTPCvsPtafterPID;//TPC nsigma vs pT after PID cuts (QA histogram)
268   TH2D *fHistNSigmaTOFvsPtafterPID;//TOF nsigma vs pT after PID cuts (QA histogram)
269   TH2D *fHistBetaVsdEdXafterPID;//TPCTOF  before PID cuts (QA histogram)//+++++++++++++++++++++
270   TH2D *fHistNSigmaTPCTOFvsPtafterPID;//TPCTOF  before PID cuts (QA histogram)//+++++++++++++++++++++
271
272   TH2D *fHistdEdxVsPTPCbeforePIDelectron; //+++++++
273   TH2D *fHistNSigmaTPCvsPtbeforePIDelectron; //+++++++
274   TH2D *fHistdEdxVsPTPCafterPIDelectron; //+++++++
275   TH2D *fHistNSigmaTPCvsPtafterPIDelectron; //+++++++
276   
277   TH3F *fHistCorrectionPlus[kCENTRALITY]; //====correction
278   TH3F *fHistCorrectionMinus[kCENTRALITY]; //===correction
279   Double_t fCentralityArrayForCorrections[kCENTRALITY];
280   Int_t fCentralityArrayBinsForCorrections;
281
282   TH1* fCentralityWeights;                   // for centrality flattening
283
284   AliPIDResponse *fPIDResponse;     //! PID response object
285   AliPIDCombined       *fPIDCombined;     //! combined PID object
286   
287   kParticleOfInterest  fParticleOfInterest;//analyzed particle
288   kDetectorUsedForPID   fPidDetectorConfig;//used detector for PID
289
290   Bool_t fUsePID; //flag to use PID 
291   Bool_t fUsePIDnSigma;//flag to use nsigma method for PID
292   Bool_t fUsePIDPropabilities;//flag to use probability method for PID
293   Double_t fPIDNSigma;//nsigma cut for PID
294   Double_t fMinAcceptedPIDProbability;//probability cut for PID
295
296   Bool_t   fElectronRejection;//flag to use electron rejection
297   Bool_t   fElectronOnlyRejection;//flag to use electron rejection with exclusive electron PID (no other particle in nsigma range)
298   Double_t fElectronRejectionNSigma;//nsigma cut for electron rejection
299   Double_t fElectronRejectionMinPt;//minimum pt for electron rejection (default = 0.)
300   Double_t fElectronRejectionMaxPt;//maximum pt for electron rejection (default = 1000.)
301   //============PID============//
302
303   AliESDtrackCuts *fESDtrackCuts; //ESD track cuts
304
305   TString fCentralityEstimator;      //"V0M","TRK","TKL","ZDC","FMD"
306   Bool_t fUseCentrality;//use the centrality (PbPb) or not (pp)
307   Double_t fCentralityPercentileMin;//centrality percentile min
308   Double_t fCentralityPercentileMax;//centrality percentile max
309   Double_t fImpactParameterMin;//impact parameter min (used for MC)
310   Double_t fImpactParameterMax;//impact parameter max (used for MC)
311
312   TString fMultiplicityEstimator;//"V0M","V0A","V0C","TPC"
313   Bool_t fUseMultiplicity;//use the multiplicity cuts
314   Double_t fNumberOfAcceptedTracksMin;//min. number of number of accepted tracks (used for the multiplicity dependence study - pp)
315   Double_t fNumberOfAcceptedTracksMax;//max. number of number of accepted tracks (used for the multiplicity dependence study - pp)
316   TH2F *fHistNumberOfAcceptedTracks;//hisot to store the number of accepted tracks
317   TH1F *fHistMultiplicity;//hisot to store the number of accepted tracks //++++++++++++++++++
318
319   Bool_t fUseOfflineTrigger;//Usage of the offline trigger selection
320   Bool_t fCheckFirstEventInChunk;//Usage of the "First Event in Chunk" check (not needed for new productions)
321   Bool_t fCheckPileUp;//Usage of the "Pile-Up" event check
322   Bool_t fCheckPrimaryFlagAOD;// Usage of check on AliAODtrack::kPrimary (default = OFF)
323   Bool_t fUseMCforKinematics;//Usage of MC information for filling the kinematics information of particles (only in MCAODrec mode)
324
325   Double_t fVxMax;//vxmax
326   Double_t fVyMax;//vymax
327   Double_t fVzMax;//vzmax
328
329   Int_t fnAODtrackCutBit;//track cut bit from track selection (only used for AODs)
330
331   Double_t fPtMin;//only used for AODs
332   Double_t fPtMax;//only used for AODs
333   Double_t fPtMinForCorrections;//only used for AODs
334   Double_t fPtMaxForCorrections;//only used for AODs
335   Double_t fPtBinForCorrections; //=================================correction
336   Double_t fEtaMin;//only used for AODs
337   Double_t fEtaMax;//only used for AODs
338   Double_t fEtaMinForCorrections;//only used for AODs
339   Double_t fEtaMaxForCorrections;//only used for AODs
340   Double_t fEtaBinForCorrections; //=================================correction
341   Double_t fPhiMin; //=================================correction 
342   Double_t fPhiMax; //=================================correction
343   Double_t fPhiMinForCorrections;//only used for AODs
344   Double_t fPhiMaxForCorrections;//only used for AODs
345   Double_t fPhiBinForCorrections; //=================================correction
346
347   Double_t fDCAxyCut;//only used for AODs
348   Double_t fDCAzCut;//only used for AODs
349
350   Double_t fTPCchi2Cut;//only used for AODs
351   Int_t fNClustersTPCCut;//only used for AODs
352
353   TF1 *fAcceptanceParameterization;//acceptance filter used for MC
354
355   TF1 *fDifferentialV2;//pt-differential v2 (from real data)
356   Bool_t fUseFlowAfterBurner;//Usage of a flow after burner
357
358   Bool_t fExcludeResonancesInMC;//flag to exclude the resonances' decay products (and conversion) from the MC analysis
359   Bool_t fExcludeElectronsInMC;//flag to exclude the electrons from the MC analysis
360   Bool_t fUseMCPdgCode; //Boolean to analyze a set of particles in MC
361   Int_t fPDGCodeToBeAnalyzed; //Analyze a set of particles in MC
362   TString fEventClass; //Can be "EventPlane", "Centrality", "Multiplicity"
363   TString fCustomBinning;//for setting customized binning (for output AliTHn of AliBalancePsi)
364   
365   //VZERO calibration
366   TH1F *fHistVZEROAGainEqualizationMap;//VZERO calibration map
367   TH1F *fHistVZEROCGainEqualizationMap;//VZERO calibration map
368   TH2F *fHistVZEROChannelGainEqualizationMap; //VZERO calibration map
369
370   AliAnalysisTaskBFPsi(const AliAnalysisTaskBFPsi&); // not implemented
371   AliAnalysisTaskBFPsi& operator=(const AliAnalysisTaskBFPsi&); // not implemented
372   
373   ClassDef(AliAnalysisTaskBFPsi, 6); // example of analysis
374 };
375
376
377
378 #endif