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