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