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