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