]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliAnalysisTaskCombinHF.h
Merge branch 'master' into TPCdev
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / AliAnalysisTaskCombinHF.h
1 #ifndef ALIANALYSISTASKCOMBINHF_H
2 #define ALIANALYSISTASKCOMBINHF_H
3
4 /* Copyright(c) 1998-2018, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 /* $Id: $ */
8
9 //*************************************************************************
10 // Class AliAnalysisTaskCombinHF
11 // AliAnalysisTaskSE to build D meson candidates by combining tracks
12 //  background is computed LS and track rotations is
13 // Authors: F. Prino, A. Rossi
14 /////////////////////////////////////////////////////////////
15
16 #include <TH1F.h>
17 #include <TH3F.h>
18 #include <THnSparse.h>
19 #include "AliAnalysisTaskSE.h"
20 #include "AliAODTrack.h"
21 #include "AliNormalizationCounter.h"
22 #include "AliRDHFCuts.h"
23
24 class AliAnalysisTaskCombinHF : public AliAnalysisTaskSE
25 {
26 public:
27   
28   AliAnalysisTaskCombinHF();
29   AliAnalysisTaskCombinHF(Int_t meson, AliRDHFCuts* analysiscuts);
30   virtual ~AliAnalysisTaskCombinHF();
31   
32   virtual void UserCreateOutputObjects();
33   virtual void Init(){};
34   virtual void LocalInit() {Init();}
35   virtual void UserExec(Option_t *option);
36   virtual void Terminate(Option_t *option);
37   virtual void FinishTaskOutput();
38   
39   void SetReadMC(Bool_t read){fReadMC=read;}
40
41   void SetEventMixingOn(){fDoEventMixing=kTRUE;}
42   void SetEventMixingOff(){fDoEventMixing=kFALSE;}
43   void SetMaxNumberOfEventsForMixing(Int_t maxn){fMaxNumberOfEventsForMixing=maxn;}
44   void SetMaxzVertDistForMix(Double_t dist){
45     fUsePoolsZ=kFALSE; 
46     fMaxzVertDistForMix=dist;
47   }
48   void SetMaxMultDiffForMix(Double_t maxd){
49     fUsePoolsM=kFALSE; 
50     fMaxMultDiffForMix=maxd;
51   }
52   void ConfigureZVertPools(Int_t nPools, Double_t*  zVertLimits);
53   void ConfigureMultiplicityPools(Int_t nPools, Double_t*  multLimits);
54   void SelectPromptD(){fPromptFeeddown=kPrompt;}
55   void SelectFeeddownD(){fPromptFeeddown=kFeeddown;}
56   void SelectPromptAndFeeddownD(){fPromptFeeddown=kBoth;}
57   void SetGoUpToQuark(Bool_t opt){fGoUpToQuark=opt;}
58   void SetKeepNegIDtracks(Bool_t nid){fKeepNegID=nid;}//set it to kTRUE only if you know what you are doing
59   void SetTrackCuts(AliESDtrackCuts* cuts){
60     if(fTrackCutsAll) delete fTrackCutsAll;
61     fTrackCutsAll=new AliESDtrackCuts(*cuts);
62   }
63   void SetPionTrackCuts(AliESDtrackCuts* cuts){
64     if(fTrackCutsPion) delete fTrackCutsPion;
65     fTrackCutsPion=new AliESDtrackCuts(*cuts);
66   }
67   void SetKaonTrackCuts(AliESDtrackCuts* cuts){
68     if(fTrackCutsKaon) delete fTrackCutsKaon;
69     fTrackCutsKaon=new AliESDtrackCuts(*cuts);
70   }
71   void SetPIDHF(AliAODPidHF* pid){
72     if(fPidHF) delete fPidHF;
73     fPidHF=new AliAODPidHF(*pid);
74   }
75   void SetRDHFCuts(AliRDHFCuts* cuts){
76     fAnalysisCuts=cuts;
77   }
78   void SetFilterMask(UInt_t mask=16){fFilterMask=mask;}
79   void SetAnalysisLevel(Int_t level){fFullAnalysis=level;}
80   void ConfigureRotation(Int_t n, Double_t phimin, Double_t phimax){
81     fNRotations=n;
82     fMinAngleForRot=phimin;
83     fMaxAngleForRot=phimax;
84   }
85   void SetMassWindow(Double_t minMass, Double_t maxMass){fMinMass=minMass; fMaxMass=maxMass;}
86   void SetMaxPt(Double_t maxPt){fMaxPt=maxPt;}
87   void SetPtBinWidth(Double_t binw){fPtBinWidth=binw;}
88   void SetEtaAccCut(Double_t etacut){fEtaAccCut=etacut;}
89   void SetPtAccCut(Double_t ptcut){fPtAccCut=ptcut;}
90   
91   void SetPIDstrategy(Int_t strat){fPIDstrategy=strat;}
92   void SetMaxPforIDPion(Double_t maxpIdPion){fmaxPforIDPion=maxpIdPion;}
93   void SetMaxPforIDKaon(Double_t maxpIdKaon){fmaxPforIDKaon=maxpIdKaon;}
94   void SetPIDselCaseZero(Int_t strat){fPIDselCaseZero=strat;}
95   void SetBayesThres(Double_t thresKaon, Double_t thresPion){
96     fBayesThresKaon=thresKaon;
97     fBayesThresPion=thresPion;
98   }
99   
100   Bool_t IsTrackSelected(AliAODTrack* track);
101   Bool_t IsKaon(AliAODTrack* track);
102   Bool_t IsPion(AliAODTrack* track);
103   Bool_t SelectAODTrack(AliAODTrack *track, AliESDtrackCuts *cuts);
104   
105   Bool_t FillHistos(Int_t pdgD,Int_t nProngs, AliAODRecoDecay* tmpRD, Double_t* px, Double_t* py, Double_t* pz, UInt_t *pdgdau, TClonesArray *arrayMC, Int_t* dgLabels);
106   void FillLSHistos(Int_t pdgD,Int_t nProngs, AliAODRecoDecay* tmpRD, Double_t* px, Double_t* py, Double_t* pz, UInt_t *pdgdau, Int_t charge);
107   void FillMEHistos(Int_t pdgD,Int_t nProngs, AliAODRecoDecay* tmpRD, Double_t* px, Double_t* py, Double_t* pz, UInt_t *pdgdau);
108   void FillGenHistos(TClonesArray* arrayMC);
109   Bool_t CheckAcceptance(TClonesArray* arrayMC, Int_t nProng, Int_t *labDau); 
110   Bool_t CanBeMixed(Double_t zv1, Double_t zv2, Double_t mult1, Double_t mult2); 
111
112   enum EMesonSpecies {kDzero, kDplus, kDstar, kDs};
113   enum EPrompFd {kNone,kPrompt,kFeeddown,kBoth};
114   enum EPIDstrategy {knSigma, kBayesianMaxProb, kBayesianThres};
115   
116 private:
117   
118   AliAnalysisTaskCombinHF(const AliAnalysisTaskCombinHF &source);
119   AliAnalysisTaskCombinHF& operator=(const AliAnalysisTaskCombinHF& source);
120   
121   TList   *fOutput; //! list send on output slot 0
122   TH1F *fHistNEvents;         //!hist. for No. of events
123   TH1F *fHistTrackStatus;     //!hist. of status of tracks
124   TH1F *fHistCheckOrigin;     //!hist. of origin (c/b) of D meson
125   TH1F *fHistCheckOriginSel;  //!hist. of origin (c/b) of D meson
126   TH1F *fHistCheckDecChan;    //!hist. of decay channel of D meson
127   TH1F *fHistCheckDecChanAcc; //!hist. of decay channel of D meson in acc.
128   TH2F *fPtVsYGen;        //! hist. of Y vs. Pt generated (all D)
129   TH2F *fPtVsYGenLargeAcc; //! hist. of Y vs. Pt generated (|y|<0.9)
130   TH2F *fPtVsYGenLimAcc;  //! hist. of Y vs. Pt generated (|y|<0.5)
131   TH2F *fPtVsYGenAcc;     //! hist. of Y vs. Pt generated (D in acc)
132   TH2F *fPtVsYReco;       //! hist. of Y vs. Pt generated (Reco D)
133   TH3F *fMassVsPtVsY;     //! hist. of Y vs. Pt vs. Mass (all cand)
134   TH3F *fMassVsPtVsYRot;   //! hist. of Y vs. Pt vs. Mass (rotations)
135   TH3F *fMassVsPtVsYLSpp;  //! hist. of Y vs. Pt vs. Mass (like sign ++)
136   TH3F *fMassVsPtVsYLSmm;  //! hist. of Y vs. Pt vs. Mass (like sign --)
137   TH3F *fMassVsPtVsYSig;   //! hist. of Y vs. Pt vs. Mass (signal)
138   TH3F *fMassVsPtVsYRefl;  //! hist. of Y vs. Pt vs. Mass (reflections)
139   TH3F *fMassVsPtVsYBkg;   //! hist. of Y vs. Pt vs. Mass (background)
140   TH1F *fNSelected;        //! hist. of n. of selected D+
141   TH1F *fNormRotated;      //! hist. rotated/selected D+
142   TH1F *fDeltaMass;        //! hist. mass difference after rotations
143   THnSparse *fDeltaMassFullAnalysis; //! hist. mass difference after rotations with more details
144   TH3F *fMassVsPtVsYME;   //! hist. of Y vs. Pt vs. Mass (mixedevents)
145   
146   UInt_t fFilterMask; // FilterMask
147   AliESDtrackCuts* fTrackCutsAll; // track selection
148   AliESDtrackCuts* fTrackCutsPion; // pion track selection
149   AliESDtrackCuts* fTrackCutsKaon; // kaon track selection
150   AliAODPidHF* fPidHF; // PID configuration
151   AliRDHFCuts *fAnalysisCuts; // Cuts for candidates
152   
153   Double_t fMinMass; // minimum value of invariant mass
154   Double_t fMaxMass; // maximum value of invariant mass
155   Double_t fMaxPt;   // maximum pT value for inv. mass histograms
156   Double_t fPtBinWidth; // width of pt bin (GeV/c)
157   Double_t fEtaAccCut; // eta limits for acceptance step
158   Double_t fPtAccCut; // pt limits for acceptance step
159   
160   Int_t fNRotations; // number of rotations
161   Double_t fMinAngleForRot; // minimum angle for track rotation
162   Double_t fMaxAngleForRot; // maximum angle for track rotation
163   Double_t fMinAngleForRot3; // minimum angle for track rotation (3rd prong)
164   Double_t fMaxAngleForRot3; // maximum angle for track rotation (3rd prong)
165   
166   AliNormalizationCounter *fCounter;//!Counter for normalization
167   
168   Int_t fMeson;          // mesonSpecies (see enum)
169   Bool_t  fReadMC;       //  flag for access to MC
170   Int_t fPromptFeeddown; // flag to select prompt (1), feeddown (2) or all (3)
171   Bool_t fGoUpToQuark;   // flag for definition of c,b origin
172   Int_t fFullAnalysis;   // flag to set analysis level (0 is the fastest)
173   
174   Int_t    fPIDstrategy;     // knSigma, kBayesianMaxProb, kBayesianThres
175   Double_t fmaxPforIDPion; // flag for upper p limit for id band for pion
176   Double_t fmaxPforIDKaon; // flag for upper p limit for id band for kaon
177   Bool_t   fKeepNegID;    // flag to keep also track with negative ID (default kFALSE, change it only if you know what you are doing)
178   Int_t    fPIDselCaseZero;  // flag to change PID strategy
179   Double_t fBayesThresKaon;  // threshold for kaon identification via Bayesian PID
180   Double_t fBayesThresPion;  // threshold for pion identification via Bayesian PID
181
182   Bool_t fDoEventMixing; // flag for event mixing
183   Int_t  fMaxNumberOfEventsForMixing; // maximum number of events to be used in event mixing
184   Double_t fMaxzVertDistForMix; // cut on zvertex distance for event mixing
185   Double_t fMaxMultDiffForMix; // cut on multiplicity difference for event mixing 
186   Bool_t fUsePoolsZ; // flag for using pools in z vertex
187   Bool_t fUsePoolsM; // flag for using pools in multiplicity
188   Int_t fNzVertPools; // number of pools in z vertex for event mixing
189   Int_t fNzVertPoolsLimSize; // number of pools in z vertex for event mixing +1
190   Double_t* fzVertPoolLims; //[fNzVertPoolsLimSize] limits of the pools in zVertex
191   Int_t fNMultPools; // number of pools in multiplicity for event mixing
192   Int_t fNMultPoolsLimSize; // number of pools in multiplicity for event mixing +1
193   Double_t* fMultPoolLims; //[fNMultPoolsLimSize] limits of the pools in multiplicity
194
195   TTree* fEventBuffer;    // structure for event mixing
196   Double_t fVtxZ;         // zVertex
197   Double_t fMultiplicity; // multiplicity
198   TObjArray* fKaonTracks; // array of kaon-compatible tracks (TLorentzVectors)
199   TObjArray* fPionTracks; // array of pion-compatible tracks (TLorentzVectors)  
200   ClassDef(AliAnalysisTaskCombinHF,7); // D0D+ task from AOD tracks
201 };
202
203 #endif