]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/TriggerPID/AliTwoParticlePIDCorr.h
Update on Two Particle Pid Corr: Debojit
[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 class TProfile;
30
31
32 #include <TObject.h> //LRCParticlePID is a derived class from"TObject"
33 #include "TMath.h"
34 #include "TNamed.h"
35 #include "AliUEHist.h"
36 #include "AliPID.h"
37 #include "AliAnalysisTask.h"
38 #include "AliUEHist.h"
39 #include "TString.h"
40 #include "AliVParticle.h"
41 #include "TParticle.h"
42 #include "AliLog.h"
43 #include "AliTHn.h"
44 #include "TBits.h"
45
46
47 #ifndef ALIANALYSISTASKSE_H
48 #include "AliAnalysisTaskSE.h"
49 #endif
50
51 namespace AliPIDNameSpace {
52   
53   enum PIDType
54   {
55     NSigmaTPC = 0,
56     NSigmaTOF,
57     NSigmaTPCTOF,  // squared sum
58     Bayes    
59   };
60   const Int_t NSigmaPIDType=NSigmaTPCTOF;//number of Nsigma PID types
61     
62   enum AliDetectorType
63   {
64     fITS = 0,
65     fTPC,
66     fTOF,
67     fNDetectors
68   };
69   
70   
71   enum AliParticleSpecies
72   {
73     SpPion = 0,
74     SpKaon,
75     SpProton,
76     unidentified,
77     NSpecies=unidentified,
78     SpUndefined=999
79   }; // Particle species used in plotting
80   
81   
82   enum AliCharge
83   {
84     Posch = 0,
85     Negch,
86     NCharge
87   };
88 }
89
90
91 using namespace AliPIDNameSpace;
92
93 class AliTwoParticlePIDCorr : public AliAnalysisTaskSE {
94  public:
95     AliTwoParticlePIDCorr();
96     AliTwoParticlePIDCorr(const char *name);
97     virtual ~AliTwoParticlePIDCorr();
98     
99     virtual void     UserCreateOutputObjects();
100     virtual void     UserExec(Option_t *option);
101     virtual void    doAODevent();
102     virtual void    doMCAODevent();
103     virtual void     Terminate(Option_t *);
104   void     SetSharedClusterCut(Double_t value) { fSharedClusterCut = value; }
105   void     SetSharedTPCmapCut(Double_t value1) { fSharedTPCmapCut = value1; }
106   void     SetSharedfraction_Pair_cut(Double_t value2) { fSharedfraction_Pair_cut = value2; }
107
108
109   void SettwoTrackEfficiencyCutDataReco(Bool_t twoTrackEfficiencyCutDataReco,Float_t twoTrackEfficiencyCutValue1,Float_t TwoTrackCutMinRadius,Float_t TwoTrackCutMaxRadius)
110   {
111     ftwoTrackEfficiencyCutDataReco=twoTrackEfficiencyCutDataReco;
112     twoTrackEfficiencyCutValue=twoTrackEfficiencyCutValue1;
113     fTwoTrackCutMinRadius=TwoTrackCutMinRadius;
114     fTwoTrackCutMaxRadius=TwoTrackCutMaxRadius;
115   }
116   void SetVertextype(Int_t Vertextype){fVertextype=Vertextype;}                                                 //Check it every time
117     void SetZvtxcut(Double_t zvtxcut) {fzvtxcut=zvtxcut;}
118     void SetCustomBinning(TString receivedCustomBinning) { fCustomBinning = receivedCustomBinning; }
119     void SetMaxNofMixingTracks(Int_t MaxNofMixingTracks) {fMaxNofMixingTracks=MaxNofMixingTracks;}               //Check it every time
120   void SetCentralityEstimator(TString CentralityMethod) { fCentralityMethod = CentralityMethod;}
121   void SetPPVsMultUtils(Bool_t val)  {fPPVsMultUtils = val;}
122   void SetSampleType(TString SampleType) {fSampleType=SampleType;}
123   void SetRequestEventPlane(Bool_t RequestEventPlane,Bool_t V2,Bool_t V3,TString EPdetector,Bool_t IsAfter2011){
124 fRequestEventPlane=RequestEventPlane;
125 fV2=V2;
126 fV3=V3;
127 fEPdet=EPdetector;
128 fIsAfter2011=IsAfter2011;
129 }
130   void SetAnalysisType(TString AnalysisType){fAnalysisType=AnalysisType;}
131   void SetFilterBit(Int_t FilterBit) {fFilterBit=FilterBit;}
132   void SetTrackStatus(UInt_t status) { fTrackStatus = status; }
133   void SetfilltrigassoUNID(Bool_t filltrigassoUNID){ffilltrigassoUNID=filltrigassoUNID;}
134   void SetfilltrigUNIDassoID(Bool_t filltrigUNIDassoID){ffilltrigUNIDassoID=filltrigUNIDassoID;}
135   void SetfilltrigIDassoUNID(Bool_t filltrigIDassoUNID){ffilltrigIDassoUNID=filltrigIDassoUNID;}
136   void SetfilltrigIDassoID(Bool_t filltrigIDassoID){ ffilltrigIDassoID=filltrigIDassoID;}
137   void SetfilltrigIDassoIDMCTRUTH(Bool_t filltrigIDassoIDMCTRUTH){ffilltrigIDassoIDMCTRUTH=filltrigIDassoIDMCTRUTH;}
138   void SetSelectHighestPtTrig(Bool_t SelectHighestPtTrig){fSelectHighestPtTrig=SelectHighestPtTrig;}
139   void SetTriggerSpeciesSelection(Bool_t TriggerSpeciesSelection,Int_t TriggerSpecies,Bool_t containPIDtrig){
140     fTriggerSpeciesSelection=TriggerSpeciesSelection;//if it is KTRUE then Set containPIDtrig=kFALSE
141     fTriggerSpecies=TriggerSpecies;
142     fcontainPIDtrig=containPIDtrig;
143   }
144   void SetAssociatedSpeciesSelection(Bool_t AssociatedSpeciesSelection,Int_t AssociatedSpecies, Bool_t containPIDasso)
145   {
146     fAssociatedSpeciesSelection=AssociatedSpeciesSelection;//if it is KTRUE then Set containPIDasso=kFALSE
147     fAssociatedSpecies=AssociatedSpecies;
148     fcontainPIDasso=containPIDasso;
149   }
150
151   void SettingChargeCounting(Int_t val) {SetChargeAxis=val;}
152
153  void SetFIllPIDQAHistos(Bool_t FIllPIDQAHistos){fFIllPIDQAHistos=FIllPIDQAHistos;}
154   void SetRejectPileUp(Bool_t rejectPileUp) {frejectPileUp=rejectPileUp;}
155   void SetCheckFirstEventInChunk(Bool_t CheckFirstEventInChunk) {fCheckFirstEventInChunk=CheckFirstEventInChunk;}
156
157   void SetKinematicCuts(Float_t minPt, Float_t maxPt,Float_t mineta,Float_t maxeta)
158   {
159     fminPt=minPt;
160     fmaxPt=maxPt;
161     fmineta=mineta;
162     fmaxeta=maxeta;
163   }
164   void SetDcaCut(Bool_t dcacut,Double_t dcacutvalue)
165   {
166     fdcacut=dcacut;
167     fdcacutvalue=dcacutvalue;
168   }
169   void SetfillHistQA(Bool_t fillhistQAReco,Bool_t fillhistQATruth)
170   {
171     ffillhistQAReco=fillhistQAReco;
172     ffillhistQATruth=fillhistQATruth;
173   }
174   void SetPtordering(Bool_t PtOrderDataReco,Bool_t PtOrderMCTruth)
175   {
176     fPtOrderDataReco=PtOrderDataReco;
177     fPtOrderMCTruth=PtOrderMCTruth;
178   }
179   void SetWeightPerEvent(Bool_t flag) {  fWeightPerEvent = flag;}
180   void Setselectprimarydatareco(Bool_t onlyprimarydatareco) {fonlyprimarydatareco=onlyprimarydatareco;}
181   void SetselectprimaryTruth(Bool_t selectprimaryTruth) {fselectprimaryTruth=selectprimaryTruth;}
182   void SetCombinedNSigmaCut(Double_t NSigmaPID) {fNSigmaPID=NSigmaPID;}
183   void SetHighPtKaonNSigmaPID(Float_t HighPtKaonSigma,Float_t HighPtKaonNSigmaPID)
184   {
185     fHighPtKaonSigma=HighPtKaonSigma;
186     fHighPtKaonNSigmaPID=HighPtKaonNSigmaPID;
187   }
188   void IgnoreoverlappedTracks(Bool_t UseExclusiveNSigma){fUseExclusiveNSigma=UseExclusiveNSigma;}
189   void SetRemoveTracksT0Fill( Bool_t RemoveTracksT0Fill){fRemoveTracksT0Fill=RemoveTracksT0Fill;}
190   void SetPairSelectCharge(Int_t SelectCharge){fSelectCharge=SelectCharge;}
191   void SetTrigAssoSelectcharge( Int_t TriggerSelectCharge,Int_t AssociatedSelectCharge)
192   {
193     fTriggerSelectCharge=TriggerSelectCharge;
194     fAssociatedSelectCharge=AssociatedSelectCharge;
195   }
196   void SetEtaOrdering(Bool_t EtaOrdering){fEtaOrdering=EtaOrdering;}
197   void SetCutConversionsResonances( Bool_t CutConversions,Bool_t CutResonances)
198   {
199      fCutConversions=CutConversions;
200      fCutResonances=CutResonances;
201   }
202   void SetCleanUp(Bool_t InjectedSignals,Bool_t RemoveWeakDecays,Bool_t RemoveDuplicates)
203   {
204     fInjectedSignals=InjectedSignals;
205     fRemoveWeakDecays=RemoveWeakDecays;
206     fRemoveDuplicates=RemoveDuplicates;
207   }
208   void SetEfficiency(Bool_t fillefficiency,Bool_t applyTrigefficiency,Bool_t applyAssoefficiency)
209     {
210       ffillefficiency=fillefficiency;
211       fapplyTrigefficiency=applyTrigefficiency;
212       fapplyAssoefficiency=applyAssoefficiency;
213
214     }
215   void SetComboeffPtRange(Double_t minPtComboeff,Double_t maxPtComboeff) {
216     fminPtComboeff=minPtComboeff;
217     fmaxPtComboeff=maxPtComboeff;}
218   //only one can be kTRUE at a time(for the next two Setters)
219   void Setmesoneffrequired(Bool_t mesoneffrequired) {fmesoneffrequired=mesoneffrequired;}
220   void Setkaonprotoneffrequired(Bool_t kaonprotoneffrequired){fkaonprotoneffrequired=kaonprotoneffrequired;}
221    void SetOnlyOneEtaSide(Int_t OnlyOneEtaSide){fOnlyOneEtaSide=OnlyOneEtaSide;}
222    void SetRejectResonanceDaughters(Int_t RejectResonanceDaughters){fRejectResonanceDaughters=RejectResonanceDaughters;}
223
224 void SetTOFPIDVal(Bool_t RequestTOFPID,Float_t PtTOFPIDmin,Float_t PtTOFPIDmax)
225    {
226 fRequestTOFPID=RequestTOFPID;
227 fPtTOFPIDmin=PtTOFPIDmin;
228 fPtTOFPIDmax=PtTOFPIDmax;
229 }
230
231  void SetEffcorectionfilePathName(TString efffilename) {fefffilename=efffilename;}
232  
233
234   //PID Type
235   void SetPIDType(PIDType PIDmethod) { fPIDType = PIDmethod; }
236   PIDType GetPIDType() {return fPIDType; }
237   //NSigma cut
238   //set cut on beyesian probability
239   void SetBayesCut(Double_t cut){fBayesCut=cut;}
240   void SetdiffPIDcutvalues(Bool_t diffPIDcutvalues,Double_t PIDCutval1, Double_t PIDCutval2, Double_t PIDCutval3,Double_t PIDCutval4){
241     fdiffPIDcutvalues=diffPIDcutvalues;
242     fPIDCutval1=PIDCutval1;
243     fPIDCutval2=PIDCutval2;
244     fPIDCutval3=PIDCutval3;
245     fPIDCutval4=PIDCutval4;
246   }
247  void  SetRandomizeReactionPlane(Bool_t RandomizeReactionPlane){fRandomizeReactionPlane=RandomizeReactionPlane;}
248  
249    //****************************************************************************************EP related part
250   void OpenInfoCalbration(Int_t run);
251   void SetTPCclusterN(Int_t ncl){fNcluster=ncl;};
252    //****************************************************************************************EP related part
253
254
255  private:                                                                                      
256  //histograms
257     TList *fOutput;        //! Output list
258     TList *fOutputList;        //! Output list
259     TList *fList;              //! List for output objects
260
261
262     TString    fCentralityMethod;     // Method to determine centrality
263     Bool_t fPPVsMultUtils;//switch to ON quantile information for pp 7 TeV case
264     TString    fSampleType;     // pp,p-Pb,Pb-Pb
265     Bool_t fRequestEventPlane; //only for PbPb
266     Int_t    fnTracksVertex;        // QA tracks pointing to principal vertex
267     AliAODVertex* trkVtx;//!
268     Float_t zvtx;
269     Int_t    fFilterBit;         // track selection cuts
270      UInt_t         fTrackStatus;       // if non-0, the bits set in this variable are required for each track
271     Double_t       fSharedClusterCut;  // cut on shared clusters (only for AOD, give the actual cut value)
272     Double_t fSharedTPCmapCut;//cut on TPC shared map(set any non negative value to implement this cut automatically, no meaning of the value itself)
273     Double_t fSharedfraction_Pair_cut;//cut on pairs at the correlation level to check whether the correlating pair has large shared clusters(set fraction percentage to be set as cut off)
274     Int_t fVertextype;
275     Int_t skipParticlesAbove;
276     Double_t fzvtxcut;
277     Bool_t ffilltrigassoUNID;
278     Bool_t ffilltrigUNIDassoID;
279     Bool_t ffilltrigIDassoUNID;
280     Bool_t ffilltrigIDassoID;
281     Bool_t ffilltrigIDassoIDMCTRUTH;
282     Int_t fMaxNofMixingTracks;
283     Bool_t fPtOrderMCTruth;
284     Bool_t fPtOrderDataReco;
285     Bool_t fWeightPerEvent;
286     Bool_t fTriggerSpeciesSelection;
287     Bool_t fAssociatedSpeciesSelection;
288     Bool_t fRandomizeReactionPlane;
289     Int_t fTriggerSpecies;
290     Int_t fAssociatedSpecies;
291     TString fCustomBinning;//for setting customized binning
292     TString fBinningString;//final binning string
293     Bool_t fSelectHighestPtTrig;
294     Bool_t fcontainPIDtrig;
295     Bool_t fcontainPIDasso;
296     Int_t SetChargeAxis;
297      Bool_t frejectPileUp;
298      Bool_t fCheckFirstEventInChunk;
299     Float_t fminPt;
300     Float_t fmaxPt;
301     Float_t fmineta;
302     Float_t fmaxeta;
303     Bool_t fselectprimaryTruth;
304     Bool_t fonlyprimarydatareco;
305     Bool_t fdcacut;
306     Double_t fdcacutvalue;
307     Bool_t ffillhistQAReco;
308     Bool_t ffillhistQATruth;
309     Int_t kTrackVariablesPair ;
310     Double_t fminPtTrig;
311     Double_t fmaxPtTrig;
312     Double_t fminPtComboeff;
313     Double_t fmaxPtComboeff;
314     Double_t fminPtAsso;
315     Double_t fmaxPtAsso;
316     Double_t fmincentmult;
317     Double_t fmaxcentmult;
318     TH1F *fPriHistShare;//!
319     TH1F *fhistcentrality;//!
320     TH1F *fEventCounter; //!
321     TH2F *fEtaSpectrasso;//!
322     TH2F *fphiSpectraasso;//!
323     TH1F *MCtruthpt;//! 
324     TH1F *MCtrutheta;//! 
325     TH1F *MCtruthphi;//!
326     TH1F *MCtruthpionpt;//!
327     TH1F *MCtruthpioneta;//!
328     TH1F *MCtruthpionphi;//!
329     TH1F *MCtruthkaonpt;//!
330     TH1F *MCtruthkaoneta;//!
331     TH1F *MCtruthkaonphi;//!
332     TH1F *MCtruthprotonpt;//!
333     TH1F *MCtruthprotoneta;//!
334     TH1F *MCtruthprotonphi;//! 
335     TH2F *fPioncont;//!
336     TH2F *fKaoncont;//!
337     TH2F *fProtoncont;//!
338     TH2F *fUNIDcont;//!
339     TH2F *fEventno;//!
340     TH2F *fEventnobaryon;//!
341     TH2F *fEventnomeson;//!
342     TH2F *fhistJetTrigestimate;//!
343     TH3F* fTwoTrackDistancePtdip;//!
344     TH3F* fTwoTrackDistancePtdipmix;//!
345     TH3F* fTwoTrackDistancePt[2];    //! control histograms for two-track efficiency study: dphi*_min vs deta (0 = before cut, 1 = after cut)
346     TH3F* fTwoTrackDistancePtmix[2];    //! control histograms for two-track efficiency study: dphi*_min vs deta (0 = before cut, 1 = after cut)
347
348     TH2D* fCentralityCorrelation;  //! centrality vs Tracks multiplicity
349  //VZERO calibration
350   TH1F *fHistVZEROAGainEqualizationMap;//VZERO calibration map
351   TH1F *fHistVZEROCGainEqualizationMap;//VZERO calibration map
352   TH2F *fHistVZEROChannelGainEqualizationMap; //VZERO calibration map
353   TH1* fCentralityWeights;                   // for centrality flattening
354
355     TH2F *fHistCentStats; //!centrality stats
356     TH2F *fHistRefmult;//!
357     TH2F *fHistEQVZEROvsTPCmultiplicity;//!
358     TH2F *fHistEQVZEROAvsTPCmultiplicity;//!
359     TH2F *fHistEQVZEROCvsTPCmultiplicity;//!
360     TH2F *fHistVZEROCvsEQVZEROCmultiplicity;//!
361     TH2F *fHistVZEROAvsEQVZEROAmultiplicity;//!
362     TH2F *fHistVZEROCvsVZEROAmultiplicity;//!
363     TH2F *fHistEQVZEROCvsEQVZEROAmultiplicity;//!
364     TH2F *fHistVZEROSignal;//!
365     TH2F *fHistEventPlaneTruth;//!
366     TH2D *fHistPsiMinusPhi;//! psi - phi QA histogram
367     TH3F *fEventPlanePID;//!
368    //****************************************************************************************EP related part
369
370     Float_t evplaneMC,fgPsi2v0a,fgPsi2v0c,fgPsi2tpc; // current Psi2
371    Float_t fgPsi3v0a,fgPsi3v0c,fgPsi3tpc; // current Psi3
372    Float_t fgPsi2v0aMC,fgPsi2v0cMC,fgPsi2tpcMC; // current Psi2
373    Float_t fgPsi3v0aMC,fgPsi3v0cMC,fgPsi3tpcMC,gReactionPlane; // current Psi3
374   Bool_t fV2; // switch to set the harmonics
375   Bool_t fV3; // switch to set the harmonics
376   Bool_t fIsAfter2011; // switch for 2011 and later runs
377
378   //    Int_t nCentrBin = 9;          //  cenrality bins
379
380   //
381   // Cuts and options
382   //
383
384   Int_t fRun;                       // current run checked to load VZERO calibrations
385
386   Int_t fNcluster;           // Numer of TPC cluster required
387
388   TString fEPdet; //Set the name of the event plane to be used to reconstruct the event plane
389     
390   // Output objects
391   TProfile *fMultV0;                //! object containing VZERO calibration information
392   Float_t fV0Cpol;          //! loaded by OADB
393   Float_t fV0Apol;          //! loaded by OADB
394   Float_t fMeanQ[9][2][2];           // and recentering
395   Float_t fWidthQ[9][2][2];          // ...
396   Float_t fMeanQv3[9][2][2];         // also for v3
397   Float_t fWidthQv3[9][2][2];        // ...
398
399   TProfile *fHResTPCv0A2;   //! TProfile for subevent resolution (output)
400   TProfile *fHResTPCv0C2;   //! TProfile for subevent resolution (output)
401   TProfile *fHResv0Cv0A2;   //! TProfile for subevent resolution (output)
402   TProfile *fHResTPCv0A3;    //! also for v3
403   TProfile *fHResTPCv0C3;   //! also for v3
404   TProfile *fHResv0Cv0A3;   //! also for v3
405
406  TProfile *fHResMA2;   //! TProfile for subevent resolution (output)
407   TProfile *fHResMC2;   //! TProfile for subevent resolution (output)
408   TProfile *fHResAC2;   //! TProfile for subevent resolution (output)
409   TProfile *fHResMA3;    //! also for v3
410   TProfile *fHResMC3;   //! also for v3
411   TProfile *fHResAC3;   //! also for v3
412
413   TH2F *fPhiRPTPC;          //! EP distribution vs. centrality (v2)
414   TH2F *fPhiRPTPCv3;          //! EP distribution vs. centrality (v2)
415   TH2F *fPhiRPv0A;          //! EP distribution vs. centrality (v2)
416   TH2F *fPhiRPv0C;          //! EP distribution vs. centrality (v2)
417   TH2F *fPhiRPv0Av3;      //! EP distribution vs. centrality (v3)
418   TH2F *fPhiRPv0Cv3;      //! EP distribution vs. centrality (v3)
419     //****************************************************************************************EP related part
420
421     TH2F* fControlConvResoncances; //! control histograms for cuts on conversions and resonances
422
423     TH2F *fHistoTPCdEdx;//!
424     TH2F *fHistoTOFbeta;//!
425     TH3F *fTPCTOFPion3d;//!
426     TH3F *fTPCTOFKaon3d;//!
427     TH3F *fTPCTOFProton3d;//!
428     TH1F *fPionPt;//!
429     TH1F *fPionEta;//!
430     TH1F *fPionPhi;//!
431     TH1F *fKaonPt;//!
432     TH1F *fKaonEta;//!
433     TH1F *fKaonPhi;//!
434     TH1F *fProtonPt;//!
435     TH1F *fProtonEta;//!
436     TH1F *fProtonPhi;//!
437     // TH3F *fHistocentNSigmaTPC;//! nsigma TPC
438     // TH3F *fHistocentNSigmaTOF;//! nsigma TOF 
439  
440     AliTHn *fCorrelatonTruthPrimary;//!
441     AliTHn *fCorrelatonTruthPrimarymix;//!
442     AliTHn *fTHnCorrUNID;//!
443     AliTHn *fTHnCorrUNIDmix;//!
444     AliTHn *fTHnCorrID;//!
445     AliTHn *fTHnCorrIDmix;//!
446     AliTHn *fTHnCorrIDUNID;//!
447     AliTHn *fTHnCorrIDUNIDmix;//!
448     AliTHn *fTHnTrigcount;//!
449     AliTHn *fTHnTrigcountMCTruthPrim;//!
450     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
451
452     
453     TH1F *fHistQA[16]; //!                  
454      
455    
456     THnSparse *effcorection[6];//!
457     // THnF *effmap[6];  
458
459     Int_t ClassifyTrack(AliAODTrack* track,AliAODVertex* vertex,Float_t magfield,Bool_t fill);
460   Double_t* GetBinning(const char* configuration, const char* tag, Int_t& nBins);
461
462
463   void Fillcorrelation(Float_t ReactionPlane,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; 
464  Bool_t CalculateSharedFraction(const TBits *triggerPadMap,const TBits *assocPadMap,const TBits *triggerShareMap,const TBits *assocShareMap);
465  Float_t GetTrackbyTrackeffvalue(AliAODTrack* track,Double_t cent,Float_t evzvtx, Int_t parpid);
466
467  //Fill PID and Event planes
468  void FillPIDEventPlane(Double_t centrality,Int_t par,Float_t trigphi,Float_t fReactionPlane);
469
470 //Mixing functions
471  // void DefineEventPool();
472   AliEventPoolManager    *fPoolMgr;//! 
473   TClonesArray          *fArrayMC;//!
474   TString          fAnalysisType;          // "MC", "ESD", "AOD"
475   TString fefffilename;
476     //PID part histograms
477
478   //PID functions
479     Bool_t HasTPCPID(AliAODTrack *track) const; // has TPC PID
480     Bool_t HasTOFPID(AliAODTrack *track) const; // has TOF PID
481     Double_t GetBeta(AliAODTrack *track);
482     void CalculateNSigmas(AliAODTrack *track, Bool_t FIllQAHistos);
483     Int_t FindMinNSigma(AliAODTrack *track, Bool_t FIllQAHistos);
484     Bool_t* GetDoubleCounting(AliAODTrack * trk, Bool_t FIllQAHistos);
485     Int_t GetParticle(AliAODTrack * trk, Bool_t FIllQAHistos);  
486  
487      TH2F* GetHistogram2D(const char * name);//return histogram "name" from fOutputList
488
489      Bool_t ftwoTrackEfficiencyCutDataReco; 
490     Float_t fTwoTrackCutMinRadius;
491     Float_t fTwoTrackCutMaxRadius;         
492    Float_t twoTrackEfficiencyCutValue;
493   //Pid objects
494   AliPIDResponse *fPID; //! PID
495   AliPIDCombined   *fPIDCombined;     //! PIDCombined
496
497   //set PID Combined
498   void SetPIDCombined(AliPIDCombined *obj){fPIDCombined=obj;}
499   AliPIDCombined *GetPIDCombined(){return fPIDCombined;}
500  
501   Double_t GetBayesCut(){return fBayesCut;}
502   Int_t GetIDBayes(AliAODTrack *trk, Bool_t FIllQAHistos);//calculate the PID according to bayesian PID
503   UInt_t CalcPIDCombined(AliAODTrack *track,Int_t detMask, Double_t* prob) const;
504   Bool_t* GetAllCompatibleIdentitiesNSigma(AliAODTrack * trk, Bool_t FIllQAHistos);//All the identities are true
505
506
507   Int_t eventno;
508   Float_t fPtTOFPIDmin; //lower pt bound for the TOCTOF combined circular pid
509   Float_t fPtTOFPIDmax; //uper pt bound for the TOCTOF combined circular pid
510   Bool_t fRequestTOFPID;//if true returns kSpUndefined if the TOF signal is missing
511   PIDType fPIDType; // PID type  Double_t fNSigmaPID; // number of sigma for PID cut
512   Bool_t fFIllPIDQAHistos; //Switch for filling the nSigma histos
513   Double_t fNSigmaPID; // number of sigma for PID cut
514   Double_t fBayesCut; // Cut on Bayesian probability
515  Bool_t fdiffPIDcutvalues;
516  Double_t fPIDCutval1;
517  Double_t fPIDCutval2;
518  Double_t fPIDCutval3;
519  Double_t fPIDCutval4;
520
521   Float_t fHighPtKaonNSigmaPID;// number of sigma for PID cut for Kaons above fHighPtKaonSigma(-1 default, no cut applied)
522   Float_t fHighPtKaonSigma;//lower pt bound for the fHighPtKaonNSigmaPID to be set >0(i.e. to make it applicable)
523   Bool_t fUseExclusiveNSigma;//if true returns the identity only if no double counting(i.e not in the overlap area)
524   Bool_t fRemoveTracksT0Fill;//if true remove tracks for which only StartTime from To-Fill is available (worst resolution)
525  Int_t fSelectCharge;           // (un)like sign selection when building correlations: 0: no selection; 1: unlike sign; 2: like sign
526     Int_t fTriggerSelectCharge;    // select charge of trigger particle: 1: positive; -1 negative
527     Int_t fAssociatedSelectCharge; // select charge of associated particle: 1: positive; -1 negative
528     Float_t fTriggerRestrictEta;   // restrict eta range for trigger particle (default: -1 [off])
529     Bool_t fEtaOrdering;           // eta ordering, see AliUEHistograms.h for documentation
530     Bool_t fCutConversions;        // cut on conversions (inv mass)
531     Bool_t fCutResonances;         // cut on resonances (inv mass)
532     Int_t fRejectResonanceDaughters; // reject all daughters of all resonance candidates (1: test method (cut at m_inv=0.9); 2: k0; 3: lambda)
533     Int_t fOnlyOneEtaSide;       // decides that only trigger particle from one eta side are considered (0 = all; -1 = negative, 1 = positive)
534     Bool_t fInjectedSignals;      // check header to skip injected signals in MC
535     Bool_t fRemoveWeakDecays;      // remove secondaries from weak decays from tracks and particles
536     Bool_t fRemoveDuplicates;// remove particles with the same label (double reconstruction)
537     Bool_t fapplyTrigefficiency;//if kTRUE then eff correction calculation starts
538     Bool_t fapplyAssoefficiency;//if kTRUE then eff correction calculation starts
539     Bool_t ffillefficiency;  //if kTRUE then THNsparses used for eff. calculation are filled up
540     Bool_t fmesoneffrequired;
541     Bool_t fkaonprotoneffrequired;
542    AliAnalysisUtils*     fAnalysisUtils;      // points to class with common analysis utilities
543   TFormula*      fDCAXYCut;          // additional pt dependent cut on DCA XY (only for AOD)
544
545
546   Float_t fnsigmas[NSpecies][NSigmaPIDType+1]; //nsigma values
547   Bool_t fHasDoubleCounting[NSpecies];//array with compatible identities
548
549   //Int_t fPIDMethod; // PID method
550
551  //functions
552   Bool_t CheckTrack(AliAODTrack * part);
553   Float_t PhiRange(Float_t DPhi);
554   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);
555 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);
556   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);
557   TObjArray* CloneAndReduceTrackList(TObjArray* tracks);
558
559
560   void ShiftTracks(TObjArray* tracks, Double_t angle);
561   Bool_t AcceptEventCentralityWeight(Double_t centrality);
562
563   //get event plane
564   Float_t GetEventPlane(AliAODEvent *event,Bool_t truth,Double_t v0Centr);
565   Double_t GetAcceptedEventMultiplicity(AliAODEvent *aod,Bool_t truth);//returns centrality after event(mainly vertex) selection IsEventAccepted  GetAcceptedEventMultiplicity
566   
567   //get vzero equalization
568   Double_t GetEqualizationFactor(Int_t run, const char* side);
569   Double_t GetChannelEqualizationFactor(Int_t run,Int_t channel);
570   void SetVZEROCalibrationFile(const char* filename,const char* lhcPeriod);
571   void SetCentralityWeights(TH1* hist) { fCentralityWeights = hist; }
572
573   Double_t GetRefMultiOrCentrality(AliAODEvent *event, Bool_t truth);
574   Double_t GetReferenceMultiplicityVZEROFromAOD(AliAODEvent *event);//mainly important for pp 7 TeV
575
576
577
578
579     
580     AliTwoParticlePIDCorr(const AliTwoParticlePIDCorr&); // not implemented
581     AliTwoParticlePIDCorr& operator=(const AliTwoParticlePIDCorr&); // not implemented
582     
583     ClassDef(AliTwoParticlePIDCorr, 1); // example of analysis
584 };
585 class LRCParticlePID : public TObject {
586 public:
587  LRCParticlePID(Int_t par,Short_t icharge,Float_t pt,Float_t eta, Float_t phi,Float_t effcorrectionval,const TBits *clustermap,const TBits *sharemap)
588    :fparticle(par),fcharge(icharge),fPt(pt), fEta(eta), fPhi(phi),feffcorrectionval(effcorrectionval),fTPCClusterMap(clustermap),fTPCHitShareMap(sharemap) {}
589   virtual ~LRCParticlePID() {}
590
591   
592     virtual Float_t Eta()        const { return fEta; }
593     virtual Float_t Phi()        const { return fPhi; }
594     virtual Float_t Pt() const { return fPt; }
595     Int_t getparticle() const {return fparticle;}
596     virtual Short_t Charge()      const { return fcharge; }
597     Float_t geteffcorrectionval() const {return feffcorrectionval;}
598     virtual Bool_t IsEqual(const TObject* obj) const { return (obj->GetUniqueID() == GetUniqueID()); }
599     virtual void SetPhi(Double_t phiv) { fPhi = phiv; }
600     virtual const TBits * GetTPCPadMap() {return fTPCClusterMap; }
601     virtual const TBits * GetTPCSharedMap() {return fTPCHitShareMap; }
602
603 private:
604   LRCParticlePID(const LRCParticlePID&);  // not implemented
605    LRCParticlePID& operator=(const LRCParticlePID&);  // not implemented
606   
607   Int_t fparticle;
608   Short_t fcharge;
609   Float_t fPt;
610   Float_t fEta;
611   Float_t fPhi;
612   Float_t feffcorrectionval;
613    const TBits   *fTPCClusterMap;
614    const TBits   *fTPCHitShareMap;
615   ClassDef(LRCParticlePID, 1);
616 } ;
617
618 #endif
619
620 //(fSampleType=="pp_2_76" || fCentralityMethod.EndsWith("_MANUAL"))
621 //(fSampleType=="pp_2_76" || fCentralityMethod.EndsWith("_MANUAL") || (fSampleType=="pp_7" && fPPVsMultUtils==kFALSE))
622 //(fCentralityMethod.EndsWith("_MANUAL"))