]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliAnalysisTaskCombinHF.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[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 SetMinNumberOfEventsForMixing(Int_t minn){fMinNumberOfEventsForMixing=minn;}
44
45   void ConfigureZVertPools(Int_t nPools, Double_t*  zVertLimits);
46   void ConfigureMultiplicityPools(Int_t nPools, Double_t*  multLimits);
47   void SelectPromptD(){fPromptFeeddown=kPrompt;}
48   void SelectFeeddownD(){fPromptFeeddown=kFeeddown;}
49   void SelectPromptAndFeeddownD(){fPromptFeeddown=kBoth;}
50   void SetGoUpToQuark(Bool_t opt){fGoUpToQuark=opt;}
51   void SetKeepNegIDtracks(Bool_t nid){fKeepNegID=nid;}//set it to kTRUE only if you know what you are doing
52   void SetTrackCuts(AliESDtrackCuts* cuts){
53     if(fTrackCutsAll) delete fTrackCutsAll;
54     fTrackCutsAll=new AliESDtrackCuts(*cuts);
55   }
56   void SetPionTrackCuts(AliESDtrackCuts* cuts){
57     if(fTrackCutsPion) delete fTrackCutsPion;
58     fTrackCutsPion=new AliESDtrackCuts(*cuts);
59   }
60   void SetKaonTrackCuts(AliESDtrackCuts* cuts){
61     if(fTrackCutsKaon) delete fTrackCutsKaon;
62     fTrackCutsKaon=new AliESDtrackCuts(*cuts);
63   }
64   void SetPIDHF(AliAODPidHF* pid){
65     if(fPidHF) delete fPidHF;
66     fPidHF=new AliAODPidHF(*pid);
67   }
68   void SetRDHFCuts(AliRDHFCuts* cuts){
69     fAnalysisCuts=cuts;
70   }
71   void SetFilterMask(UInt_t mask=16){fFilterMask=mask;}
72   void SetAnalysisLevel(Int_t level){fFullAnalysis=level;}
73   void ConfigureRotation(Int_t n, Double_t phimin, Double_t phimax){
74     fNRotations=n;
75     fMinAngleForRot=phimin;
76     fMaxAngleForRot=phimax;
77   }
78   void SetMassWindow(Double_t minMass, Double_t maxMass){fMinMass=minMass; fMaxMass=maxMass;}
79   void SetMaxPt(Double_t maxPt){fMaxPt=maxPt;}
80   void SetPtBinWidth(Double_t binw){fPtBinWidth=binw;}
81   void SetEtaAccCut(Double_t etacut){fEtaAccCut=etacut;}
82   void SetPtAccCut(Double_t ptcut){fPtAccCut=ptcut;}
83   
84   void SetPIDstrategy(Int_t strat){fPIDstrategy=strat;}
85   void SetMaxPforIDPion(Double_t maxpIdPion){fmaxPforIDPion=maxpIdPion;}
86   void SetMaxPforIDKaon(Double_t maxpIdKaon){fmaxPforIDKaon=maxpIdKaon;}
87   void SetPIDselCaseZero(Int_t strat){fPIDselCaseZero=strat;}
88   void SetBayesThres(Double_t thresKaon, Double_t thresPion){
89     fBayesThresKaon=thresKaon;
90     fBayesThresPion=thresPion;
91   }
92   
93   Bool_t IsTrackSelected(AliAODTrack* track);
94   Bool_t IsKaon(AliAODTrack* track);
95   Bool_t IsPion(AliAODTrack* track);
96   Bool_t SelectAODTrack(AliAODTrack *track, AliESDtrackCuts *cuts);
97   
98   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);
99   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);
100   void FillMEHistos(Int_t pdgD,Int_t nProngs, AliAODRecoDecay* tmpRD, Double_t* px, Double_t* py, Double_t* pz, UInt_t *pdgdau);
101   void FillGenHistos(TClonesArray* arrayMC);
102   Bool_t CheckAcceptance(TClonesArray* arrayMC, Int_t nProng, Int_t *labDau); 
103   Int_t GetPoolIndex(Double_t zvert, Double_t mult);
104   void ResetPool(Int_t poolIndex);
105   void DoMixing(Int_t poolIndex);
106
107   enum EMesonSpecies {kDzero, kDplus, kDstar, kDs};
108   enum EPrompFd {kNone,kPrompt,kFeeddown,kBoth};
109   enum EPIDstrategy {knSigma, kBayesianMaxProb, kBayesianThres};
110   
111 private:
112   
113   AliAnalysisTaskCombinHF(const AliAnalysisTaskCombinHF &source);
114   AliAnalysisTaskCombinHF& operator=(const AliAnalysisTaskCombinHF& source);
115   
116   TList   *fOutput; //! list send on output slot 0
117   TH1F *fHistNEvents;         //!hist. for No. of events
118   TH1F *fHistTrackStatus;     //!hist. of status of tracks
119   TH1F *fHistCheckOrigin;     //!hist. of origin (c/b) of D meson
120   TH1F *fHistCheckOriginSel;  //!hist. of origin (c/b) of D meson
121   TH1F *fHistCheckDecChan;    //!hist. of decay channel of D meson
122   TH1F *fHistCheckDecChanAcc; //!hist. of decay channel of D meson in acc.
123   TH2F *fPtVsYGen;        //! hist. of Y vs. Pt generated (all D)
124   TH2F *fPtVsYGenLargeAcc; //! hist. of Y vs. Pt generated (|y|<0.9)
125   TH2F *fPtVsYGenLimAcc;  //! hist. of Y vs. Pt generated (|y|<0.5)
126   TH2F *fPtVsYGenAcc;     //! hist. of Y vs. Pt generated (D in acc)
127   TH2F *fPtVsYReco;       //! hist. of Y vs. Pt generated (Reco D)
128   TH3F *fMassVsPtVsY;     //! hist. of Y vs. Pt vs. Mass (all cand)
129   TH3F *fMassVsPtVsYRot;   //! hist. of Y vs. Pt vs. Mass (rotations)
130   TH3F *fMassVsPtVsYLSpp;  //! hist. of Y vs. Pt vs. Mass (like sign ++)
131   TH3F *fMassVsPtVsYLSmm;  //! hist. of Y vs. Pt vs. Mass (like sign --)
132   TH3F *fMassVsPtVsYSig;   //! hist. of Y vs. Pt vs. Mass (signal)
133   TH3F *fMassVsPtVsYRefl;  //! hist. of Y vs. Pt vs. Mass (reflections)
134   TH3F *fMassVsPtVsYBkg;   //! hist. of Y vs. Pt vs. Mass (background)
135   TH1F *fNSelected;        //! hist. of n. of selected D+
136   TH1F *fNormRotated;      //! hist. rotated/selected D+
137   TH1F *fDeltaMass;        //! hist. mass difference after rotations
138   THnSparse *fDeltaMassFullAnalysis; //! hist. mass difference after rotations with more details
139   TH3F *fMassVsPtVsYME;   //! hist. of Y vs. Pt vs. Mass (mixedevents)
140   TH2F* fEventsPerPool;   //! hist with number of events per pool  
141   UInt_t fFilterMask; // FilterMask
142   AliESDtrackCuts* fTrackCutsAll; // track selection
143   AliESDtrackCuts* fTrackCutsPion; // pion track selection
144   AliESDtrackCuts* fTrackCutsKaon; // kaon track selection
145   AliAODPidHF* fPidHF; // PID configuration
146   AliRDHFCuts *fAnalysisCuts; // Cuts for candidates
147   
148   Double_t fMinMass; // minimum value of invariant mass
149   Double_t fMaxMass; // maximum value of invariant mass
150   Double_t fMaxPt;   // maximum pT value for inv. mass histograms
151   Double_t fPtBinWidth; // width of pt bin (GeV/c)
152   Double_t fEtaAccCut; // eta limits for acceptance step
153   Double_t fPtAccCut; // pt limits for acceptance step
154   
155   Int_t fNRotations; // number of rotations
156   Double_t fMinAngleForRot; // minimum angle for track rotation
157   Double_t fMaxAngleForRot; // maximum angle for track rotation
158   Double_t fMinAngleForRot3; // minimum angle for track rotation (3rd prong)
159   Double_t fMaxAngleForRot3; // maximum angle for track rotation (3rd prong)
160   
161   AliNormalizationCounter *fCounter;//!Counter for normalization
162   
163   Int_t fMeson;          // mesonSpecies (see enum)
164   Bool_t  fReadMC;       //  flag for access to MC
165   Int_t fPromptFeeddown; // flag to select prompt (1), feeddown (2) or all (3)
166   Bool_t fGoUpToQuark;   // flag for definition of c,b origin
167   Int_t fFullAnalysis;   // flag to set analysis level (0 is the fastest)
168   
169   Int_t    fPIDstrategy;     // knSigma, kBayesianMaxProb, kBayesianThres
170   Double_t fmaxPforIDPion; // flag for upper p limit for id band for pion
171   Double_t fmaxPforIDKaon; // flag for upper p limit for id band for kaon
172   Bool_t   fKeepNegID;    // flag to keep also track with negative ID (default kFALSE, change it only if you know what you are doing)
173   Int_t    fPIDselCaseZero;  // flag to change PID strategy
174   Double_t fBayesThresKaon;  // threshold for kaon identification via Bayesian PID
175   Double_t fBayesThresPion;  // threshold for pion identification via Bayesian PID
176
177   Bool_t fDoEventMixing; // flag for event mixing
178   Int_t  fMinNumberOfEventsForMixing; // maximum number of events to be used in event mixing
179   Int_t fNzVertPools; // number of pools in z vertex for event mixing
180   Int_t fNzVertPoolsLimSize; // number of pools in z vertex for event mixing +1
181   Double_t* fzVertPoolLims; //[fNzVertPoolsLimSize] limits of the pools in zVertex
182   Int_t fNMultPools; // number of pools in multiplicity for event mixing
183   Int_t fNMultPoolsLimSize; // number of pools in multiplicity for event mixing +1
184   Double_t* fMultPoolLims; //[fNMultPoolsLimSize] limits of the pools in multiplicity
185
186   TTree** fEventBuffer;   //! structure for event mixing
187   Double_t fVtxZ;         // zVertex
188   Double_t fMultiplicity; // multiplicity
189   TObjArray* fKaonTracks; // array of kaon-compatible tracks (TLorentzVectors)
190   TObjArray* fPionTracks; // array of pion-compatible tracks (TLorentzVectors)  
191   ClassDef(AliAnalysisTaskCombinHF,8); // D0D+ task from AOD tracks
192 };
193
194 #endif