]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/hfe/AliAnalysisTaskEMCalHFEpA.h
updated the task
[u/mrichter/AliRoot.git] / PWGHF / hfe / AliAnalysisTaskEMCalHFEpA.h
1 #ifndef AliAnalysisTaskEMCalHFEpA_cxx
2 #define AliAnalysisTaskEMCalHFEpA_cxx
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 ////////////////////////////////////////////////////////////////////////
8 //                                                                    //
9 //      Task for Heavy-flavour electron analysis in pPb collisions    //
10 //      (+ Electron-Hadron Jetlike Azimuthal Correlation)             //
11 //                                                                                                                                        //
12 //              v1.0                                                                                                              //
13 //                                                                    //
14 //          Authors                                                                               //
15 //              Elienos Pereira de Oliveira Filho (epereira@cern.ch)          //
16 //              Cristiane Jahnke                (cristiane.jahnke@cern.ch)                    //
17 //                                                                    //
18 ////////////////////////////////////////////////////////////////////////
19
20 class TH1F;
21 class TH2F;
22 class AliESDEvent;
23 class AliESDtrackCuts;
24 class AliESDtrack;
25 class AliHFEcontainer;
26 class AliHFEcuts;
27 class AliHFEpid;
28 class AliHFEpidQAmanager;
29 class AliCFManager;
30 class AliPIDResponse;
31 class AliCentrality;
32 class AliAODEvent;
33 class AliVEvent;
34 class AliAODMCHeader;
35 class AliSelectNonHFE;
36 class AliEventPoolManager;
37 class AliEventPool;
38 class TObjArray;
39
40 //______________________________________________________________________
41 //Library
42 #include "AliAnalysisTaskSE.h"
43 #include "AliHFEpid.h"
44 #include "AliLog.h"
45 //______________________________________________________________________
46
47 //______________________________________________________________________
48 class AliAnalysisTaskEMCalHFEpA : public AliAnalysisTaskSE 
49 {
50 //______________________________________________________________________
51         public:
52         AliAnalysisTaskEMCalHFEpA();
53         AliAnalysisTaskEMCalHFEpA(const char *name);
54         virtual ~AliAnalysisTaskEMCalHFEpA();
55   
56         virtual void   UserCreateOutputObjects();
57         virtual void   UserExec(Option_t *option);
58         virtual void   Terminate(Option_t *);
59
60         //Setters
61         void SetHFECuts(AliHFEcuts * const cuts) {fCuts = cuts;};
62         void SetRejectKinkMother(Bool_t rejectKinkMother = kFALSE) {fRejectKinkMother = rejectKinkMother;};
63         void SetCorrelationAnalysis(Bool_t CorrelationFlag=kTRUE) {fCorrelationFlag = CorrelationFlag;};
64         void SetMCanalysis() {fIsMC = kTRUE;};
65         void SetCentrality(Double_t CentralityMin, Double_t CentralityMax) { fCentralityMin = CentralityMin; fCentralityMax = CentralityMax; fHasCentralitySelection = kTRUE; };
66         void SetAODanalysis(Bool_t IsAOD) {fIsAOD = IsAOD;};
67         void SetEventMixing(Bool_t EventMixingFlag) { fEventMixingFlag = EventMixingFlag;};
68         void SetNonHFEmassCut(Double_t MassCut) { fMassCut = MassCut; fMassCutFlag = kTRUE;};
69         void SetEtaCut(Double_t EtaCutMin,Double_t EtaCutMax ) { fEtaCutMin = EtaCutMin; fEtaCutMax = EtaCutMax; };
70         
71         void SetdPhidEtaCut(Double_t dPhiCut, Double_t dEtaCut ) { fdPhiCut = dPhiCut;fdEtaCut = dEtaCut ;};
72         
73         void SetEoverPCut(Double_t EoverPCutMin,Double_t EoverPCutMax ) { fEoverPCutMin = EoverPCutMin; fEoverPCutMax = EoverPCutMax; };
74         
75         void SetM02Cut(Double_t M02CutMin,Double_t M02CutMax ) { fM02CutMin = M02CutMin; fM02CutMax = M02CutMax; };
76         void SetM20Cut(Double_t M20CutMin,Double_t M20CutMax ) { fM20CutMin = M20CutMin; fM20CutMax = M20CutMax; };
77         
78         
79         void SetNonHFEangleCut(Double_t AngleCut) { fAngleCut = AngleCut; fAngleCutFlag = kTRUE;};
80         void SetNonHFEchi2Cut(Double_t Chi2Cut) { fChi2Cut = Chi2Cut; fChi2CutFlag = kTRUE;};
81         void SetNonHFEdcaCut(Double_t DCAcut) { fDCAcut = DCAcut; fDCAcutFlag = kTRUE;};
82         void SetUseEMCal() { fUseEMCal=kTRUE;};
83         void SetUseShowerShapeCut(Bool_t UseShowerShapeCut=kFALSE) { fUseShowerShapeCut=UseShowerShapeCut;};
84         void SetBackground(Bool_t FillBackground=kFALSE) { fFillBackground=FillBackground;};
85         void SetEMCalTriggerEG1() { fEMCEG1=kTRUE; };
86         void SetEMCalTriggerEG2() { fEMCEG2=kTRUE; };
87         void SetCentralityEstimator(Int_t Estimator) { fEstimator=Estimator; }; //0 = V0A, 1 = Other
88         
89         //Getters
90         AliHFEpid *GetPID() const {return fPID;};
91 //______________________________________________________________________
92   
93 //______________________________________________________________________
94         private:
95         
96 //Function to process track cuts
97         Bool_t ProcessCutStep(Int_t cutStep, AliVParticle *track);
98 //Function to process eh analysis
99         void ElectronHadronCorrelation(AliVTrack *track, Int_t trackIndex, AliVParticle *vtrack);
100 //Function to find non-HFE and fill histos
101         void Background(AliVTrack *track, Int_t trackIndex, AliVParticle *vtrack);
102 //Selected Hadrons, for mixed event analysis
103         TObjArray* SelectedHadrons();
104 //DiHadron Correlation Background
105         void DiHadronCorrelation(AliVTrack *track, Int_t trackIndex);
106 //Find Mothers (Finde HFE and NonHFE from MC information)
107         Bool_t FindMother(Int_t mcIndex);
108
109 //Flags for specifics analysis
110         Bool_t                          fCorrelationFlag;
111         Bool_t                          fIsMC;
112         Bool_t                          fUseEMCal;
113         Bool_t                          fUseShowerShapeCut;
114         Bool_t                          fFillBackground;
115
116         Bool_t                          fEMCEG1;
117         Bool_t                          fEMCEG2;
118
119 //Used in the function FindMother
120         Bool_t                          fIsHFE1;
121         Bool_t                          fIsHFE2;
122         Bool_t                          fIsNonHFE;
123         Bool_t                          fIsFromD;
124         Bool_t                          fIsFromB;
125         Bool_t                          fIsFromPi0;
126         Bool_t                          fIsFromEta;
127         Bool_t                          fIsFromGamma;
128         
129 //General variables
130         AliESDEvent                     *fESD;
131         AliAODEvent                     *fAOD;                          /// new
132         AliVEvent                       *fVevent;                       /// new
133         AliESDtrackCuts         *fPartnerCuts;
134         TList                           *fOutputList;
135         AliPIDResponse                  *fPidResponse;
136         AliSelectNonHFE                 *fNonHFE;
137         
138 //For the case of AOD analysis
139         Bool_t                                  fIsAOD;                                 //flag for AOD analysis
140         
141 //For Centrality Selection
142         AliCentrality                   *fCentrality;
143         Double_t                                fCentralityMin;
144         Double_t                                fCentralityMax;
145         Bool_t                                  fHasCentralitySelection;
146         TH1F                                    *fCentralityHist;
147         TH1F                                    *fCentralityHistPass;
148         Float_t                                 fZvtx;  
149         Int_t                                   fEstimator;
150         
151 //EMCal
152         
153         AliVCluster                             *fClus;
154         //AliESDCaloCluster             *fClusESD;
155         
156 //Histograms
157         TH1F                            *fNevent;
158         TH1F                            *fPtElec_Inc;
159
160         
161         TH1F                            *fCharge_n;
162         TH1F                            *fCharge_p;
163         
164         TH2D                            *fTime;
165         TH2D                            *fTime2;
166         TH2D                            *ftimingEle;
167         TH2D                            *ftimingEle2;   
168         
169         TH1F                            *fPtElec_ULS;
170         TH1F                            *fPtElec_LS;
171         
172         //PID Histograms
173         TH1F                            *fpid;          
174         
175         TH2F                            **fEoverP_pt;
176         TH2F                            **fEoverP_tpc;
177         
178         TH1F                            **fTPC_pt;
179         TH2F                            **fTPC_p;
180         
181         TH1F                            **fTPCnsigma_pt;
182         TH2F                            **fTPCnsigma_p;
183         TH2F                            *fTPCnsigma_pt_2D;
184         TH2F                            *fShowerShapeCut;
185
186         
187         TH2F                            *fTPCnsigma_eta;
188         TH2F                            *fTPCnsigma_phi;
189         
190
191         TH1F                            **fECluster;
192         TH2F                            **fEtaPhi;
193         TH1F                            **fVtxZ;
194         TH1F                            **fNTracks;
195         TH1F                            **fNClusters;
196         TH2F                            **fTPCNcls_EoverP;
197         
198         TH1F                            **fEta;
199         TH1F                            **fPhi;
200         TH1F                            **fR;
201         TH2F                            **fR_EoverP;
202         TH1F                            **fNcells;
203         TH2F                            **fNcells_EoverP;
204         TH1F                            **fNcells_electrons;
205         TH1F                            **fNcells_hadrons;
206         TH1F                            **fECluster_ptbins;
207         TH1F                            **fEoverP_ptbins;
208         TH1F                            **fEoverP_wSSCut;
209         TH2F                            **fM02_EoverP;
210         TH2F                            **fM20_EoverP;
211         TH2F                            **fTPCnsigma_eta_electrons;
212         TH2F                            **fTPCnsigma_eta_hadrons;
213         
214         TH2F                            *fEoverP_pt_pions;
215         
216         TH2F                            *ftpc_p_EoverPcut;
217         TH2F                            *fnsigma_p_EoverPcut;
218         
219         TH2F                            *fEoverP_pt_pions2;
220         TH2F                            *fNcells_pt;
221         TH2F                            *fEoverP_pt_hadrons;
222         
223         //Electron-Hadron Correlation Histograms
224         TH2F                            **fCEtaPhi_Inc;
225         
226         TH2F                            **fCEtaPhi_ULS;
227         TH2F                            **fCEtaPhi_LS;
228         TH2F                            **fCEtaPhi_ULS_NoP;
229         TH2F                            **fCEtaPhi_LS_NoP;
230         
231         TH2F                            **fCEtaPhi_ULS_Weight;
232         TH2F                            **fCEtaPhi_LS_Weight;
233         TH2F                            **fCEtaPhi_ULS_NoP_Weight;
234         TH2F                            **fCEtaPhi_LS_NoP_Weight;
235         
236         TH1F                            *fInvMass;
237         TH1F                            *fInvMassBack;
238         TH1F                            *fDCA;
239         TH1F                            *fDCABack;
240         TH1F                            *fOpAngle;
241         TH1F                            *fOpAngleBack;
242         
243         Double_t                        fMassCut;
244         Double_t                        fEtaCutMin;
245         Double_t                        fEtaCutMax;
246         
247         Double_t                        fdPhiCut;
248         Double_t                        fdEtaCut;
249         
250         Double_t                        fEoverPCutMin;
251         Double_t                        fEoverPCutMax;
252         Double_t                        fM20CutMin;
253         Double_t                        fM20CutMax;
254         Double_t                        fM02CutMin;
255         Double_t                        fM02CutMax;
256         
257         Double_t                        fAngleCut;
258         Double_t                        fChi2Cut;
259         Double_t                        fDCAcut;
260         Bool_t                          fMassCutFlag;
261         Bool_t                          fAngleCutFlag;
262         Bool_t                          fChi2CutFlag;
263         Bool_t                          fDCAcutFlag;
264         
265         //Non-HFE reconstruction efficiency
266         TH1F                            *fPtBackgroundBeforeReco;
267         TH1F                            *fPtBackgroundAfterReco;        
268         //Tracking Efficiency
269         TH1F                            *fPtMCparticleAll;
270         TH1F                            *fPtMCparticleReco;
271         TH1F                            *fPtMCparticleAllHfe1;
272         TH1F                            *fPtMCparticleRecoHfe1;
273         TH1F                            *fPtMCparticleAllHfe2;
274         TH1F                            *fPtMCparticleRecoHfe2;
275         TH1F                            *fPtMCelectronAfterAll;
276         
277         TH1F                            *fPtMCpi0;
278         
279         TH1F                            *fPtMC_EMCal_All;
280         TH1F                            *fPtMC_EMCal_Selected;
281         TH1F                            *fPtMC_TPC_All;
282         TH1F                            *fPtMC_TPC_Selected;
283         
284         TH1F                            *fPtMCWithLabel;
285         TH1F                            *fPtMCWithoutLabel;
286         TH1F                            *fPtIsPhysicaPrimary;
287         
288 //For the HFE package
289         AliHFEcuts                      *fCuts;                                 // Cut Collection for HFE
290         AliCFManager            *fCFM;                                  // Correction Framework Manager
291         AliHFEpid                       *fPID;                                  // PID
292         AliHFEpidQAmanager      *fPIDqa;                                                // PID QA manager
293         
294 //Others
295         AliStack                        *fMCstack;                                              //
296         Bool_t              fRejectKinkMother;                          //
297         TParticle                       *fMCtrack;
298         TParticle                       *fMCtrackMother;
299         TParticle                       *fMCtrackGMother;
300         TParticle                       *fMCtrackGGMother;
301         TParticle                       *fMCtrackGGGMother;
302         TClonesArray            *fMCarray;
303         AliAODMCHeader          *fMCheader;
304         AliAODMCParticle        *fMCparticle;
305         AliAODMCParticle        *fMCparticleMother;
306         AliAODMCParticle        *fMCparticleGMother;
307         AliAODMCParticle        *fMCparticleGGMother;
308         AliAODMCParticle        *fMCparticleGGGMother;
309         AliMCEventHandler       *fEventHandler;
310         AliMCEvent                      *fMCevent;
311
312 //______________________________________________________________________
313 //Mixed event analysis
314         AliEventPoolManager *fPoolMgr;
315         AliEventPool            *fPool;
316         TObjArray                       *fTracksClone;
317         TObjArray                       *fTracks;
318         
319         TH2F                            **fCEtaPhi_Inc_EM;
320         
321         TH2F                            **fCEtaPhi_ULS_EM;
322         TH2F                            **fCEtaPhi_LS_EM;
323         
324         TH2F                            **fCEtaPhi_ULS_Weight_EM;
325         TH2F                            **fCEtaPhi_LS_Weight_EM;
326         
327         TH1F                            *fPoolNevents;
328         
329         Bool_t                          fEventMixingFlag;
330 //______________________________________________________________________
331
332 //______________________________________________________________________
333 //Di-hadron correlation
334         TH2F                            **fCEtaPhi_Inc_DiHadron;
335         TH1F                            *fPtTrigger_Inc;
336 //______________________________________________________________________
337
338         AliAnalysisTaskEMCalHFEpA(const AliAnalysisTaskEMCalHFEpA&);                    // not implemented
339         AliAnalysisTaskEMCalHFEpA& operator=(const AliAnalysisTaskEMCalHFEpA&);                 // not implemented
340   
341         ClassDef(AliAnalysisTaskEMCalHFEpA, 1);                                                                 // example of analysis
342 //______________________________________________________________________
343 };
344
345 ///_________________________________________________________________________________________________
346 ///Class copied from : $ALICE_ROOT/PWGCF/Correlations/DPhi/AliAnalysisTaskLongRangeCorrelations.h
347 ///Author: Christoph Mayer
348 class AliEHCParticle : public TObject {
349 public:
350   AliEHCParticle(Double_t eta=0, Double_t phi=0, Double_t pt=0)
351     : fEta(eta), fPhi(phi), fPt(pt) {}
352   virtual ~AliEHCParticle() {}
353
354   Double_t Eta() const { return fEta; }
355   Double_t Phi() const { return fPhi; }
356   Double_t Pt() const { return fPt; }
357
358 protected:
359 private:
360   AliEHCParticle(const AliEHCParticle&);
361   AliEHCParticle& operator=(const AliEHCParticle&);
362
363   Double_t fEta;
364   Double_t fPhi;
365   Double_t fPt;
366   
367   ClassDef(AliEHCParticle, 1);
368 } ;
369 ///_________________________________________________________________________________________________
370
371 #endif