]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliConversionCuts.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[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 #include "TF1.h"
18 #include "AliAnalysisUtils.h"
19 #include "AliAnalysisManager.h"
20
21 class AliESDEvent;
22 class AliAODEvent;
23 class AliConversionPhotonBase;
24 class AliKFVertex;
25 class TH1F;
26 class TH2F;
27 class TF1;
28 class AliPIDResponse;
29 class AliAnalysisCuts;
30 class iostream;
31 class TList;
32 class AliAnalysisManager;
33 class AliAODMCParticle;
34
35 using namespace std;
36
37 class AliConversionCuts : public AliAnalysisCuts {
38       
39    public: 
40       
41
42    enum cutIds {
43       kisHeavyIon,                  
44       kCentralityMin,               
45       kCentralityMax,               
46       kselectV0AND,                 
47       kmultiplicityMethod,             
48       kremovePileUp,                
49       kExtraSignals, 
50       kv0FinderType,                
51       ketaCut,                                     
52       kRCut,                     
53       ksinglePtCut,                 
54       kclsTPCCut,                   
55       kededxSigmaCut,               
56       kpidedxSigmaCut,              
57       kpiMomdedxSigmaCut,        
58       kpiMaxMomdedxSigmaCut,        
59       kLowPRejectionSigmaCut,       
60       kTOFelectronPID,              
61       kQtMaxCut,                    
62       kchi2GammaCut,                
63       kPsiPair, 
64       kdoPhotonAsymmetryCut,
65       kCosPAngle,
66       kElecShare,
67       kToCloseV0s,
68       kDcaRPrimVtx,
69       kDcaZPrimVtx,
70       kInPlaneOutOfPlane,
71       kNCuts
72    };
73
74    enum photonCuts {
75          kPhotonIn=0,
76          kOnFly,
77          kNoTracks,
78          kTrackCuts,
79          kdEdxCuts,
80          kConvPointFail,
81          kPhotonCuts,
82          kEventPlane,
83          kPhotonOut
84    };
85
86
87    Bool_t SetCutIds(TString cutString); 
88    Int_t fCuts[kNCuts];
89    Bool_t SetCut(cutIds cutID, Int_t cut);
90    Bool_t UpdateCutString();
91
92
93    static const char * fgkCutNames[kNCuts];
94
95    Double_t GetCosineOfPointingAngle(const AliConversionPhotonBase * photon, AliVEvent * event) const; 
96
97
98    Bool_t InitializeCutsFromCutString(const TString analysisCutSelection);
99    void SelectCollisionCandidates(UInt_t offlineTriggerMask = AliVEvent::kAny) {
100       fOfflineTriggerMask = offlineTriggerMask;
101       fTriggerSelectedManually = kTRUE;
102    }
103    void SelectSpecialTrigger(UInt_t offlineTriggerMask = AliVEvent::kAny, TString TriggerClassName = "AliVEvent::kAny" ) {
104       fOfflineTriggerMask = offlineTriggerMask;
105       fSpecialTriggerName = TriggerClassName;
106       cout << fSpecialTriggerName.Data() << endl;
107       
108    }   
109    void FillElectonLabelArray(AliAODConversionPhoton* photon, Int_t nV0);
110    void SetAcceptedHeader(TList *HeaderList){fHeaderList = HeaderList;}   
111    void SetPreSelectionCutFlag(Bool_t preSelFlag){fPreSelCut = preSelFlag;}   
112    TString *GetFoundHeader(){return fGeneratorNames;}
113
114    Int_t GetEventQuality(){return fEventQuality;}
115    Bool_t GetIsFromPileup(){return fRemovePileUp;}
116    
117    AliConversionCuts(const char *name="V0Cuts", const char * title="V0 Cuts");
118    AliConversionCuts(const AliConversionCuts&);
119    AliConversionCuts& operator=(const AliConversionCuts&);
120
121    virtual ~AliConversionCuts();                            //virtual destructor
122
123    static AliConversionCuts * GetStandardCuts2010PbPb();
124    static AliConversionCuts * GetStandardCuts2010pp();
125
126    virtual Bool_t IsSelected(TObject* /*obj*/){return kTRUE;}
127    virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
128
129    TString GetCutNumber();
130    
131    void GetCentralityRange(Double_t range[2]){range[0]=10*fCentralityMin;range[1]=10*fCentralityMax;};
132    
133    // Cut Selection
134    Bool_t EventIsSelected(AliVEvent *fInputEvent, AliVEvent *fMCEvent);
135    Int_t IsEventAcceptedByConversionCut(AliConversionCuts *ReaderCuts, AliVEvent *InputEvent, AliMCEvent *MCEvent, Bool_t isHeavyIon);
136    Bool_t PhotonIsSelected(AliConversionPhotonBase * photon, AliVEvent  * event);
137    Bool_t PhotonIsSelectedMC(TParticle *particle,AliStack *fMCStack,Bool_t checkForConvertedGamma=kTRUE);
138    Bool_t PhotonIsSelectedAODMC(AliAODMCParticle *particle,TClonesArray *aodmcArray,Bool_t checkForConvertedGamma=kTRUE);
139    Bool_t ElectronIsSelectedMC(TParticle *particle,AliStack *fMCStack);
140    Bool_t TracksAreSelected(AliVTrack * negTrack, AliVTrack * posTrack);
141    Bool_t MesonIsSelected(AliAODConversionMother *pi0,Bool_t IsSignal=kTRUE);
142    Bool_t MesonIsSelectedMC(TParticle *fMCMother,AliStack *fMCStack, Bool_t bMCDaughtersInAcceptance=kFALSE);
143
144    void InitAODpidUtil(Int_t type);
145    Bool_t InitPIDResponse();
146    
147    void SetPIDResponse(AliPIDResponse * pidResponse) {fPIDResponse = pidResponse;}
148    AliPIDResponse * GetPIDResponse() { return fPIDResponse;}
149    
150    void PrintCuts();
151    void PrintCutsWithValues();
152
153    void InitCutHistograms(TString name="",Bool_t preCut = kTRUE);
154    void SetFillCutHistograms(TString name="",Bool_t preCut = kTRUE){if(!fHistograms){InitCutHistograms(name,preCut);};}
155    TList *GetCutHistograms(){return fHistograms;}
156    void FillPhotonCutIndex(Int_t photoncut){if(hCutIndex)hCutIndex->Fill(photoncut);}
157    void SetEtaShift(Double_t etaShift) {
158       fEtaShift = etaShift;
159       fLineCutZRSlope = tan(2*atan(exp(-fEtaCut + etaShift)));
160       if(fEtaCutMin > -0.1)
161          fLineCutZRSlopeMin = tan(2*atan(exp(-fEtaCutMin + etaShift)));
162    }
163    void SetEtaShift(TString pPbOrPbp) {
164       Double_t etaShift = 0.0;
165       if(!pPbOrPbp.CompareTo("pPb"))      etaShift = -0.465;
166       else if(!pPbOrPbp.CompareTo("Pbp")) etaShift =  0.465;
167       
168       fEtaShift = etaShift;
169       fLineCutZRSlope = tan(2*atan(exp(-fEtaCut + etaShift)));
170       if(fEtaCutMin > -0.1)
171          fLineCutZRSlopeMin = tan(2*atan(exp(-fEtaCutMin + etaShift)));
172    }
173    Double_t GetEtaShift() {return fEtaShift;}
174    Bool_t GetDoEtaShift(){return fDoEtaShift;}
175    void DoEtaShift(Bool_t doEtaShift){fDoEtaShift = doEtaShift;}
176    void GetCorrectEtaShiftFromPeriod(TString periodName);
177    
178    static AliVTrack * GetTrack(AliVEvent * event, Int_t label);
179    static AliESDtrack *GetESDTrack(AliESDEvent * event, Int_t label);
180    
181    ///Cut functions
182    Bool_t SpecificTrackCuts(AliAODTrack * negTrack, AliAODTrack * posTrack,Int_t &cutIndex);
183    Bool_t SpecificTrackCuts(AliESDtrack * negTrack, AliESDtrack * posTrack,Int_t &cutIndex);
184    Bool_t AcceptanceCuts(AliConversionPhotonBase *photon);
185    Bool_t AcceptanceCut(TParticle *particle, TParticle * ePos,TParticle* eNeg);
186    Bool_t dEdxCuts(AliVTrack * track);
187    Bool_t ArmenterosQtCut(AliConversionPhotonBase *photon);
188    Bool_t AsymmetryCut(AliConversionPhotonBase *photon,AliVEvent *event);
189    Bool_t PIDProbabilityCut(AliConversionPhotonBase *photon, AliVEvent * event);
190    Bool_t SelectV0Finder(Bool_t onfly){
191       if(onfly == fUseOnFlyV0Finder) return kTRUE;
192       else return kFALSE;
193    }
194    Bool_t PhotonCuts(AliConversionPhotonBase *photon,AliVEvent *event);
195    Bool_t CorrectedTPCClusterCut(AliConversionPhotonBase *photon, AliVEvent * event);
196    Bool_t PsiPairCut(const AliConversionPhotonBase * photon) const;
197    Bool_t CosinePAngleCut(const AliConversionPhotonBase * photon, AliVEvent * event) const;
198    Bool_t RejectSharedElectronV0s(AliAODConversionPhoton* photon, Int_t nV0, Int_t nV0s);
199    Bool_t RejectToCloseV0s(AliAODConversionPhoton* photon, TList *photons, Int_t nV0);
200    Int_t IsParticleFromBGEvent(Int_t index, AliStack *MCStack, AliVEvent *InputEvent = 0x0);
201    void GetNotRejectedParticles(Int_t rejection, TList *HeaderList, AliVEvent *MCEvent);
202    void SetUseReweightingWithHistogramFromFile( Bool_t pi0reweight=kTRUE, Bool_t etareweight=kFALSE, Bool_t k0sreweight=kFALSE, TString path="$ALICE_ROOT/PWGGA/GammaConv/MCSpectraInput.root", 
203                                                 TString histoNamePi0 = "", TString histoNameEta = "", TString histoNameK0s = "",
204                                                 TString fitNamePi0 = "", TString fitNameEta = "", TString fitNameK0s ="" ) {
205       AliInfo(Form("enabled reweighting for: pi0 : %i, eta: %i, K0s: %i",pi0reweight, etareweight, k0sreweight));
206       fDoReweightHistoMCPi0 = pi0reweight; 
207       fDoReweightHistoMCEta = etareweight; 
208       fDoReweightHistoMCK0s = k0sreweight; 
209       fPathTrFReweighting=path;
210       fNameHistoReweightingPi0 =histoNamePi0;
211       fNameHistoReweightingEta =histoNameEta;
212       fNameHistoReweightingK0s =histoNameK0s; 
213       fNameFitDataPi0 =fitNamePi0;
214       fNameFitDataEta =fitNameEta;
215       fNameFitDataK0s =fitNameK0s; 
216       
217    }
218    void  LoadReweightingHistosMCFromFile ();
219    // Event Cuts
220    Bool_t IsCentralitySelected(AliVEvent *fInputEvent, AliVEvent *fMCEvent = NULL);
221    Double_t GetCentrality(AliVEvent *event);
222    Int_t GetNumberOfContributorsVtx(AliVEvent *event);
223    Bool_t VertexZCut(AliVEvent *fInputEvent);
224    Bool_t IsTriggerSelected(AliVEvent *fInputEvent);
225    Bool_t HasV0AND(){return fHasV0AND;}
226    Bool_t IsSDDFired(){return fIsSDDFired;}
227    Int_t IsSpecialTrigger(){return fSpecialTrigger;}
228    TString GetSpecialTriggerName(){return fSpecialTriggerName;}
229    Bool_t InPlaneOutOfPlaneCut(Double_t photonPhi, Double_t eventPlaneAngle = -100, Bool_t fill = kTRUE);
230    Int_t GetInPlaneOutOfPlaneCut(){return fInPlaneOutOfPlane;}
231
232
233    // Set Individual Cuts
234    Bool_t SetRCut(Int_t RCut);
235    Bool_t SetV0Finder(Int_t v0FinderType);
236    Bool_t SetChi2GammaCut(Int_t chi2GammaCut);
237    Bool_t SetTPCdEdxCutPionLine(Int_t pidedxSigmaCut);
238    Bool_t SetTPCdEdxCutElectronLine(Int_t ededxSigmaCut);
239    Bool_t SetSinglePtCut(Int_t singlePtCut);
240    Bool_t SetTPCClusterCut(Int_t clsTPCCut);
241    Bool_t SetEtaCut(Int_t etaCut);
242    Bool_t SetMinMomPiondEdxCut(Int_t piMinMomdedxSigmaCut);
243    Bool_t SetMaxMomPiondEdxCut(Int_t piMaxMomdedxSigmaCut);
244    Bool_t SetLowPRejectionCuts(Int_t LowPRejectionSigmaCut);
245    Bool_t SetQtMaxCut(Int_t QtMaxCut);
246    Bool_t SetTOFElectronPIDCut(Int_t TOFelectronPID);
247    Bool_t SetTRDElectronCut(Int_t TRDElectronCut);
248    Bool_t SetCentralityMin(Int_t useCentrality);
249    Bool_t SetIsHeavyIon(Int_t isHeavyIon);
250    Bool_t SetCentralityMax(Int_t centralityBin);
251    Bool_t SetPhotonAsymmetryCut(Int_t doPhotonAsymmetryCut);
252    Bool_t SetRemovePileUp(Int_t removePileUp);  
253    Bool_t SetMultiplicityMethod(Int_t multiplicityMethod);
254    Int_t SetSelectSpecialTrigger(Int_t selectSpecialTrigger);
255    Bool_t SetCosPAngleCut(Int_t cosCut);
256    Bool_t SetPsiPairCut(Int_t psiCut);
257    Bool_t SetSharedElectronCut(Int_t sharedElec);
258    Bool_t SetToCloseV0sCut(Int_t toClose);
259    Bool_t SetRejectExtraSignalsCut(Int_t extraSignal);
260    Bool_t SetDCARPhotonPrimVtxCut(Int_t DCARPhotonPrimVtx);
261    Bool_t SetDCAZPhotonPrimVtxCut(Int_t DCAZPhotonPrimVtx);
262    Bool_t SetInPlaneOutOfPlane(Int_t inOutPlane);
263    
264    // Request Flags
265
266    Int_t IsHeavyIon(){return fIsHeavyIon;}
267    Int_t GetFirstTPCRow(Double_t radius);
268    Float_t GetWeightForMeson(TString period, Int_t index, AliStack *MCStack, AliVEvent *InputEvent = 0x0);
269
270    Bool_t UseElecSharingCut(){return fDoSharedElecCut;}
271    Bool_t UseToCloseV0sCut(){return fDoToCloseV0sCut;}
272    Int_t GetMultiplicityMethod(){return fMultiplicityMethod;}
273    Double_t GetEtaCut(){return fEtaCut;}
274    Int_t GetSignalRejection(){return fRejectExtraSignals;}
275    Int_t GetNAcceptedHeaders(){return fnHeaders; }
276    TString * GetAcceptedHeaderNames(){return fGeneratorNames;}
277    Int_t * GetAcceptedHeaderStart(){return fNotRejectedStart;}
278    Int_t * GetAcceptedHeaderEnd(){return fNotRejectedEnd;}
279    TList* GetAcceptedHeader(){return fHeaderList;}
280    
281    
282    protected:
283    TList *fHistograms;
284    TList *fHeaderList;
285    AliPIDResponse *fPIDResponse;
286
287
288    Int_t fEventQuality; // EventQuality
289    //cuts
290    Double_t fMaxR; //r cut
291    Double_t fMinR; //r cut
292    Double_t fEtaCut; //eta cut
293    Double_t fEtaCutMin; //eta cut
294    Double_t fPtCut; // pt cut
295    Double_t fSinglePtCut; // pt cut for electron/positron
296    Double_t fMaxZ; //z cut
297    Double_t fMinClsTPC; // minimum clusters in the TPC
298    Double_t fMinClsTPCToF; // minimum clusters to findable clusters
299    Double_t fLineCutZRSlope; //linecut
300    Double_t fLineCutZValue; //linecut
301    Double_t fLineCutZRSlopeMin; //linecut
302    Double_t fLineCutZValueMin; //linecut
303    Double_t fChi2CutConversion; //chi2cut
304    Double_t fPIDProbabilityCutNegativeParticle;
305    Double_t fPIDProbabilityCutPositiveParticle;
306    Bool_t   fDodEdxSigmaCut; // flag to use the dEdxCut based on sigmas
307    Bool_t   fDoTOFsigmaCut; // flag to use TOF pid cut RRnewTOF
308    Double_t fPIDTRDEfficiency; // required electron efficiency for TRD PID
309    Bool_t   fDoTRDPID; // flag to use TRD pid
310    Double_t fPIDnSigmaAboveElectronLine; // sigma cut
311    Double_t fPIDnSigmaBelowElectronLine; // sigma cut
312    Double_t fTofPIDnSigmaAboveElectronLine; // sigma cut RRnewTOF
313    Double_t fTofPIDnSigmaBelowElectronLine; // sigma cut RRnewTOF 
314    Double_t fPIDnSigmaAbovePionLine;     // sigma cut
315    Double_t fPIDnSigmaAbovePionLineHighPt;     // sigma cut
316    Double_t fPIDMinPnSigmaAbovePionLine; // sigma cut
317    Double_t fPIDMaxPnSigmaAbovePionLine; // sigma cut
318    Double_t fDoKaonRejectionLowP;   // Kaon rejection at low p
319    Double_t fDoProtonRejectionLowP; // Proton rejection at low p
320    Double_t fDoPionRejectionLowP;   // Pion rejection at low p
321    Double_t fPIDnSigmaAtLowPAroundKaonLine; // sigma cut
322    Double_t fPIDnSigmaAtLowPAroundProtonLine; // sigma cut
323    Double_t fPIDnSigmaAtLowPAroundPionLine; // sigma cut
324    Double_t fPIDMinPKaonRejectionLowP; // Momentum limit to apply kaon rejection
325    Double_t fPIDMinPProtonRejectionLowP; // Momentum limit to apply proton rejection
326    Double_t fPIDMinPPionRejectionLowP; // Momentum limit to apply proton rejection
327    Bool_t   fDoQtGammaSelection; // Select gammas using qtMax
328    Bool_t   fDo2DQt; // Select gammas using ellipse cut
329    Double_t fQtMax; // Maximum Qt from Armenteros to select Gammas
330    Double_t fXVertexCut; //vertex cut
331    Double_t fYVertexCut; //vertex cut
332    Double_t fZVertexCut; // vertexcut
333    Double_t fNSigmaMass; //nsigma cut
334    Bool_t fUseEtaMinCut; //flag
335    Bool_t fUseOnFlyV0Finder; //flag
336    Bool_t   fDoPhotonAsymmetryCut; // flag to use the PhotonAsymetryCut
337    Double_t fMinPPhotonAsymmetryCut; // Min Momentum for Asymmetry Cut
338    Double_t fMinPhotonAsymmetry;  // Asymmetry Cut
339    Int_t  fIsHeavyIon;               // flag for heavy ion
340    Int_t fDetectorCentrality;  // centrality detecotor V0M or CL1
341    Int_t fModCentralityClass; // allows to select smaller centrality classes
342    Double_t fMaxVertexZ;    // max z offset of vertex
343    Int_t fCentralityMin;  // centrality selection lower bin value
344    Int_t fCentralityMax; // centrality selection upper bin value
345    Bool_t fUseCorrectedTPCClsInfo; // flag to use corrected tpc cl info
346    Bool_t fUseTOFpid; // flag to use tof pid
347    Int_t fMultiplicityMethod; // selected multiplicity method
348    Int_t fSpecialTrigger; // flag
349    Bool_t fRemovePileUp; //flag
350    Float_t fOpeningAngle; // min opening angle for meson
351    Float_t fPsiPairCut;
352    Bool_t fDo2DPsiPairChi2;
353    Float_t fCosPAngleCut;
354    Bool_t fDoToCloseV0sCut; //
355    Int_t fRejectExtraSignals;//
356    Double_t fminV0Dist; //
357    Bool_t fDoSharedElecCut; //
358    UInt_t fOfflineTriggerMask;   //  Task processes collision candidates only
359    Bool_t fHasV0AND; // V0AND Offline Trigger
360    Bool_t fIsSDDFired; // SDD FIRED to select with SDD events
361    TRandom3 fRandom; //
362    Int_t fElectronArraySize; // Size of electron array
363    Int_t *fElectronLabelArray; //[fElectronArraySize]
364    Double_t fDCAZPrimVtxCut; // cut value for the maximum distance in Z between the photon & the primary vertex [cm]
365    Double_t fDCARPrimVtxCut; // cut value for the maximum distance in R between the photon & the primary vertex [cm]
366    Int_t fInPlaneOutOfPlane; // In-Plane Out-Of Plane Analysis
367    Float_t fConversionPointXArray; // Array with conversion Point x
368    Float_t fConversionPointYArray; // Array with conversion Point y
369    Float_t fConversionPointZArray; // Array with conversion Point z
370    Int_t fnHeaders; // Number of Headers
371    Int_t *fNotRejectedStart; //[fnHeaders]
372    Int_t *fNotRejectedEnd; //[fnHeaders]
373    TString *fGeneratorNames; //[fnHeaders]
374    TObjString *fCutString; // cut number used for analysis
375    AliAnalysisUtils *fUtils;
376    Double_t fEtaShift;
377    Bool_t fDoEtaShift;            // Flag for Etashift
378    Bool_t fDoReweightHistoMCPi0; // Flag for reweighting Pi0 input with histogram
379    Bool_t fDoReweightHistoMCEta; // Flag for reweighting Eta input with histogram
380    Bool_t fDoReweightHistoMCK0s; // Flag for reweighting K0s input with histogram
381    TString fPathTrFReweighting; // Path for file used in reweighting
382    TString fNameHistoReweightingPi0; //Histogram name for reweighting Pi0
383    TString fNameHistoReweightingEta; //Histogram name for reweighting Eta
384    TString fNameHistoReweightingK0s; //Histogram name for reweighting K0s
385    TString fNameFitDataPi0; //Fit name for fit to spectrum of pi0s in Data
386    TString fNameFitDataEta; //Fit name for fit to spectrum of etas in Data
387    TString fNameFitDataK0s; //Fit name for fit to spectrum of k0s in Data
388    // Histograms
389    TH1F *hdEdxCuts;  // bookkeeping for dEdx cuts
390    TH2F *hTPCdEdxbefore; // TPC dEdx before cuts
391    TH2F *hTPCdEdxafter; // TPC dEdx after cuts
392    TH2F *hTPCdEdxSigbefore; // TPC Sigma dEdx before cuts
393    TH2F *hTPCdEdxSigafter; // TPC Sigm dEdx after cuts
394    TH2F *hTOFbefore; // TOF before cuts
395    TH2F *hTOFSigbefore; // TOF Sigma before cuts
396    TH2F *hTOFSigafter; // TOF Sigma after cuts
397    TH2F *hPsiPairDeltaPhiafter; // TOF Sigma after cuts
398    TH1F *hTrackCuts; // bookkeeping for track cuts
399    TH1F *hPhotonCuts; // bookkeeping for photon specific cuts
400    TH1F *hInvMassbefore; // e+e- inv mass distribution before cuts
401    TH2F *hArmenterosbefore; // armenteros podolanski plot before cuts
402    TH1F *hInvMassafter; // e+e- inv mass distribution after cuts
403    TH2F *hArmenterosafter;  // armenteros podolanski plot after cuts
404    TH1F *hAcceptanceCuts; // bookkeeping for acceptance cuts
405    TH1F *hCutIndex; // bookkeeping for cuts
406    TH1F *hV0EventCuts; // bookkeeping for event selection cuts
407    TH1F *hCentrality; // centrality distribution for selected events
408    TH2F *hCentralityVsNumberOfPrimaryTracks; // centrality distribution for selected events
409    TH1F *hVertexZ; // vertex z distribution for selected events
410    TH1F *hEventPlanePhi; //EventPlaneAngle Minus Photon Angle
411    TH1F *hTriggerClass; //fired offline trigger class
412    TH1F *hTriggerClassSelected; //selected fired offline trigger class
413    TH1D *hReweightMCHistPi0; //histogram input for reweighting Pi0
414    TH1D *hReweightMCHistEta; //histogram input for reweighting Eta
415    TH1D *hReweightMCHistK0s; //histogram input for reweighting K0s
416    TF1  *fFitDataPi0; //fit to pi0 spectrum in Data
417    TF1  *fFitDataEta; //fit to eta spectrum in Data
418    TF1  *fFitDataK0s; //fit to K0s spectrum in Data
419    Bool_t fPreSelCut; // Flag for preselection cut used in V0Reader
420    Bool_t fTriggerSelectedManually; // Flag for manual trigger selection
421    TString fSpecialTriggerName; // Name of the Special Triggers
422 private:
423
424    ClassDef(AliConversionCuts,7)
425 };
426
427
428 inline void AliConversionCuts::InitAODpidUtil(Int_t type) {
429   if (!fPIDResponse) fPIDResponse = new AliAODpidUtil();
430   Double_t alephParameters[5];
431   // simulation
432   alephParameters[0] = 2.15898e+00/50.;
433   alephParameters[1] = 1.75295e+01;
434   alephParameters[2] = 3.40030e-09;
435   alephParameters[3] = 1.96178e+00;
436   alephParameters[4] = 3.91720e+00;
437   fPIDResponse->GetTOFResponse().SetTimeResolution(80.);
438   
439   // data
440   if (type==1){
441     alephParameters[0] = 0.0283086/0.97;
442     alephParameters[1] = 2.63394e+01;
443     alephParameters[2] = 5.04114e-11;
444     alephParameters[3] = 2.12543e+00;
445     alephParameters[4] = 4.88663e+00;
446     fPIDResponse->GetTOFResponse().SetTimeResolution(130.);
447     fPIDResponse->GetTPCResponse().SetMip(50.);
448   }
449   
450   fPIDResponse->GetTPCResponse().SetBetheBlochParameters(
451     alephParameters[0],alephParameters[1],alephParameters[2],
452     alephParameters[3],alephParameters[4]);
453   
454   fPIDResponse->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);
455 }
456
457
458 #endif