]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliConversionCuts.h
Transition PWG4 --> PWGGA
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliConversionCuts.h
1 #ifndef ALICONVERSIONCUTS_H
2 #define ALICONVERSIONCUTS_H
3
4 // Class handling all kinds of selection cuts for Gamma Conversion analysis
5 // Authors: Svein Lindal, Daniel Lohner                                                                                         *
6
7 #include "AliAODpidUtil.h"
8 #include "AliConversionPhotonBase.h"
9 #include "AliAODConversionMother.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
17 class AliESDEvent;
18 class AliAODEvent;
19 class AliConversionPhotonBase;
20 class AliKFVertex;
21 class TH1F;
22 class TH2F;
23 class AliPIDResponse;
24 class AliAnalysisCuts;
25 class iostream;
26 class TList;
27 class AliAnalysisManager;
28
29
30 using namespace std;
31
32 class AliConversionCuts : public AliAnalysisCuts {
33         
34  public: 
35
36
37   enum cutIds {
38         kgoodId=0, 
39         kv0FinderType,                
40         keProbCut,                    
41         kededxSigmaCut,               
42         kpidedxSigmaCut,              
43         kpiMomdedxSigmaCut,           
44         kchi2GammaCut,                
45         ksinglePtCut,                 
46         kclsTPCCut,                   
47         ketaCut,                      
48         kchi2MesonCut,                
49         kLowPRejectionSigmaCut,       
50         kQtMaxCut,                    
51         kpiMaxMomdedxSigmaCut,        
52         kalphaMesonCut,               
53         kminRCut,                     
54         kRapidityMesonCut,            
55         kBackgroundScheme,            
56         kDegreesForRotationMethod,    
57         kNumberOfRotations,           
58         kremovePileUp,                
59         kselectV0AND,                 
60         kmultiplicityBin,             
61         kisHeavyIon,                  
62         kuseCentrality,               
63         kcentralityBin,               
64         kTOFelectronPID,              
65         kuseMCPSmearing,              
66         kdoPhotonAsymmetryCut,
67         kPsiPair, 
68         kCosPAngle,
69      //   kHBTmultiplicityBin,
70      //   kprimaryCutNumber,
71      //   kuseBayesPID,
72      //   kdalitzelectronsPID,
73      //   kpsiCutNumber,
74      //   kdalitzBackgroundType,
75      //   keleclsTPCCut,
76         kNCuts
77   };
78
79   Bool_t SetCutIds(TString cutString); 
80   Int_t fCuts[kNCuts];
81   Bool_t SetCut(cutIds cutID, Int_t cut);
82   Bool_t UpdateCutString(cutIds cutID, Int_t value);
83
84
85   static const char * fgkCutNames[kNCuts];
86
87   Double_t GetPsiPair(const AliESDv0* v0, const AliESDEvent * event) const;
88   Double_t GetCosineOfPointingAngle(const AliConversionPhotonBase * photon, const AliVEvent * event) const; 
89
90
91   Bool_t InitializeCutsFromCutString(const TString analysisCutSelection);
92
93   AliConversionCuts(const char *name="V0Cuts", const char * title="V0 Cuts");
94   virtual ~AliConversionCuts();                            //virtual destructor
95
96   virtual Bool_t IsSelected(TObject* /*obj*/){return kTRUE;}
97   virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
98
99   TString GetCutNumber();
100
101   // Cut Selection
102   Bool_t EventIsSelected(AliVEvent *fInputEvent);
103   Bool_t PhotonIsSelected(AliConversionPhotonBase * photon, AliVEvent  * event,Bool_t DoOnlyPhotonCuts=kFALSE);
104   Bool_t PhotonIsSelectedMC(TParticle *particle,AliStack *fMCStack);
105   Bool_t TracksAreSelected(AliVTrack * negTrack, AliVTrack * posTrack);
106   Bool_t MesonIsSelected(AliAODConversionMother *pi0,Bool_t IsSignal=kTRUE);
107   Bool_t MesonIsSelectedMC(TParticle *fMCMother,AliStack *fMCStack,Bool_t bMCDaughtersInAcceptance=kFALSE);
108
109   void InitAODpidUtil(Int_t type);
110   static AliConversionCuts * GetStandardCuts2010PbPb();
111   Bool_t InitPIDResponse();
112   
113   void SetPIDResponse(AliPIDResponse * pidResponse) {fPIDResponse = pidResponse;}
114   AliPIDResponse * GetPIDResponse() { return fPIDResponse;}
115   
116   void PrintCuts();
117
118   void InitCutHistograms();
119   void SetFillCutHistograms(){if(!fHistograms){InitCutHistograms();};}
120   TList *GetCutHistograms(){return fHistograms;}
121
122   AliVTrack * GetTrack(AliVEvent * event, Int_t label) const;
123
124   ///Cut functions
125   Bool_t SpecificTrackCuts(AliAODTrack * negTrack, AliAODTrack * posTrack,Int_t &cutIndex);
126   Bool_t SpecificTrackCuts(AliESDtrack * negTrack, AliESDtrack * posTrack,Int_t &cutIndex);
127   Bool_t AcceptanceCuts(AliConversionPhotonBase *photon);
128   Bool_t AcceptanceCut(TParticle *particle, TParticle * ePos,TParticle* eNeg);
129   Bool_t dEdxCuts(AliVTrack * track);
130   Bool_t ArmenterosQtCut(AliConversionPhotonBase *photon);
131   Bool_t AsymmetryCut(AliConversionPhotonBase *photon,AliVEvent *event);
132   Bool_t PIDProbabilityCut(AliConversionPhotonBase *photon, AliVEvent * event);
133   Bool_t SelectV0Finder(Bool_t onfly){return onfly&&fUseOnFlyV0Finder;}
134   Bool_t PhotonCuts(AliConversionPhotonBase *photon,AliVEvent *event);
135   Bool_t CorrectedTPCClusterCut(AliConversionPhotonBase *photon, AliVEvent * event);
136   Bool_t PsiPairCut(const AliConversionPhotonBase * photon, const AliVEvent * event) const;
137   Bool_t CosinePAngleCut(const AliConversionPhotonBase * photon, const AliVEvent * event) const;
138   // Event Cuts
139   Bool_t IsCentralitySelected(AliVEvent *fInputEvent);
140   Double_t GetCentrality(AliVEvent *event);
141   Int_t GetNumberOfContributorsVtx(AliVEvent *event);
142   Bool_t VertexZCut(AliVEvent *fInputEvent);
143
144   // Set Individual Cuts
145   Bool_t SetRCut(Int_t RCut);
146   Bool_t SetV0Finder(Int_t v0FinderType);
147   Bool_t SetElectronProbCut(Int_t eProbCut);
148   Bool_t SetChi2GammaCut(Int_t chi2GammaCut);
149   Bool_t SetTPCdEdxCutPionLine(Int_t pidedxSigmaCut);
150   Bool_t SetTPCdEdxCutElectronLine(Int_t ededxSigmaCut);
151   Bool_t SetSinglePtCut(Int_t singlePtCut);
152   Bool_t SetTPCClusterCut(Int_t clsTPCCut);
153   Bool_t SetEtaCut(Int_t etaCut);
154   Bool_t SetChi2MesonCut(Int_t chi2MesonCut);
155   Bool_t SetMinMomPiondEdxCut(Int_t piMinMomdedxSigmaCut);
156   Bool_t SetMaxMomPiondEdxCut(Int_t piMaxMomdedxSigmaCut);
157   Bool_t SetLowPRejectionCuts(Int_t LowPRejectionSigmaCut);
158   Bool_t SetQtMaxCut(Int_t QtMaxCut);
159   Bool_t SetAlphaMesonCut(Int_t alphaMesonCut);
160   Bool_t SetTOFElectronPIDCut(Int_t TOFelectronPID);
161   Bool_t SetTRDElectronCut(Int_t TRDElectronCut);
162   Bool_t SetRapidityMesonCut(Int_t RapidityMesonCut);
163   Bool_t SetUseCentrality(Int_t useCentrality);
164   Bool_t SetIsHeavyIon(Int_t isHeavyIon);
165   Bool_t SetCentralityBin(Int_t centralityBin);
166   Bool_t SetPhotonAsymmetryCut(Int_t doPhotonAsymmetryCut);
167   Bool_t SetRemovePileUp(Int_t removePileUp);
168   Bool_t SetBackgroundScheme(Int_t BackgroundScheme);
169   Bool_t SetNDegreesForRotationMethod(Int_t DegreesForRotationMethod);
170   Bool_t SetNumberOfRotations(Int_t NumberOfRotations);
171   Bool_t SetMCPSmearing(Int_t useMCPSmearing);
172   Bool_t SetMultiplicityBin(Int_t multiplicityBin);
173   Bool_t SetSelectV0AND(Int_t selectV0AND);
174   Bool_t SetCosPAngleCut(Int_t cosCut);
175   Bool_t SetPsiPairCut(Int_t psiCut);
176
177   // Request Flags
178
179   Bool_t UseRotationMethod(){return fUseRotationMethodInBG;}
180   Int_t NumberOfRotationEvents(){return fnumberOfRotationEventsForBG;}
181   Int_t NDegreesRotation(){return fnDegreeRotationPMForBG;}
182   Bool_t IsHeavyIon(){return fIsHeavyIon;}
183   Bool_t UseTrackMultiplicity(){return fUseTrackMultiplicityForBG;}
184
185
186   Int_t GetFirstTPCRow(Double_t radius);
187   
188
189   protected:
190   TList *fHistograms;
191   AliPIDResponse *fPIDResponse;
192
193   //cuts
194   Double_t fMaxR; //r cut
195   Double_t fMinR; //r cut
196   Double_t fEtaCut; //eta cut
197   Double_t fEtaCutMin; //eta cut
198   Double_t fPtCut; // pt cut
199   Double_t fSinglePtCut; // pt cut for electron/positron
200   Double_t fMaxZ; //z cut
201   Double_t fMinClsTPC; // minimum clusters in the TPC
202   Double_t fMinClsTPCToF; // minimum clusters to findable clusters
203   Double_t fLineCutZRSlope; //linecut
204   Double_t fLineCutZValue; //linecut
205   Double_t fLineCutZRSlopeMin; //linecut
206   Double_t fLineCutZValueMin; //linecut
207   Double_t fChi2CutConversion; //chi2cut
208   Double_t fChi2CutMeson; //chicut meson
209   Double_t fPIDProbabilityCutNegativeParticle;
210   Double_t fPIDProbabilityCutPositiveParticle;
211   Bool_t   fDodEdxSigmaCut; // flag to use the dEdxCut based on sigmas
212   Bool_t   fDoTOFsigmaCut; // flag to use TOF pid cut RRnewTOF
213   Double_t fPIDTRDEfficiency; // required electron efficiency for TRD PID
214   Bool_t   fDoTRDPID; // flag to use TRD pid
215   Double_t fPIDnSigmaAboveElectronLine; // sigma cut
216   Double_t fPIDnSigmaBelowElectronLine; // sigma cut
217   Double_t fTofPIDnSigmaAboveElectronLine; // sigma cut RRnewTOF
218   Double_t fTofPIDnSigmaBelowElectronLine; // sigma cut RRnewTOF 
219   Double_t fPIDnSigmaAbovePionLine;     // sigma cut
220   Double_t fPIDnSigmaAbovePionLineHighPt;     // sigma cut
221   Double_t fPIDMinPnSigmaAbovePionLine; // sigma cut
222   Double_t fPIDMaxPnSigmaAbovePionLine; // sigma cut
223   Double_t fDoKaonRejectionLowP;   // Kaon rejection at low p
224   Double_t fDoProtonRejectionLowP; // Proton rejection at low p
225   Double_t fDoPionRejectionLowP;   // Pion rejection at low p
226   Double_t fPIDnSigmaAtLowPAroundKaonLine; // sigma cut
227   Double_t fPIDnSigmaAtLowPAroundProtonLine; // sigma cut
228   Double_t fPIDnSigmaAtLowPAroundPionLine; // sigma cut
229   Double_t fPIDMinPKaonRejectionLowP; // Momentum limit to apply kaon rejection
230   Double_t fPIDMinPProtonRejectionLowP; // Momentum limit to apply proton rejection
231   Double_t fPIDMinPPionRejectionLowP; // Momentum limit to apply proton rejection
232   Bool_t   fDoQtGammaSelection; // Select gammas using qtMax
233   Bool_t   fDoHighPtQtGammaSelection; // RRnew Select gammas using qtMax for high pT
234   Double_t fQtMax; // Maximum Qt from Armenteros to select Gammas
235   Double_t fHighPtQtMax; // RRnew Maximum Qt for High pT from Armenteros to select Gammas
236   Double_t fPtBorderForQt; // RRnew 
237   Double_t fXVertexCut; //vertex cut
238   Double_t fYVertexCut; //vertex cut
239   Double_t fZVertexCut; // vertexcut
240   Double_t fNSigmaMass; //nsigma cut
241   Bool_t fUseEtaMinCut; //flag
242   Bool_t fUseOnFlyV0Finder; //flag
243   Bool_t   fDoPhotonAsymmetryCut; // flag to use the PhotonAsymetryCut
244   Double_t fMinPPhotonAsymmetryCut; // Min Momentum for Asymmetry Cut
245   Double_t fMinPhotonAsymmetry;  // Asymmetry Cut
246   Bool_t fIsHeavyIon;               // flag for heavy ion
247   Double_t fMaxVertexZ;    // max z offset of vertex
248   Bool_t fUseCentrality;  // centrality selection
249   Bool_t fUseCentralityBin; // centrality selection with individual bins
250   Bool_t fUseCorrectedTPCClsInfo; // flag to use corrected tpc cl info
251   Bool_t fUseTOFpid; // flag to use tof pid
252   Double_t fAlphaMinCutMeson; // min value for meson alpha cut
253   Double_t fAlphaCutMeson; // max value for meson alpha cut
254   Double_t fRapidityCutMeson; // max value for meson rapidity
255   Bool_t fUseRotationMethodInBG; // flag to apply rotation method for meson bg estimation
256   Bool_t fdoBGProbability; // flag to use probability method for meson bg estimation
257   Bool_t fUseTrackMultiplicityForBG; // flag to use track multiplicity for meson bg estimation (else V0 mult)
258   Int_t fnDegreeRotationPMForBG; //
259   Int_t fnumberOfRotationEventsForBG; //
260   Bool_t fUseMCPSmearing; // flag
261   Double_t fPBremSmearing;//
262   Double_t fPSigSmearing; //
263   Double_t fPSigSmearingCte; //
264   Bool_t fUseMultiplicity; // flag
265   Int_t fUseMultiplicityBin; // selected multiplicity bin
266   Bool_t fSelectV0AND; // flag
267   Bool_t fRemovePileUp; //flag
268   Float_t fOpeningAngle; // min opening angle for meson
269   Float_t fPsiPairCut;
270   Float_t fCosPAngleCut;
271
272   // Histograms
273   TObjString *fCutString; // cut number used for analysis
274   TH1F *hdEdxCuts;  // bookkeeping for dEdx cuts
275   TH2F *hTPCdEdxbefore; // TPC dEdx before cuts
276   TH2F *hTPCdEdxafter; // TPC dEdx after cuts
277   TH1F *hTrackCuts; // bookkeeping for track cuts
278   TH1F *hPhotonCuts; // bookkeeping for photon specific cuts
279   TH1F *hInvMassbefore; // e+e- inv mass distribution before cuts
280   TH2F *hArmenterosbefore; // armenteros podolanski plot before cuts
281   TH1F *hInvMassafter; // e+e- inv mass distribution after cuts
282   TH2F *hArmenterosafter;  // armenteros podolanski plot after cuts
283   TH1F *hAcceptanceCuts; // bookkeeping for acceptance cuts
284   TH1F *hCutIndex; // bookkeeping for cuts
285   TH1F *hV0EventCuts; // bookkeeping for event selection cuts
286   TH1F *hCentrality; // centrality distribution for selected events
287   TH1F *hVertexZ; // vertex z distribution for selected events
288   TH1F *hMesonCuts; // bookkeeping for meson cuts
289   TH1F *hMesonBGCuts; // bookkeeping for meson bg cuts
290
291
292 private:
293
294   AliConversionCuts(const AliConversionCuts&); // not implemented
295   AliConversionCuts& operator=(const AliConversionCuts&); // not implemented
296
297
298   ClassDef(AliConversionCuts,2)
299 };
300
301
302 inline void AliConversionCuts::InitAODpidUtil(Int_t type) {
303   if (!fPIDResponse) fPIDResponse = new AliAODpidUtil();
304   Double_t alephParameters[5];
305   // simulation
306   alephParameters[0] = 2.15898e+00/50.;
307   alephParameters[1] = 1.75295e+01;
308   alephParameters[2] = 3.40030e-09;
309   alephParameters[3] = 1.96178e+00;
310   alephParameters[4] = 3.91720e+00;
311   fPIDResponse->GetTOFResponse().SetTimeResolution(80.);
312   
313   // data
314   if (type==1){
315     alephParameters[0] = 0.0283086/0.97;
316     alephParameters[1] = 2.63394e+01;
317     alephParameters[2] = 5.04114e-11;
318     alephParameters[3] = 2.12543e+00;
319     alephParameters[4] = 4.88663e+00;
320     fPIDResponse->GetTOFResponse().SetTimeResolution(130.);
321     fPIDResponse->GetTPCResponse().SetMip(50.);
322   }
323   
324   fPIDResponse->GetTPCResponse().SetBetheBlochParameters(
325     alephParameters[0],alephParameters[1],alephParameters[2],
326     alephParameters[3],alephParameters[4]);
327   
328   fPIDResponse->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);
329 }
330
331
332 #endif