]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/TriggerPID/AliTwoParticlePIDCorr.h
3836d6074fa5bdd24ce06d7206c80b3bc38c525b
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / TriggerPID / AliTwoParticlePIDCorr.h
1 #ifndef ALITWOPARTICLEPIDCORR_H
2 #define ALITWOPARTICLEPIDCORR_H
3
4 #include "THn.h" // in cxx file causes .../THn.h:257: error: conflicting declaration ‘typedef class THnT<float> THnF’
5
6
7 class TH1F;
8 class TH2F;
9 class TH3F;
10 class THnSparse;
11 class TString;
12 class TList;
13 //class AliESDtrackCuts;
14 class TSeqCollection;
15 class AliPIDResponse;
16 class AliPIDCombined;  
17 class AliAODEvent;
18 class AliAODTrack;
19 class AliAODVertex;
20 class AliEventPoolManager;
21 class TFormula;
22 class AliAnalysisUtils;
23 class LRCParticlePID;
24 class AliVParticle;
25 class AliCFContainer;
26 class AliCFGridSparse;
27 class THnBase;
28 class AliTHn;
29
30
31 #include <TObject.h> //LRCParticlePID is a derived class from"TObject"
32 #include "TMath.h"
33 #include "TNamed.h"
34 #include "AliUEHist.h"
35 #include "AliPID.h"
36 #include "AliAnalysisTask.h"
37 #include "AliUEHist.h"
38 #include "TString.h"
39 #include "AliVParticle.h"
40 #include "TParticle.h"
41 #include "AliLog.h"
42 #include "AliTHn.h"
43
44
45
46 #ifndef ALIANALYSISTASKSE_H
47 #include "AliAnalysisTaskSE.h"
48 #endif
49
50 namespace AliPIDNameSpace {
51   
52   enum PIDType
53   {
54     NSigmaTPC = 0,
55     NSigmaTOF,
56     NSigmaTPCTOF,  // squared sum
57     Bayes    
58   };
59   const Int_t NSigmaPIDType=NSigmaTPCTOF;//number of Nsigma PID types
60     
61   enum AliDetectorType
62   {
63     fITS = 0,
64     fTPC,
65     fTOF,
66     fNDetectors
67   };
68   
69   
70   enum AliParticleSpecies
71   {
72     SpPion = 0,
73     SpKaon,
74     SpProton,
75     unidentified,
76     NSpecies=unidentified,
77     SpUndefined=999
78   }; // Particle species used in plotting
79   
80   
81   enum AliCharge
82   {
83     Posch = 0,
84     Negch,
85     NCharge
86   };
87 }
88
89
90 using namespace AliPIDNameSpace;
91
92 class AliTwoParticlePIDCorr : public AliAnalysisTaskSE {
93  public:
94     AliTwoParticlePIDCorr();
95     AliTwoParticlePIDCorr(const char *name);
96     virtual ~AliTwoParticlePIDCorr();
97     
98     virtual void     UserCreateOutputObjects();
99     virtual void     UserExec(Option_t *option);
100     virtual void    doAODevent();
101     virtual void    doMCAODevent();
102     virtual void     Terminate(Option_t *);
103   void     SetSharedClusterCut(Double_t value) { fSharedClusterCut = value; }
104
105   void SettwoTrackEfficiencyCutDataReco(Bool_t twoTrackEfficiencyCutDataReco,Float_t twoTrackEfficiencyCutValue1)
106   {
107     ftwoTrackEfficiencyCutDataReco=twoTrackEfficiencyCutDataReco;
108     twoTrackEfficiencyCutValue=twoTrackEfficiencyCutValue1;
109   }
110   void SetVertextype(Int_t Vertextype){fVertextype=Vertextype;}                                                 //Check it every time
111     void SetZvtxcut(Double_t zvtxcut) {fzvtxcut=zvtxcut;}
112     void SetCustomBinning(TString receivedCustomBinning) { fCustomBinning = receivedCustomBinning; }
113     void SetMaxNofMixingTracks(Int_t MaxNofMixingTracks) {fMaxNofMixingTracks=MaxNofMixingTracks;}               //Check it every time
114   void SetCentralityEstimator(TString CentralityMethod) { fCentralityMethod = CentralityMethod;}
115   void SetSampleType(TString SampleType) {fSampleType=SampleType;}
116   void SetRequestEventPlane(Bool_t RequestEventPlane){fRequestEventPlane=RequestEventPlane;}
117   void SetAnalysisType(TString AnalysisType){fAnalysisType=AnalysisType;}
118   void SetFilterBit(Int_t FilterBit) {fFilterBit=FilterBit;}
119   void SetTrackStatus(UInt_t status) { fTrackStatus = status; }
120   void SetfilltrigassoUNID(Bool_t filltrigassoUNID){ffilltrigassoUNID=filltrigassoUNID;}
121   void SetfilltrigUNIDassoID(Bool_t filltrigUNIDassoID){ffilltrigUNIDassoID=filltrigUNIDassoID;}
122   void SetfilltrigIDassoUNID(Bool_t filltrigIDassoUNID){ffilltrigIDassoUNID=filltrigIDassoUNID;}
123   void SetfilltrigIDassoID(Bool_t filltrigIDassoID){ ffilltrigIDassoID=filltrigIDassoID;}
124   void SetfilltrigIDassoIDMCTRUTH(Bool_t filltrigIDassoIDMCTRUTH){ffilltrigIDassoIDMCTRUTH=filltrigIDassoIDMCTRUTH;}
125   void SetSelectHighestPtTrig(Bool_t SelectHighestPtTrig){fSelectHighestPtTrig=SelectHighestPtTrig;}
126   void SetTriggerSpeciesSelection(Bool_t TriggerSpeciesSelection,Int_t TriggerSpecies,Bool_t containPIDtrig){
127     fTriggerSpeciesSelection=TriggerSpeciesSelection;//if it is KTRUE then Set containPIDtrig=kFALSE
128     fTriggerSpecies=TriggerSpecies;
129     fcontainPIDtrig=containPIDtrig;
130   }
131   void SetAssociatedSpeciesSelection(Bool_t AssociatedSpeciesSelection,Int_t AssociatedSpecies, Bool_t containPIDasso)
132   {
133     fAssociatedSpeciesSelection=AssociatedSpeciesSelection;//if it is KTRUE then Set containPIDasso=kFALSE
134     fAssociatedSpecies=AssociatedSpecies;
135     fcontainPIDasso=containPIDasso;
136   }
137
138   void SettingChargeCounting(Int_t val) {SetChargeAxis=val;}
139
140  void SetFIllPIDQAHistos(Bool_t FIllPIDQAHistos){fFIllPIDQAHistos=FIllPIDQAHistos;}
141   void SetRejectPileUp(Bool_t rejectPileUp) {frejectPileUp=rejectPileUp;}
142   void SetCheckFirstEventInChunk(Bool_t CheckFirstEventInChunk) {fCheckFirstEventInChunk=CheckFirstEventInChunk;}
143
144   void SetKinematicCuts(Float_t minPt, Float_t maxPt,Float_t mineta,Float_t maxeta)
145   {
146     fminPt=minPt;
147     fmaxPt=maxPt;
148     fmineta=mineta;
149     fmaxeta=maxeta;
150   }
151   void SetDcaCut(Bool_t dcacut,Double_t dcacutvalue)
152   {
153     fdcacut=dcacut;
154     fdcacutvalue=dcacutvalue;
155   }
156   void SetfillHistQA(Bool_t fillhistQAReco,Bool_t fillhistQATruth)
157   {
158     ffillhistQAReco=fillhistQAReco;
159     ffillhistQATruth=fillhistQATruth;
160   }
161   void SetPtordering(Bool_t PtOrderDataReco,Bool_t PtOrderMCTruth)
162   {
163     fPtOrderDataReco=PtOrderDataReco;
164     fPtOrderMCTruth=PtOrderMCTruth;
165   }
166   void SetWeightPerEvent(Bool_t flag) {  fWeightPerEvent = flag;}
167   void Setselectprimarydatareco(Bool_t onlyprimarydatareco) {fonlyprimarydatareco=onlyprimarydatareco;}
168   void SetselectprimaryTruth(Bool_t selectprimaryTruth) {fselectprimaryTruth=selectprimaryTruth;}
169   void SetCombinedNSigmaCut(Double_t NSigmaPID) {fNSigmaPID=NSigmaPID;}
170   void SetHighPtKaonNSigmaPID(Float_t HighPtKaonSigma,Float_t HighPtKaonNSigmaPID)
171   {
172     fHighPtKaonSigma=HighPtKaonSigma;
173     fHighPtKaonNSigmaPID=HighPtKaonNSigmaPID;
174   }
175   void IgnoreoverlappedTracks(Bool_t UseExclusiveNSigma){fUseExclusiveNSigma=UseExclusiveNSigma;}
176   void SetRemoveTracksT0Fill( Bool_t RemoveTracksT0Fill){fRemoveTracksT0Fill=RemoveTracksT0Fill;}
177   void SetPairSelectCharge(Int_t SelectCharge){fSelectCharge=SelectCharge;}
178   void SetTrigAssoSelectcharge( Int_t TriggerSelectCharge,Int_t AssociatedSelectCharge)
179   {
180     fTriggerSelectCharge=TriggerSelectCharge;
181     fAssociatedSelectCharge=AssociatedSelectCharge;
182   }
183   void SetEtaOrdering(Bool_t EtaOrdering){fEtaOrdering=EtaOrdering;}
184   void SetCutConversionsResonances( Bool_t CutConversions,Bool_t CutResonances)
185   {
186      fCutConversions=CutConversions;
187      fCutResonances=CutResonances;
188   }
189   void SetCleanUp(Bool_t InjectedSignals,Bool_t RemoveWeakDecays,Bool_t RemoveDuplicates)
190   {
191     fInjectedSignals=InjectedSignals;
192     fRemoveWeakDecays=RemoveWeakDecays;
193     fRemoveDuplicates=RemoveDuplicates;
194   }
195   void SetEfficiency(Bool_t fillefficiency,Bool_t applyTrigefficiency,Bool_t applyAssoefficiency)
196     {
197       ffillefficiency=fillefficiency;
198       fapplyTrigefficiency=applyTrigefficiency;
199       fapplyAssoefficiency=applyAssoefficiency;
200
201     }
202   void SetComboeffPtRange(Double_t minPtComboeff,Double_t maxPtComboeff) {
203     fminPtComboeff=minPtComboeff;
204     fmaxPtComboeff=maxPtComboeff;}
205   //only one can be kTRUE at a time(for the next two Setters)
206   void Setmesoneffrequired(Bool_t mesoneffrequired) {fmesoneffrequired=mesoneffrequired;}
207   void Setkaonprotoneffrequired(Bool_t kaonprotoneffrequired){fkaonprotoneffrequired=kaonprotoneffrequired;}
208    void SetOnlyOneEtaSide(Int_t OnlyOneEtaSide){fOnlyOneEtaSide=OnlyOneEtaSide;}
209    void SetRejectResonanceDaughters(Int_t RejectResonanceDaughters){fRejectResonanceDaughters=RejectResonanceDaughters;}
210
211 void SetTOFPIDVal(Bool_t RequestTOFPID,Float_t PtTOFPIDmin,Float_t PtTOFPIDmax)
212    {
213 fRequestTOFPID=RequestTOFPID;
214 fPtTOFPIDmin=PtTOFPIDmin;
215 fPtTOFPIDmax=PtTOFPIDmax;
216 }
217
218  void SetEffcorectionfilePathName(TString efffilename) {fefffilename=efffilename;}
219  
220  private:
221  //histograms
222     TList *fOutput;        //! Output list
223     TList *fOutputList;        //! Output list
224
225     TString    fCentralityMethod;     // Method to determine centrality
226     TString    fSampleType;     // pp,p-Pb,Pb-Pb
227     Bool_t fRequestEventPlane; //only for PbPb
228     Int_t    fnTracksVertex;        // QA tracks pointing to principal vertex
229     AliAODVertex* trkVtx;//!
230     Float_t zvtx;
231     Int_t    fFilterBit;         // track selection cuts
232      UInt_t         fTrackStatus;       // if non-0, the bits set in this variable are required for each track
233     Double_t       fSharedClusterCut;  // cut on shared clusters (only for AOD)
234     Int_t fVertextype;
235     Int_t skipParticlesAbove;
236     Double_t fzvtxcut;
237     Bool_t ffilltrigassoUNID;
238     Bool_t ffilltrigUNIDassoID;
239     Bool_t ffilltrigIDassoUNID;
240     Bool_t ffilltrigIDassoID;
241     Bool_t ffilltrigIDassoIDMCTRUTH;
242     Int_t fMaxNofMixingTracks;
243     Bool_t fPtOrderMCTruth;
244     Bool_t fPtOrderDataReco;
245     Bool_t fWeightPerEvent;
246     Bool_t fTriggerSpeciesSelection;
247     Bool_t fAssociatedSpeciesSelection;
248     Bool_t fRandomizeReactionPlane;
249     Int_t fTriggerSpecies;
250     Int_t fAssociatedSpecies;
251     TString fCustomBinning;//for setting customized binning
252     TString fBinningString;//final binning string
253     Bool_t fSelectHighestPtTrig;
254     Bool_t fcontainPIDtrig;
255     Bool_t fcontainPIDasso;
256     Int_t SetChargeAxis;
257      Bool_t frejectPileUp;
258      Bool_t fCheckFirstEventInChunk;
259     Float_t fminPt;
260     Float_t fmaxPt;
261     Float_t fmineta;
262     Float_t fmaxeta;
263     Bool_t fselectprimaryTruth;
264     Bool_t fonlyprimarydatareco;
265     Bool_t fdcacut;
266     Double_t fdcacutvalue;
267     Bool_t ffillhistQAReco;
268     Bool_t ffillhistQATruth;
269     Int_t kTrackVariablesPair ;
270     Double_t fminPtTrig;
271     Double_t fmaxPtTrig;
272     Double_t fminPtComboeff;
273     Double_t fmaxPtComboeff;
274     Double_t fminPtAsso;
275     Double_t fmaxPtAsso;
276     Double_t fmincentmult;
277     Double_t fmaxcentmult;
278     TH1F *fhistcentrality;//!
279     TH1F *fEventCounter; //!
280     TH2F *fEtaSpectrasso;//!
281     TH2F *fphiSpectraasso;//!
282     TH1F *MCtruthpt;//! 
283     TH1F *MCtrutheta;//! 
284     TH1F *MCtruthphi;//!
285     TH1F *MCtruthpionpt;//!
286     TH1F *MCtruthpioneta;//!
287     TH1F *MCtruthpionphi;//!
288     TH1F *MCtruthkaonpt;//!
289     TH1F *MCtruthkaoneta;//!
290     TH1F *MCtruthkaonphi;//!
291     TH1F *MCtruthprotonpt;//!
292     TH1F *MCtruthprotoneta;//!
293     TH1F *MCtruthprotonphi;//! 
294     TH2F *fPioncont;//!
295     TH2F *fKaoncont;//!
296     TH2F *fProtoncont;//!
297     TH2F *fUNIDcont;//!
298     TH2F *fEventno;//!
299     TH2F *fEventnobaryon;//!
300     TH2F *fEventnomeson;//!
301     TH2F *fhistJetTrigestimate;//!
302
303     TH2D* fCentralityCorrelation;  //! centrality vs Tracks multiplicity
304  //VZERO calibration
305   TH1F *fHistVZEROAGainEqualizationMap;//VZERO calibration map
306   TH1F *fHistVZEROCGainEqualizationMap;//VZERO calibration map
307   TH2F *fHistVZEROChannelGainEqualizationMap; //VZERO calibration map
308   TH1* fCentralityWeights;                   // for centrality flattening
309
310     TH2F *fHistCentStats; //!centrality stats
311     TH2F *fHistRefmult;//!
312     TH2F *fHistEQVZEROvsTPCmultiplicity;//!
313     TH2F *fHistEQVZEROAvsTPCmultiplicity;//!
314     TH2F *fHistEQVZEROCvsTPCmultiplicity;//!
315     TH2F *fHistVZEROCvsEQVZEROCmultiplicity;//!
316     TH2F *fHistVZEROAvsEQVZEROAmultiplicity;//!
317     TH2F *fHistVZEROCvsVZEROAmultiplicity;//!
318     TH2F *fHistEQVZEROCvsEQVZEROAmultiplicity;//!
319     TH2F *fHistVZEROSignal;//!
320     TH2F *fHistEventPlaneReco;//!
321     TH2F *fHistEventPlaneTruth;//!
322     TH2D *fHistPsiMinusPhi;//! psi - phi QA histogram
323
324     TH2F* fControlConvResoncances; //! control histograms for cuts on conversions and resonances
325
326     TH2F *fHistoTPCdEdx;//!
327     TH2F *fHistoTOFbeta;//!
328     TH3F *fTPCTOFPion3d;//!
329     TH3F *fTPCTOFKaon3d;//!
330     TH3F *fTPCTOFProton3d;//!
331     TH1F *fPionPt;//!
332     TH1F *fPionEta;//!
333     TH1F *fPionPhi;//!
334     TH1F *fKaonPt;//!
335     TH1F *fKaonEta;//!
336     TH1F *fKaonPhi;//!
337     TH1F *fProtonPt;//!
338     TH1F *fProtonEta;//!
339     TH1F *fProtonPhi;//!
340     // TH3F *fHistocentNSigmaTPC;//! nsigma TPC
341     // TH3F *fHistocentNSigmaTOF;//! nsigma TOF 
342     
343     AliTHn *fCorrelatonTruthPrimary;//!
344     AliTHn *fCorrelatonTruthPrimarymix;//!
345     AliTHn *fTHnCorrUNID;//!
346     AliTHn *fTHnCorrUNIDmix;//!
347     AliTHn *fTHnCorrID;//!
348     AliTHn *fTHnCorrIDmix;//!
349     AliTHn *fTHnCorrIDUNID;//!
350     AliTHn *fTHnCorrIDUNIDmix;//!
351     AliTHn *fTHnTrigcount;//!
352     AliTHn *fTHnTrigcountMCTruthPrim;//!
353     AliTHn* fTrackHistEfficiency[6]; //! container for tracking efficiency and contamination (all particles filled including leading one): axes: eta, pT, particle species:::::::::0 pion, 1 kaon,2 proton,3 mesons,4 kaons+protons,5 all
354
355     
356     TH1F *fHistQA[16]; //!                  gReactionPlane
357      
358    
359     THnSparse *effcorection[6];//!
360     // THnF *effmap[6];  
361
362     Int_t ClassifyTrack(AliAODTrack* track,AliAODVertex* vertex,Float_t magfield,Bool_t fill);
363   Double_t* GetBinning(const char* configuration, const char* tag, Int_t& nBins);
364
365
366   void Fillcorrelation(Double_t gReactionPlane,TObjArray *trackstrig,TObjArray *tracksasso,Double_t cent,Float_t vtx,Float_t weight,Bool_t firstTime,Float_t bSign,Bool_t fPtOrder,Bool_t twoTrackEfficiencyCut,Bool_t mixcase,TString fillup);//mixcase=kTRUE in case of mixing; 
367  Float_t GetTrackbyTrackeffvalue(AliAODTrack* track,Double_t cent,Float_t evzvtx, Int_t parpid);
368
369 //Mixing functions
370   void DefineEventPool();
371   AliEventPoolManager    *fPoolMgr;//! 
372   TClonesArray          *fArrayMC;//!
373   TString          fAnalysisType;          // "MC", "ESD", "AOD"
374   TString fefffilename;
375
376     //PID part histograms
377
378   //PID functions
379     Bool_t HasTPCPID(AliAODTrack *track) const; // has TPC PID
380     Bool_t HasTOFPID(AliAODTrack *track) const; // has TOF PID
381     Double_t GetBeta(AliAODTrack *track);
382     void CalculateNSigmas(AliAODTrack *track, Bool_t FIllQAHistos);
383     Int_t FindMinNSigma(AliAODTrack *track, Bool_t FIllQAHistos);
384     Bool_t* GetDoubleCounting(AliAODTrack * trk, Bool_t FIllQAHistos);
385     Int_t GetParticle(AliAODTrack * trk, Bool_t FIllQAHistos);  
386  
387      TH2F* GetHistogram2D(const char * name);//return histogram "name" from fOutputList
388
389      Bool_t ftwoTrackEfficiencyCutDataReco;        
390    Float_t twoTrackEfficiencyCutValue;
391   //Pid objects
392   AliPIDResponse *fPID; //! PID
393   AliPIDCombined   *fPIDCombined;     //! PIDCombined
394
395   //set PID Combined
396   void SetPIDCombined(AliPIDCombined *obj){fPIDCombined=obj;}
397   AliPIDCombined *GetPIDCombined(){return fPIDCombined;}
398   //set cut on beyesian probability
399   void SetBayesCut(Double_t cut){fBayesCut=cut;}
400   void SetdiffPIDcutvalues(Bool_t diffPIDcutvalues,Double_t PIDCutval1, Double_t PIDCutval2, Double_t PIDCutval3,Double_t PIDCutval4){
401     fdiffPIDcutvalues=diffPIDcutvalues;
402     fPIDCutval1=PIDCutval1;
403     fPIDCutval2=PIDCutval2;
404     fPIDCutval3=PIDCutval3;
405     fPIDCutval4=PIDCutval4;
406   }
407  void  SetRandomizeReactionPlane(Bool_t RandomizeReactionPlane){fRandomizeReactionPlane=RandomizeReactionPlane;}
408   Double_t GetBayesCut(){return fBayesCut;}
409   Int_t GetIDBayes(AliAODTrack *trk, Bool_t FIllQAHistos);//calculate the PID according to bayesian PID
410   UInt_t CalcPIDCombined(AliAODTrack *track,Int_t detMask, Double_t* prob) const;
411   Bool_t* GetAllCompatibleIdentitiesNSigma(AliAODTrack * trk, Bool_t FIllQAHistos);//All the identities are true
412
413
414   Int_t eventno;
415   Float_t fPtTOFPIDmin; //lower pt bound for the TOCTOF combined circular pid
416   Float_t fPtTOFPIDmax; //uper pt bound for the TOCTOF combined circular pid
417   Bool_t fRequestTOFPID;//if true returns kSpUndefined if the TOF signal is missing
418   PIDType fPIDType; // PID type  Double_t fNSigmaPID; // number of sigma for PID cut
419   Bool_t fFIllPIDQAHistos; //Switch for filling the nSigma histos
420   Double_t fNSigmaPID; // number of sigma for PID cut
421   Double_t fBayesCut; // Cut on Bayesian probability
422  Bool_t fdiffPIDcutvalues;
423  Double_t fPIDCutval1;
424  Double_t fPIDCutval2;
425  Double_t fPIDCutval3;
426  Double_t fPIDCutval4;
427
428   Float_t fHighPtKaonNSigmaPID;// number of sigma for PID cut for Kaons above fHighPtKaonSigma(-1 default, no cut applied)
429   Float_t fHighPtKaonSigma;//lower pt bound for the fHighPtKaonNSigmaPID to be set >0(i.e. to make it applicable)
430   Bool_t fUseExclusiveNSigma;//if true returns the identity only if no double counting(i.e not in the overlap area)
431   Bool_t fRemoveTracksT0Fill;//if true remove tracks for which only StartTime from To-Fill is available (worst resolution)
432  Int_t fSelectCharge;           // (un)like sign selection when building correlations: 0: no selection; 1: unlike sign; 2: like sign
433     Int_t fTriggerSelectCharge;    // select charge of trigger particle: 1: positive; -1 negative
434     Int_t fAssociatedSelectCharge; // select charge of associated particle: 1: positive; -1 negative
435     Float_t fTriggerRestrictEta;   // restrict eta range for trigger particle (default: -1 [off])
436     Bool_t fEtaOrdering;           // eta ordering, see AliUEHistograms.h for documentation
437     Bool_t fCutConversions;        // cut on conversions (inv mass)
438     Bool_t fCutResonances;         // cut on resonances (inv mass)
439     Int_t fRejectResonanceDaughters; // reject all daughters of all resonance candidates (1: test method (cut at m_inv=0.9); 2: k0; 3: lambda)
440     Int_t fOnlyOneEtaSide;       // decides that only trigger particle from one eta side are considered (0 = all; -1 = negative, 1 = positive)
441     Bool_t fInjectedSignals;      // check header to skip injected signals in MC
442     Bool_t fRemoveWeakDecays;      // remove secondaries from weak decays from tracks and particles
443     Bool_t fRemoveDuplicates;// remove particles with the same label (double reconstruction)
444     Bool_t fapplyTrigefficiency;//if kTRUE then eff correction calculation starts
445     Bool_t fapplyAssoefficiency;//if kTRUE then eff correction calculation starts
446     Bool_t ffillefficiency;  //if kTRUE then THNsparses used for eff. calculation are filled up
447     Bool_t fmesoneffrequired;
448     Bool_t fkaonprotoneffrequired;
449    AliAnalysisUtils*     fAnalysisUtils;      // points to class with common analysis utilities
450   TFormula*      fDCAXYCut;          // additional pt dependent cut on DCA XY (only for AOD)
451
452
453   Float_t fnsigmas[NSpecies][NSigmaPIDType+1]; //nsigma values
454   Bool_t fHasDoubleCounting[NSpecies];//array with compatible identities
455
456   //Int_t fPIDMethod; // PID method
457
458  //functions
459   Bool_t CheckTrack(AliAODTrack * part);
460   Float_t PhiRange(Float_t DPhi);
461   Float_t GetInvMassSquared(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2);
462 Float_t GetInvMassSquaredCheap(Float_t pt1, Float_t eta1, Float_t phi1, Float_t pt2, Float_t eta2, Float_t phi2, Float_t m0_1, Float_t m0_2);
463   Float_t GetDPhiStar(Float_t phi1, Float_t pt1, Float_t charge1, Float_t phi2, Float_t pt2, Float_t charge2, Float_t radius, Float_t bSign);
464   TObjArray* CloneAndReduceTrackList(TObjArray* tracks);
465
466
467   void ShiftTracks(TObjArray* tracks, Double_t angle);
468   Bool_t AcceptEventCentralityWeight(Double_t centrality);
469
470   //get event plane
471   Double_t GetEventPlane(AliAODEvent *event,Bool_t truth);
472   Double_t GetAcceptedEventMultiplicity(AliAODEvent *aod,Bool_t truth);//returns centrality after event(mainly vertex) selection IsEventAccepted  GetAcceptedEventMultiplicity
473   
474   //get vzero equalization
475   Double_t GetEqualizationFactor(Int_t run, const char* side);
476   Double_t GetChannelEqualizationFactor(Int_t run,Int_t channel);
477   void SetVZEROCalibrationFile(const char* filename,const char* lhcPeriod);
478   void SetCentralityWeights(TH1* hist) { fCentralityWeights = hist; }
479
480   Double_t GetRefMultiOrCentrality(AliAODEvent *event, Bool_t truth);
481   Double_t GetReferenceMultiplicityVZEROFromAOD(AliAODEvent *event);//mainly important for pp 7 TeV
482
483
484
485
486     
487     AliTwoParticlePIDCorr(const AliTwoParticlePIDCorr&); // not implemented
488     AliTwoParticlePIDCorr& operator=(const AliTwoParticlePIDCorr&); // not implemented
489     
490     ClassDef(AliTwoParticlePIDCorr, 1); // example of analysis
491 };
492 class LRCParticlePID : public TObject {
493 public:
494  LRCParticlePID(Int_t par,Short_t icharge,Float_t pt,Float_t eta, Float_t phi,Float_t effcorrectionval)
495    :fparticle(par),fcharge(icharge),fPt(pt), fEta(eta), fPhi(phi),feffcorrectionval(effcorrectionval)  {}
496   virtual ~LRCParticlePID() {}
497
498   
499     virtual Float_t Eta()        const { return fEta; }
500     virtual Float_t Phi()        const { return fPhi; }
501     virtual Float_t Pt() const { return fPt; }
502     Int_t getparticle() const {return fparticle;}
503     virtual Short_t Charge()      const { return fcharge; }
504     Float_t geteffcorrectionval() const {return feffcorrectionval;}
505     virtual Bool_t IsEqual(const TObject* obj) const { return (obj->GetUniqueID() == GetUniqueID()); }
506     virtual void SetPhi(Double_t phiv) { fPhi = phiv; }
507
508
509 private:
510   LRCParticlePID(const LRCParticlePID&);  // not implemented
511    LRCParticlePID& operator=(const LRCParticlePID&);  // not implemented
512   
513   Int_t fparticle;
514   Short_t fcharge;
515   Float_t fPt;
516   Float_t fEta;
517   Float_t fPhi;
518   Float_t feffcorrectionval;
519   ClassDef(LRCParticlePID, 1);
520 } ;
521
522 #endif
523