]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/TriggerPID/AliTwoParticlePIDCorr.h
99a328d62d076309cd10745d93160dd60b6cfa30
[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 AliAODEvent;
17 class AliAODTrack;
18 class AliAODVertex;
19 class AliEventPoolManager;
20 class TFormula;
21 //class AliAnalysisUtils;
22 class LRCParticlePID;
23 class AliVParticle;
24 class AliCFContainer;
25 class AliCFGridSparse;
26 class THnBase;
27 class AliTHn;
28
29
30 #include <TObject.h> //LRCParticlePID is a derived class from"TObject"
31 #include "TMath.h"
32 #include "TNamed.h"
33 #include "AliUEHist.h"
34 #include "AliPID.h"
35 #include "AliAnalysisTask.h"
36 #include "AliUEHist.h"
37 #include "TString.h"
38 #include "AliVParticle.h"
39 #include "TParticle.h"
40 #include "AliLog.h"
41 #include "AliTHn.h"
42
43
44
45 #ifndef ALIANALYSISTASKSE_H
46 #include "AliAnalysisTaskSE.h"
47 #endif
48
49 namespace AliPIDNameSpace {
50   
51   enum PIDType
52   {
53     NSigmaTPC = 0,
54     NSigmaTOF,
55     NSigmaTPCTOF,  // squared sum
56     NSigmaPIDType=NSigmaTPCTOF
57   };
58     
59   enum AliDetectorType
60   {
61     fITS = 0,
62     fTPC,
63     fTOF,
64     fNDetectors
65   };
66   
67   
68   enum AliParticleSpecies
69   {
70     SpPion = 0,
71     SpKaon,
72     SpProton,
73     unidentified,
74     NSpecies=unidentified,
75     SpUndefined=999
76   }; // Particle species used in plotting
77   
78   
79   enum AliCharge
80   {
81     Posch = 0,
82     Negch,
83     NCharge
84   };
85 }
86
87
88 using namespace AliPIDNameSpace;
89
90 class AliTwoParticlePIDCorr : public AliAnalysisTaskSE {
91  public:
92     AliTwoParticlePIDCorr();
93     AliTwoParticlePIDCorr(const char *name);
94     virtual ~AliTwoParticlePIDCorr();
95     
96     virtual void     UserCreateOutputObjects();
97     virtual void     UserExec(Option_t *option);
98     virtual void    doAODevent();
99     virtual void    doMCAODevent();
100     virtual void     Terminate(Option_t *);
101   void           SetSharedClusterCut(Double_t value) { fSharedClusterCut = value; }
102
103   void SetVertextype(Int_t Vertextype){fVertextype=Vertextype;}                                                 //Check it every time
104     void SetZvtxcut(Double_t zvtxcut) {fzvtxcut=zvtxcut;}
105     void SetCustomBinning(TString receivedCustomBinning) { fCustomBinning = receivedCustomBinning; }
106     void SetMaxNofMixingTracks(Int_t MaxNofMixingTracks) {fMaxNofMixingTracks=MaxNofMixingTracks;}               //Check it every time
107   void SetCentralityEstimator(TString CentralityMethod) { fCentralityMethod = CentralityMethod;}
108   void SetSampleType(TString SampleType) {fSampleType=SampleType;}
109   void SetAnalysisType(TString AnalysisType){fAnalysisType=AnalysisType;}
110   void SetFilterBit(Int_t FilterBit) {fFilterBit=FilterBit;}
111   void SetfilltrigassoUNID(Bool_t filltrigassoUNID){ffilltrigassoUNID=filltrigassoUNID;}
112   void SetfilltrigUNIDassoID(Bool_t filltrigUNIDassoID){ffilltrigUNIDassoID=filltrigUNIDassoID;}
113   void SetfilltrigIDassoUNID(Bool_t filltrigIDassoUNID){ffilltrigIDassoUNID=filltrigIDassoUNID;}
114   void SetfilltrigIDassoID(Bool_t filltrigIDassoID){ ffilltrigIDassoID=filltrigIDassoID;}
115   void SetfilltrigIDassoIDMCTRUTH(Bool_t filltrigIDassoIDMCTRUTH){ffilltrigIDassoIDMCTRUTH=filltrigIDassoIDMCTRUTH;}
116   void SetSelectHighestPtTrig(Bool_t SelectHighestPtTrig){fSelectHighestPtTrig=SelectHighestPtTrig;}
117   void SetTriggerSpeciesSelection(Bool_t TriggerSpeciesSelection,Int_t TriggerSpecies,Bool_t containPIDtrig){
118     fTriggerSpeciesSelection=TriggerSpeciesSelection;//if it is KTRUE then Set containPIDtrig=kFALSE
119     fTriggerSpecies=TriggerSpecies;
120     fcontainPIDtrig=containPIDtrig;
121   }
122   void SetAssociatedSpeciesSelection(Bool_t AssociatedSpeciesSelection,Int_t AssociatedSpecies, Bool_t containPIDasso)
123   {
124     fAssociatedSpeciesSelection=AssociatedSpeciesSelection;//if it is KTRUE then Set containPIDasso=kFALSE
125     fAssociatedSpecies=AssociatedSpecies;
126     fcontainPIDasso=containPIDasso;
127   }
128
129  void SetFIllPIDQAHistos(Bool_t FIllPIDQAHistos){fFIllPIDQAHistos=FIllPIDQAHistos;}
130   // void SetRejectPileUp(Bool_t rejectPileUp) {frejectPileUp=rejectPileUp;}
131   void SetKinematicCuts(Float_t minPt, Float_t maxPt,Float_t mineta,Float_t maxeta)
132   {
133     fminPt=minPt;
134     fmaxPt=maxPt;
135     fmineta=mineta;
136     fmaxeta=maxeta;
137   }
138   void SetAsymmetricnSigmaCut( Float_t minprotonsigmacut,Float_t maxprotonsigmacut,Float_t minpionsigmacut,Float_t maxpionsigmacut)
139   {
140     fminprotonsigmacut=minprotonsigmacut;
141     fmaxprotonsigmacut=maxprotonsigmacut;
142     fminpionsigmacut=minpionsigmacut;
143     fmaxpionsigmacut=maxpionsigmacut;
144   }
145   void SetDcaCut(Bool_t dcacut,Double_t dcacutvalue)
146   {
147     fdcacut=dcacut;
148     fdcacutvalue=dcacutvalue;
149   }
150   void SetfillHistQA(Bool_t fillhistQAReco,Bool_t fillhistQATruth)
151   {
152     ffillhistQAReco=fillhistQAReco;
153     ffillhistQATruth=fillhistQATruth;
154   }
155   void Setselectprimarydatareco(Bool_t onlyprimarydatareco) {fonlyprimarydatareco=onlyprimarydatareco;}
156   void SetselectprimaryTruth(Bool_t selectprimaryTruth) {fselectprimaryTruth=selectprimaryTruth;}
157   void SetCombinedNSigmaCut(Double_t NSigmaPID) {fNSigmaPID=NSigmaPID;}
158   void IgnoreoverlappedTracks(Bool_t UseExclusiveNSigma){fUseExclusiveNSigma=UseExclusiveNSigma;}
159   void SetRemoveTracksT0Fill( Bool_t RemoveTracksT0Fill){fRemoveTracksT0Fill=RemoveTracksT0Fill;}
160   void SetPairSelectCharge(Int_t SelectCharge){fSelectCharge=SelectCharge;}
161   void SetTrigAssoSelectcharge( Int_t TriggerSelectCharge,Int_t AssociatedSelectCharge)
162   {
163     fTriggerSelectCharge=TriggerSelectCharge;
164     fAssociatedSelectCharge=AssociatedSelectCharge;
165   }
166   void SetEtaOrdering(Bool_t EtaOrdering){fEtaOrdering=EtaOrdering;}
167   void SetCutConversionsResonances( Bool_t CutConversions,Bool_t CutResonances)
168   {
169      fCutConversions=CutConversions;
170      fCutResonances=CutResonances;
171   }
172   void SetCleanUp(Bool_t InjectedSignals,Bool_t RemoveWeakDecays,Bool_t RemoveDuplicates)
173   {
174     fInjectedSignals=InjectedSignals;
175     fRemoveWeakDecays=RemoveWeakDecays;
176     fRemoveDuplicates=RemoveDuplicates;
177   }
178   void SetEfficiency(Bool_t fillefficiency,Bool_t applyTrigefficiency,Bool_t applyAssoefficiency)
179     {
180       ffillefficiency=fillefficiency;
181       fapplyTrigefficiency=applyTrigefficiency;
182       fapplyAssoefficiency=applyAssoefficiency;
183
184     }
185   void SetComboeffPtRange(Double_t minPtComboeff,Double_t maxPtComboeff) {
186     fminPtComboeff=minPtComboeff;
187     fmaxPtComboeff=maxPtComboeff;}
188   //only one can be kTRUE at a time(for the next two Setters)
189   void Setmesoneffrequired(Bool_t mesoneffrequired) {fmesoneffrequired=mesoneffrequired;}
190   void Setkaonprotoneffrequired(Bool_t kaonprotoneffrequired){fkaonprotoneffrequired=kaonprotoneffrequired;}
191    void SetOnlyOneEtaSide(Int_t OnlyOneEtaSide){fOnlyOneEtaSide=OnlyOneEtaSide;}
192    void SetRejectResonanceDaughters(Int_t RejectResonanceDaughters){fRejectResonanceDaughters=RejectResonanceDaughters;}
193
194 void SetTOFPIDVal(Bool_t RequestTOFPID,Float_t PtTOFPIDmin,Float_t PtTOFPIDmax)
195    {
196 fRequestTOFPID=RequestTOFPID;
197 fPtTOFPIDmin=PtTOFPIDmin;
198 fPtTOFPIDmax=PtTOFPIDmax;
199 }
200
201  void SetEffcorectionfilePathName(TString efffilename) {fefffilename=efffilename;}
202  
203  private:
204  //histograms
205     TList *fOutput;        //! Output list
206     TList *fOutputList;        //! Output list
207
208     TString    fCentralityMethod;     // Method to determine centrality
209     TString    fSampleType;     // pp,p-Pb,Pb-Pb
210     Int_t    fnTracksVertex;        // QA tracks pointing to principal vertex
211     AliAODVertex* trkVtx;//!
212     Float_t zvtx;
213     Int_t    fFilterBit;         // track selection cuts
214     Double_t       fSharedClusterCut;  // cut on shared clusters (only for AOD)
215     Int_t fVertextype;
216     Double_t fzvtxcut;
217     Bool_t ffilltrigassoUNID;
218     Bool_t ffilltrigUNIDassoID;
219     Bool_t ffilltrigIDassoUNID;
220     Bool_t ffilltrigIDassoID;
221     Bool_t ffilltrigIDassoIDMCTRUTH;
222     Int_t fMaxNofMixingTracks;
223     Bool_t fPtOrderMCTruth;
224     Bool_t fTriggerSpeciesSelection;
225     Bool_t fAssociatedSpeciesSelection;
226     Int_t fTriggerSpecies;
227     Int_t fAssociatedSpecies;
228     TString fCustomBinning;//for setting customized binning
229     TString fBinningString;//final binning string
230     Bool_t fSelectHighestPtTrig;
231     Bool_t fcontainPIDtrig;
232     Bool_t fcontainPIDasso;
233     // Bool_t frejectPileUp;
234     Float_t fminPt;
235     Float_t fmaxPt;
236     Float_t fmineta;
237     Float_t fmaxeta;
238     Float_t fminprotonsigmacut;
239     Float_t fmaxprotonsigmacut;
240     Float_t fminpionsigmacut;
241     Float_t fmaxpionsigmacut;
242     Bool_t fselectprimaryTruth;
243     Bool_t fonlyprimarydatareco;
244     Bool_t fdcacut;
245     Double_t fdcacutvalue;
246     Bool_t ffillhistQAReco;
247     Bool_t ffillhistQATruth;
248     Int_t kTrackVariablesPair ;
249     Double_t fminPtTrig;
250     Double_t fmaxPtTrig;
251     Double_t fminPtComboeff;
252     Double_t fmaxPtComboeff;
253     Double_t fminPtAsso;
254     Double_t fmaxPtAsso;
255     TH1F *fhistcentrality;//!
256     TH1F *fEventCounter; //!
257     TH2F *fEtaSpectrasso;//!
258     TH2F *fphiSpectraasso;//!
259     TH1F *MCtruthpt;//! 
260     TH1F *MCtrutheta;//! 
261     TH1F *MCtruthphi;//!
262     TH1F *MCtruthpionpt;//!
263     TH1F *MCtruthpioneta;//!
264     TH1F *MCtruthpionphi;//!
265     TH1F *MCtruthkaonpt;//!
266     TH1F *MCtruthkaoneta;//!
267     TH1F *MCtruthkaonphi;//!
268     TH1F *MCtruthprotonpt;//!
269     TH1F *MCtruthprotoneta;//!
270     TH1F *MCtruthprotonphi;//! 
271     TH2F *fPioncont;//!
272     TH2F *fKaoncont;//!
273     TH2F *fProtoncont;//!
274     TH2F *fEventno;//!
275     TH2F *fEventnobaryon;//!
276     TH2F *fEventnomeson;//!
277
278     TH2D* fCentralityCorrelation;  //! centrality vs multiplicity
279
280     TH2F *fHistoTPCdEdx;//!
281     TH2F *fHistoTOFbeta;//!
282     TH3F *fTPCTOFPion3d;//!
283     TH3F *fTPCTOFKaon3d;//!
284     TH3F *fTPCTOFProton3d;//!
285     TH1F *fPionPt;//!
286     TH1F *fPionEta;//!
287     TH1F *fPionPhi;//!
288     TH1F *fKaonPt;//!
289     TH1F *fKaonEta;//!
290     TH1F *fKaonPhi;//!
291     TH1F *fProtonPt;//!
292     TH1F *fProtonEta;//!
293     TH1F *fProtonPhi;//!
294     // TH3F *fHistocentNSigmaTPC;//! nsigma TPC
295     // TH3F *fHistocentNSigmaTOF;//! nsigma TOF 
296     
297     AliTHn *fCorrelatonTruthPrimary;//!
298     AliTHn *fCorrelatonTruthPrimarymix;//!
299     AliTHn *fTHnCorrUNID;//!
300     AliTHn *fTHnCorrUNIDmix;//!
301     AliTHn *fTHnCorrID;//!
302     AliTHn *fTHnCorrIDmix;//!
303     AliTHn *fTHnCorrIDUNID;//!
304     AliTHn *fTHnCorrIDUNIDmix;//!
305     AliTHn *fTHnTrigcount;//!
306     AliTHn *fTHnTrigcountMCTruthPrim;//!
307     
308     TH1F *fHistQA[16]; //!
309      
310     THnSparse *fTHnrecomatchedallPid[6];//!                 //0 pion, 1 kaon,2 proton,3 mesons,4 kaons+protons,5 all
311     THnSparse *fTHngenprimPidTruth[6];//!
312     THnSparse *effcorection[6];//!
313     // THnF *effmap[6];  
314     //TH2F* fControlConvResoncances; //! control histograms for cuts on conversions and resonances
315
316     Int_t ClassifyTrack(AliAODTrack* track,AliAODVertex* vertex,Float_t magfield);
317   Double_t* GetBinning(const char* configuration, const char* tag, Int_t& nBins);
318
319
320   void Fillcorrelation(TObjArray *trackstrig,TObjArray *tracksasso,Double_t cent,Float_t vtx,Float_t weight,Float_t bSign,Bool_t fPtOrder,Bool_t twoTrackEfficiencyCut,Bool_t mixcase,TString fillup);//mixcase=kTRUE in case of mixing; 
321  Float_t GetTrackbyTrackeffvalue(AliAODTrack* track,Double_t cent,Float_t evzvtx, Int_t parpid);
322
323 //Mixing functions
324   void DefineEventPool();
325   // AliAnalysisUtils *fUtils;
326   AliEventPoolManager    *fPoolMgr;//! 
327   TClonesArray          *fArrayMC;//!
328   TString          fAnalysisType;          // "MC", "ESD", "AOD"
329   TString fefffilename;
330
331     //PID part histograms
332
333   //PID functions
334     Bool_t HasTPCPID(AliAODTrack *track) const; // has TPC PID
335     Bool_t HasTOFPID(AliAODTrack *track) const; // has TOF PID
336     Double_t GetBeta(AliAODTrack *track);
337     void CalculateNSigmas(AliAODTrack *track, Bool_t FIllQAHistos);
338     Int_t FindMinNSigma(AliAODTrack *track, Bool_t FIllQAHistos);
339     Bool_t* GetDoubleCounting(AliAODTrack * trk, Bool_t FIllQAHistos);
340     Int_t GetParticle(AliAODTrack * trk, Bool_t FIllQAHistos);  
341  
342      TH2F* GetHistogram2D(const char * name);//return histogram "name" from fOutputList
343
344            
345    Float_t twoTrackEfficiencyCutValue;
346   //Pid objects
347   AliPIDResponse *fPID; //! PID
348   Int_t eventno;
349   Float_t fPtTOFPIDmin; //lower pt bound for the TOCTOF combined circular pid
350   Float_t fPtTOFPIDmax; //uper pt bound for the TOCTOF combined circular pid
351   Bool_t fRequestTOFPID;//if true returns kSpUndefined if the TOF signal is missing
352   PIDType fPIDType; // PID type  Double_t fNSigmaPID; // number of sigma for PID cut
353   Bool_t fFIllPIDQAHistos; //Switch for filling the nSigma histos
354   Double_t fNSigmaPID; // number of sigma for PID cut
355   Bool_t fUseExclusiveNSigma;//if true returns the identity only if no double counting(i.e not in the overlap area)
356   Bool_t fRemoveTracksT0Fill;//if true remove tracks for which only StartTime from To-Fill is available (worst resolution)
357  Int_t fSelectCharge;           // (un)like sign selection when building correlations: 0: no selection; 1: unlike sign; 2: like sign
358     Int_t fTriggerSelectCharge;    // select charge of trigger particle: 1: positive; -1 negative
359     Int_t fAssociatedSelectCharge; // select charge of associated particle: 1: positive; -1 negative
360     Float_t fTriggerRestrictEta;   // restrict eta range for trigger particle (default: -1 [off])
361     Bool_t fEtaOrdering;           // eta ordering, see AliUEHistograms.h for documentation
362     Bool_t fCutConversions;        // cut on conversions (inv mass)
363     Bool_t fCutResonances;         // cut on resonances (inv mass)
364     Int_t fRejectResonanceDaughters; // reject all daughters of all resonance candidates (1: test method (cut at m_inv=0.9); 2: k0; 3: lambda)
365     Int_t fOnlyOneEtaSide;       // decides that only trigger particle from one eta side are considered (0 = all; -1 = negative, 1 = positive)
366     Bool_t fInjectedSignals;      // check header to skip injected signals in MC
367     Bool_t fRemoveWeakDecays;      // remove secondaries from weak decays from tracks and particles
368     Bool_t fRemoveDuplicates;// remove particles with the same label (double reconstruction)
369     Bool_t fapplyTrigefficiency;//if kTRUE then eff correction calculation starts
370     Bool_t fapplyAssoefficiency;//if kTRUE then eff correction calculation starts
371     Bool_t ffillefficiency;//if kTRUE then THNsparses used for eff. calculation are filled up
372     Bool_t fmesoneffrequired;
373     Bool_t fkaonprotoneffrequired;
374     //  AliAnalysisUtils*     fAnalysisUtils;      // points to class with common analysis utilities
375   TFormula*      fDCAXYCut;          // additional pt dependent cut on DCA XY (only for AOD)
376
377
378   Float_t fnsigmas[NSpecies][NSigmaPIDType+1]; //nsigma values
379   Bool_t fHasDoubleCounting[NSpecies];//array with compatible identities
380
381   //Int_t fPIDMethod; // PID method
382
383  //functions
384   Float_t PhiRange(Float_t DPhi);
385   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);
386 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);
387   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);
388   TObjArray* CloneAndReduceTrackList(TObjArray* tracks);
389            
390     
391     AliTwoParticlePIDCorr(const AliTwoParticlePIDCorr&); // not implemented
392     AliTwoParticlePIDCorr& operator=(const AliTwoParticlePIDCorr&); // not implemented
393     
394     ClassDef(AliTwoParticlePIDCorr, 1); // example of analysis
395 };
396 class LRCParticlePID : public TObject {
397 public:
398  LRCParticlePID(Int_t par,Short_t icharge,Float_t pt,Float_t eta, Float_t phi,Float_t effcorrectionval)
399    :fparticle(par),fcharge(icharge),fPt(pt), fEta(eta), fPhi(phi),feffcorrectionval(effcorrectionval)  {}
400   virtual ~LRCParticlePID() {}
401
402   
403     virtual Float_t Eta()        const { return fEta; }
404     virtual Float_t Phi()        const { return fPhi; }
405     virtual Float_t Pt() const { return fPt; }
406     Int_t getparticle() const {return fparticle;}
407     virtual Short_t Charge()      const { return fcharge; }
408     Float_t geteffcorrectionval() const {return feffcorrectionval;}
409     virtual Bool_t IsEqual(const TObject* obj) const { return (obj->GetUniqueID() == GetUniqueID()); }
410
411
412 private:
413   LRCParticlePID(const LRCParticlePID&);  // not implemented
414    LRCParticlePID& operator=(const LRCParticlePID&);  // not implemented
415   
416   Int_t fparticle;
417   Short_t fcharge;
418   Float_t fPt;
419   Float_t fEta;
420   Float_t fPhi;
421   Float_t feffcorrectionval;
422   ClassDef(LRCParticlePID, 1);
423 } ;
424
425 #endif
426