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