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