]> 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                 kMassCut,
50                 kNCuts
51         };
52
53
54         enum pionCuts {
55                 kPionIn=0,
56                 kNoTracks,
57                 kTrackCuts,
58                 kdEdxCuts,
59                 kPionOut
60         };
61
62
63         Bool_t SetCutIds(TString cutString); 
64         Int_t fCuts[kNCuts];
65         Bool_t SetCut(cutIds cutID, Int_t cut);
66         Bool_t UpdateCutString();
67         static const char * fgkCutNames[kNCuts];
68
69
70         Bool_t InitializeCutsFromCutString(const TString analysisCutSelection); 
71         
72
73         AliPrimaryPionCuts(const char *name="PionCuts", const char * title="Pion Cuts");
74         virtual ~AliPrimaryPionCuts();                            //virtual destructor
75
76         virtual Bool_t IsSelected(TObject* /*obj*/){return kTRUE;}
77         virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
78
79         TString GetCutNumber();
80
81                 // Cut Selection
82         Bool_t PionIsSelectedMC(Int_t labelParticle,AliStack *fMCStack);
83         Bool_t TrackIsSelected(AliESDtrack* lTrack);
84         Bool_t PionIsSelected(AliESDtrack* lTrack);
85         static AliPrimaryPionCuts * GetStandardCuts2010PbPb();
86         static AliPrimaryPionCuts * GetStandardCuts2010pp();
87         Bool_t InitPIDResponse();
88         
89         void SetPIDResponse(AliPIDResponse * pidResponse) {fPIDResponse = pidResponse;}
90         AliPIDResponse * GetPIDResponse() { return fPIDResponse;}
91         
92         void PrintCuts();
93
94         void InitCutHistograms(TString name="",Bool_t preCut = kTRUE,TString cutName="");
95         void SetFillCutHistograms(TString name="",Bool_t preCut = kTRUE,TString cutName=""){if(!fHistograms){InitCutHistograms(name,preCut,cutName);};}
96         TList *GetCutHistograms(){return fHistograms;}
97
98         static AliVTrack * GetTrack(AliVEvent * event, Int_t label);
99
100         ///Cut functions
101         Bool_t dEdxCuts(AliVTrack * track);
102
103         Bool_t SetTPCdEdxCutPionLine(Int_t pidedxSigmaCut);
104         Bool_t SetITSdEdxCutPionLine(Int_t ededxSigmaCut);
105         Bool_t SetITSClusterCut(Int_t clsITSCut);
106         Bool_t SetTPCClusterCut(Int_t clsTPCCut);
107         Bool_t SetEtaCut(Int_t etaCut);
108         Bool_t SetPtCut(Int_t ptCut);
109         Bool_t SetDCACut(Int_t dcaCut);
110         void SetEtaShift(Double_t etaShift){fEtaShift = etaShift;}
111         Bool_t SetTOFPionPIDCut(Int_t TOFelectronPID);
112         Bool_t SetMassCut(Int_t massCut);
113         Double_t GetMassCut(){return fMassCut;}
114         
115         // Request Flags
116         Double_t GetEtaCut(){ return  fEtaCut;}
117         Double_t GetNFindableClustersTPC(AliESDtrack* lTrack);
118         Bool_t   DoWeights(){return fDoWeights;}
119         Bool_t   DoMassCut(){return fDoMassCut;}
120         
121         protected:
122
123         TList *fHistograms;
124         AliPIDResponse *fPIDResponse;
125         AliESDtrackCuts *fEsdTrackCuts;
126
127         Double_t fEtaCut; //eta cutç
128         Double_t fEtaShift;
129         Bool_t   fDoEtaCut;
130         Double_t fPtCut;
131         Double_t fMinClsTPC; // minimum clusters in the TPC
132         Double_t fMinClsTPCToF; // minimum clusters to findable clusters
133         Bool_t   fDodEdxSigmaITSCut; // flag to use the dEdxCut ITS based on sigmas
134         Bool_t   fDodEdxSigmaTPCCut; // flag to use the dEdxCut TPC based on sigmas
135         Bool_t   fDoTOFsigmaCut; // flag to use TOF pid cut RRnewTOF
136         Double_t fPIDnSigmaAbovePionLineITS; // sigma cut ITS
137         Double_t fPIDnSigmaBelowPionLineITS; // sigma cut ITS
138         Double_t fPIDnSigmaAbovePionLineTPC; // sigma cut TPC
139         Double_t fPIDnSigmaBelowPionLineTPC; // sigma cut TPC
140         Double_t fPIDnSigmaAbovePionLineTOF; // sigma cut TOF
141         Double_t fPIDnSigmaBelowPionLineTOF; // sigma cut TOF 
142         Bool_t   fUseCorrectedTPCClsInfo; // flag to use corrected tpc cl info
143         Bool_t   fUseTOFpid; // flag to use tof pid
144         Bool_t   fRequireTOF; //flg to analyze only tracks with TOF signal
145         Bool_t   fDoMassCut;
146         Double_t fMassCut;      
147         Bool_t   fDoWeights;
148         
149
150
151         // Histograms
152         TObjString *fCutString; // cut number used for analysis
153         TH1F *fHistCutIndex; // bookkeeping for cuts
154         TH1F *fHistdEdxCuts;  // bookkeeping for dEdx cuts
155         TH2F *fHistITSdEdxbefore; // ITS dEdx before cuts
156         TH2F *fHistITSdEdxafter;
157         TH2F *fHistTPCdEdxbefore; // TPC dEdx before cuts
158         TH2F *fHistTPCdEdxafter; // TPC dEdx after cuts
159         TH2F *fHistTPCdEdxSignalbefore; //TPC dEdx signal before
160         TH2F *fHistTPCdEdxSignalafter; //TPC dEdx signal  after
161         TH2F *fHistTOFbefore; // TOF after cuts
162         TH2F *fHistTOFafter; // TOF after cuts
163         TH2F *fHistTrackDCAxyPtbefore;
164         TH2F *fHistTrackDCAxyPtafter;
165         TH2F *fHistTrackDCAzPtbefore;
166         TH2F *fHistTrackDCAzPtafter;
167         TH2F *fHistTrackNFindClsPtTPCbefore;
168         TH2F *fHistTrackNFindClsPtTPCafter;
169         
170         private:
171
172         AliPrimaryPionCuts(const AliPrimaryPionCuts&); // not implemented
173         AliPrimaryPionCuts& operator=(const AliPrimaryPionCuts&); // not implemented
174
175
176         ClassDef(AliPrimaryPionCuts,3)
177 };
178
179 #endif