]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliDalitzElectronCuts.h
Fix in filling the hit-map of ESDtrack
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliDalitzElectronCuts.h
1 #ifndef ALIDALITZELECTRONCUTS_H
2 #define ALIDALITZELECTRONCUTS_H
3
4 // Class handling all kinds of selection cuts for electrons
5
6 // Authors: Svein Lindal, Daniel Lohner                                                                                         *
7
8
9 #include "AliAODpidUtil.h"
10 //#include "AliConversionPhotonBase.h"
11 //#include "AliAODConversionMother.h"
12 #include "AliAODTrack.h"
13 #include "AliESDtrack.h"
14 #include "AliVTrack.h"
15 #include "AliAODTrack.h"
16 #include "AliStack.h"
17 #include "AliAnalysisCuts.h"
18 #include "AliESDtrackCuts.h"
19 #include "TH1F.h"
20
21 class AliESDEvent;
22 class AliAODEvent;
23 class AliConversionPhotonBase;
24 class AliKFVertex;
25 class AliKFParticle;
26 class TH1F;
27 class TH2F;
28 class AliPIDResponse;
29 class AliAnalysisCuts;
30 class iostream;
31 class TList;
32 class AliAnalysisManager;
33
34
35 using namespace std;
36
37 class AliDalitzElectronCuts : public AliAnalysisCuts {
38         
39  public: 
40
41
42   enum cutIds {
43         kgoodId=0, 
44         kededxSigmaITSCut,
45         kededxSigmaTPCCut,
46         kpidedxSigmaTPCCut,
47         kpiMinMomdedxSigmaTPCCut,
48         kpiMaxMomdedxSigmaTPCCut,
49         kLowPRejectionSigmaCut,
50         kTOFelectronPID,
51         kclsTPCCut,
52         ketaCut,
53         kPsiPair,
54         kRejectSharedElecGamma,
55         kBackgroundScheme,
56         kNumberOfRotations,
57         kNCuts
58   };
59
60
61  enum electronCuts {
62       kElectronIn=0,
63       kNoTracks,
64       kTrackCuts,
65       kdEdxCuts,
66       kElectronOut
67  };
68
69
70   Bool_t SetCutIds(TString cutString); 
71   Int_t fCuts[kNCuts];
72   Bool_t SetCut(cutIds cutID, Int_t cut);
73   Bool_t UpdateCutString(cutIds cutID, Int_t value);
74   static const char * fgkCutNames[kNCuts];
75
76
77   Bool_t InitializeCutsFromCutString(const TString analysisCutSelection); 
78   
79
80   AliDalitzElectronCuts(const char *name="ElectronCuts", const char * title="Electron Cuts");
81   virtual ~AliDalitzElectronCuts();                            //virtual destructor
82
83   virtual Bool_t IsSelected(TObject* /*obj*/){return kTRUE;}
84   virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
85
86   TString GetCutNumber();
87
88     // Cut Selection
89   Bool_t TrackIsSelected(AliESDtrack* lTrack);
90   Bool_t ElectronIsSelected(AliESDtrack* lTrack);
91   void InitAODpidUtil(Int_t type);
92   static AliDalitzElectronCuts * GetStandardCuts2010PbPb();
93   static AliDalitzElectronCuts * GetStandardCuts2010pp();
94   Bool_t InitPIDResponse();
95   
96   void SetPIDResponse(AliPIDResponse * pidResponse) {fPIDResponse = pidResponse;}
97   AliPIDResponse * GetPIDResponse() { return fPIDResponse;}
98   
99   void PrintCuts();
100
101   void InitCutHistograms(TString name="",Bool_t preCut = kTRUE,TString cutName="");
102   void SetFillCutHistograms(TString name="",Bool_t preCut = kTRUE,TString cutName=""){if(!fHistograms){InitCutHistograms(name,preCut,cutName);};}
103   TList *GetCutHistograms(){return fHistograms;}
104
105   static AliVTrack * GetTrack(AliVEvent * event, Int_t label);
106
107   
108
109   ///Cut functions
110   Bool_t dEdxCuts(AliVTrack * track);
111   Bool_t PIDProbabilityCut(AliConversionPhotonBase *photon, AliVEvent * event);
112   Bool_t RejectSharedElecGamma(TList *photons, Int_t indexEle);
113   Bool_t IsFromGammaConversion( Double_t psiPair, Double_t deltaPhi );
114
115   // Event Cuts
116
117   //Double_t GetPsiPair( const AliESDtrack *trackPos, const AliESDtrack *trackNeg );
118
119   Bool_t SetTPCdEdxCutPionLine(Int_t pidedxSigmaCut);
120   Bool_t SetTPCdEdxCutElectronLine(Int_t ededxSigmaCut);
121   Bool_t SetITSdEdxCutElectronLine(Int_t ededxSigmaCut);
122   Bool_t SetMinMomPiondEdxTPCCut(Int_t piMomdedxSigmaCut);
123   Bool_t SetMaxMomPiondEdxTPCCut(Int_t piMomdedxSigmaCut);
124   Bool_t SetTPCClusterCut(Int_t clsTPCCut);
125   Bool_t SetEtaCut(Int_t etaCut);
126   Bool_t SetMinMomPiondEdxCut(Int_t piMinMomdedxSigmaCut);
127   Bool_t SetMaxMomPiondEdxCut(Int_t piMaxMomdedxSigmaCut);
128   Bool_t SetLowPRejectionCuts(Int_t LowPRejectionSigmaCut);
129   Bool_t SetTOFElectronPIDCut(Int_t TOFelectronPID);
130   Bool_t SetPsiPairCut(Int_t psiCut);
131   Bool_t SetRejectSharedElecGamma(Int_t RCut);
132   Bool_t SetBackgroundScheme(Int_t BackgroundScheme);
133   Bool_t SetNumberOfRotations(Int_t NumberOfRotations);
134   
135   // Request Flags
136
137   Double_t GetEtaCut(){ return  fEtaCut;}
138   Double_t GetPsiPairCut(){ return fPsiPairCut; }
139   Double_t DoRejectSharedElecGamma(){ return fDoRejectSharedElecGamma;}
140   Double_t DoPsiPairCut(){return fDoPsiPairCut;}
141   Bool_t   UseTrackMultiplicity(){ return fUseTrackMultiplicityForBG;}
142   Int_t    GetBKGMethod(){ return fBKGMethod; }
143   Int_t    NumberOfRotationEvents(){return fnumberOfRotationEventsForBG;}
144   
145
146   
147   protected:
148
149   TList *fHistograms;
150   AliPIDResponse *fPIDResponse;
151   AliESDtrackCuts *fesdTrackCuts;
152
153   Double_t fEtaCut; //eta cutç
154   Double_t fRadiusCut; // radius cut
155   Double_t fPsiPairCut;
156   Double_t fDeltaPhiCutMin;
157   Double_t fDeltaPhiCutMax;
158   Double_t fMinClsTPC; // minimum clusters in the TPC
159   Double_t fMinClsTPCToF; // minimum clusters to findable clusters
160   Bool_t   fDodEdxSigmaITSCut; // flag to use the dEdxCut ITS based on sigmas
161   Bool_t   fDodEdxSigmaTPCCut; // flag to use the dEdxCut TPC based on sigmas
162   Bool_t   fDoTOFsigmaCut; // flag to use TOF pid cut RRnewTOF
163   Bool_t   fDoRejectSharedElecGamma; //Reject electrons from the gammas with Radius < RadiusCut
164   Bool_t   fDoPsiPairCut; // PsiPair Cut
165   Double_t fPIDnSigmaAboveElectronLineITS; // sigma cut
166   Double_t fPIDnSigmaBelowElectronLineITS; // sigma cut
167   Double_t fPIDnSigmaAboveElectronLineTPC; // sigma cut
168   Double_t fPIDnSigmaBelowElectronLineTPC; // sigma cut
169   Double_t fPIDnSigmaAbovePionLineTPC;
170   Double_t fPIDnSigmaAbovePionLineTPCHighPt;
171   Double_t fTofPIDnSigmaAboveElectronLine; // sigma cut RRnewTOF
172   Double_t fTofPIDnSigmaBelowElectronLine; // sigma cut RRnewTOF 
173   Double_t fPIDMinPnSigmaAbovePionLineTPC; // sigma cut
174   Double_t fPIDMaxPnSigmaAbovePionLineTPC; // sigma cut
175   Double_t fDoKaonRejectionLowP;   // Kaon rejection at low p
176   Double_t fDoProtonRejectionLowP; // Proton rejection at low p
177   Double_t fDoPionRejectionLowP;   // Pion rejection at low p
178   Double_t fPIDnSigmaAtLowPAroundKaonLine; // sigma cut
179   Double_t fPIDnSigmaAtLowPAroundProtonLine; // sigma cut
180   Double_t fPIDnSigmaAtLowPAroundPionLine; // sigma cut
181   Double_t fPIDMinPKaonRejectionLowP; // Momentum limit to apply kaon rejection
182   Double_t fPIDMinPProtonRejectionLowP; // Momentum limit to apply proton rejection
183   Double_t fPIDMinPPionRejectionLowP; // Momentum limit to apply proton rejection
184
185   Bool_t   fUseCorrectedTPCClsInfo; // flag to use corrected tpc cl info
186   Bool_t   fUseTOFpid; // flag to use tof pid
187   Bool_t   fUseTrackMultiplicityForBG; // use multiplicity
188   Int_t    fBKGMethod;
189   Int_t    fnumberOfRotationEventsForBG;
190
191
192   // Histograms
193   TObjString *fCutString; // cut number used for analysis
194   TH1F *hCutIndex; // bookkeeping for cuts
195   TH1F *hdEdxCuts;  // bookkeeping for dEdx cuts
196   TH2F *hITSdEdxbefore; // ITS dEdx before cuts
197   TH2F *hITSdEdxafter;
198   TH2F *hTPCdEdxbefore; // TPC dEdx before cuts
199   TH2F *hTPCdEdxafter; // TPC dEdx after cuts
200   TH2F *hTPCdEdxSignalafter; //TPC dEdx signal
201   TH2F *hTOFbefore; // TOF after cuts
202   TH2F *hTOFafter; // TOF after cuts
203   
204
205 private:
206
207   AliDalitzElectronCuts(const AliDalitzElectronCuts&); // not implemented
208   AliDalitzElectronCuts& operator=(const AliDalitzElectronCuts&); // not implemented
209
210
211   ClassDef(AliDalitzElectronCuts,2)
212 };
213
214
215 inline void AliDalitzElectronCuts::InitAODpidUtil(Int_t type) {
216   if (!fPIDResponse) fPIDResponse = new AliAODpidUtil();
217   Double_t alephParameters[5];
218   // simulation
219   alephParameters[0] = 2.15898e+00/50.;
220   alephParameters[1] = 1.75295e+01;
221   alephParameters[2] = 3.40030e-09;
222   alephParameters[3] = 1.96178e+00;
223   alephParameters[4] = 3.91720e+00;
224   fPIDResponse->GetTOFResponse().SetTimeResolution(80.);
225   
226   // data
227   if (type==1){
228     alephParameters[0] = 0.0283086/0.97;
229     alephParameters[1] = 2.63394e+01;
230     alephParameters[2] = 5.04114e-11;
231     alephParameters[3] = 2.12543e+00;
232     alephParameters[4] = 4.88663e+00;
233     fPIDResponse->GetTOFResponse().SetTimeResolution(130.);
234     fPIDResponse->GetTPCResponse().SetMip(50.);
235   }
236   
237   fPIDResponse->GetTPCResponse().SetBetheBlochParameters(
238     alephParameters[0],alephParameters[1],alephParameters[2],
239     alephParameters[3],alephParameters[4]);
240   
241   fPIDResponse->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);
242 }
243
244
245 #endif