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