]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliPrimaryPionCuts.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliPrimaryPionCuts.h
1 #ifndef ALIPRIMARYPIONCUTS_H
2 #define ALIPRIMARYPIONCUTS_H
3
4 // Class handling all kinds of selection cuts for primary
5
6 // Authors: Svein Lindal, Daniel Lohner                                                                                         *
7
8
9 #include "AliAODpidUtil.h"
10 #include "AliAODTrack.h"
11 #include "AliESDtrack.h"
12 #include "AliVTrack.h"
13 #include "AliAODTrack.h"
14 #include "AliStack.h"
15 #include "AliAnalysisCuts.h"
16 #include "AliESDtrackCuts.h"
17 #include "TH1F.h"
18
19 class AliESDEvent;
20 class AliAODEvent;
21 class AliConversionPhotonBase;
22 class AliKFVertex;
23 class AliKFParticle;
24 class TH1F;
25 class TH2F;
26 class AliPIDResponse;
27 class AliAnalysisCuts;
28 class iostream;
29 class TList;
30 class AliAnalysisManager;
31
32
33 using namespace std;
34
35 class AliPrimaryPionCuts : public AliAnalysisCuts {
36                 
37         public: 
38
39
40         enum cutIds {
41                 kEtaCut,
42                 kClsITSCut,
43                 kClsTPCCut,
44                 kDCACut,
45                 kPtCut,
46                 kPidedxSigmaITSCut,
47                 kPidedxSigmaTPCCut,
48                 kPiTOFSigmaPID,
49                 kNCuts
50         };
51
52
53         enum pionCuts {
54                 kPionIn=0,
55                 kNoTracks,
56                 kTrackCuts,
57                 kdEdxCuts,
58                 kPionOut
59         };
60
61
62         Bool_t SetCutIds(TString cutString); 
63         Int_t fCuts[kNCuts];
64         Bool_t SetCut(cutIds cutID, Int_t cut);
65         Bool_t UpdateCutString();
66         static const char * fgkCutNames[kNCuts];
67
68
69         Bool_t InitializeCutsFromCutString(const TString analysisCutSelection); 
70         
71
72         AliPrimaryPionCuts(const char *name="PionCuts", const char * title="Pion Cuts");
73         virtual ~AliPrimaryPionCuts();                            //virtual destructor
74
75         virtual Bool_t IsSelected(TObject* /*obj*/){return kTRUE;}
76         virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
77
78         TString GetCutNumber();
79
80                 // Cut Selection
81         Bool_t PionIsSelectedMC(Int_t labelParticle,AliStack *fMCStack);
82         Bool_t TrackIsSelected(AliESDtrack* lTrack);
83         Bool_t PionIsSelected(AliESDtrack* lTrack);
84         static AliPrimaryPionCuts * GetStandardCuts2010PbPb();
85         static AliPrimaryPionCuts * GetStandardCuts2010pp();
86         Bool_t InitPIDResponse();
87         
88         void SetPIDResponse(AliPIDResponse * pidResponse) {fPIDResponse = pidResponse;}
89         AliPIDResponse * GetPIDResponse() { return fPIDResponse;}
90         
91         void PrintCuts();
92
93         void InitCutHistograms(TString name="",Bool_t preCut = kTRUE,TString cutName="");
94         void SetFillCutHistograms(TString name="",Bool_t preCut = kTRUE,TString cutName=""){if(!fHistograms){InitCutHistograms(name,preCut,cutName);};}
95         TList *GetCutHistograms(){return fHistograms;}
96
97         static AliVTrack * GetTrack(AliVEvent * event, Int_t label);
98
99         ///Cut functions
100         Bool_t dEdxCuts(AliVTrack * track);
101
102         Bool_t SetTPCdEdxCutPionLine(Int_t pidedxSigmaCut);
103         Bool_t SetITSdEdxCutPionLine(Int_t ededxSigmaCut);
104         Bool_t SetITSClusterCut(Int_t clsITSCut);
105         Bool_t SetTPCClusterCut(Int_t clsTPCCut);
106         Bool_t SetEtaCut(Int_t etaCut);
107         Bool_t SetPtCut(Int_t ptCut);
108         Bool_t SetDCACut(Int_t dcaCut);
109         void SetEtaShift(Double_t etaShift){fEtaShift = etaShift;}
110         Bool_t SetTOFPionPIDCut(Int_t TOFelectronPID);
111         
112         // Request Flags
113         Double_t GetEtaCut(){ return  fEtaCut;}
114         Double_t GetNFindableClustersTPC(AliESDtrack* lTrack);
115         Bool_t   DoWeights(){return fDoWeights;}
116         
117         protected:
118
119         TList *fHistograms;
120         AliPIDResponse *fPIDResponse;
121         AliESDtrackCuts *fEsdTrackCuts;
122
123         Double_t fEtaCut; //eta cutç
124         Double_t fEtaShift;
125         Bool_t   fDoEtaCut;
126         Double_t fPtCut;
127         Double_t fMinClsTPC; // minimum clusters in the TPC
128         Double_t fMinClsTPCToF; // minimum clusters to findable clusters
129         Bool_t   fDodEdxSigmaITSCut; // flag to use the dEdxCut ITS based on sigmas
130         Bool_t   fDodEdxSigmaTPCCut; // flag to use the dEdxCut TPC based on sigmas
131         Bool_t   fDoTOFsigmaCut; // flag to use TOF pid cut RRnewTOF
132         Double_t fPIDnSigmaAbovePionLineITS; // sigma cut ITS
133         Double_t fPIDnSigmaBelowPionLineITS; // sigma cut ITS
134         Double_t fPIDnSigmaAbovePionLineTPC; // sigma cut TPC
135         Double_t fPIDnSigmaBelowPionLineTPC; // sigma cut TPC
136         Double_t fPIDnSigmaAbovePionLineTOF; // sigma cut TOF
137         Double_t fPIDnSigmaBelowPionLineTOF; // sigma cut TOF 
138         Bool_t   fUseCorrectedTPCClsInfo; // flag to use corrected tpc cl info
139         Bool_t   fUseTOFpid; // flag to use tof pid
140         Bool_t   fRequireTOF; //flg to analyze only tracks with TOF signal
141         Bool_t   fDoWeights;
142
143
144         // Histograms
145         TObjString *fCutString; // cut number used for analysis
146         TH1F *fHistCutIndex; // bookkeeping for cuts
147         TH1F *fHistdEdxCuts;  // bookkeeping for dEdx cuts
148         TH2F *fHistITSdEdxbefore; // ITS dEdx before cuts
149         TH2F *fHistITSdEdxafter;
150         TH2F *fHistTPCdEdxbefore; // TPC dEdx before cuts
151         TH2F *fHistTPCdEdxafter; // TPC dEdx after cuts
152         TH2F *fHistTPCdEdxSignalbefore; //TPC dEdx signal before
153         TH2F *fHistTPCdEdxSignalafter; //TPC dEdx signal  after
154         TH2F *fHistTOFbefore; // TOF after cuts
155         TH2F *fHistTOFafter; // TOF after cuts
156         TH2F *fHistTrackDCAxyPtbefore;
157         TH2F *fHistTrackDCAxyPtafter;
158         TH2F *fHistTrackDCAzPtbefore;
159         TH2F *fHistTrackDCAzPtafter;
160         TH2F *fHistTrackNFindClsPtTPCbefore;
161         TH2F *fHistTrackNFindClsPtTPCafter;
162         
163         private:
164
165         AliPrimaryPionCuts(const AliPrimaryPionCuts&); // not implemented
166         AliPrimaryPionCuts& operator=(const AliPrimaryPionCuts&); // not implemented
167
168
169         ClassDef(AliPrimaryPionCuts,2)
170 };
171
172 #endif