]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliAnalysisTaskCombinHF.h
Merge branch 'master' into dev
[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
38   void SetReadMC(Bool_t read){fReadMC=read;}
39   void SelectPromptD(){fPromptFeeddown=kPrompt;}
40   void SelectFeeddownD(){fPromptFeeddown=kFeeddown;}
41   void SelectPromptAndFeeddownD(){fPromptFeeddown=kBoth;}
42   void SetGoUpToQuark(Bool_t opt){fGoUpToQuark=opt;}
43   void SetKeepNegIDtracks(Bool_t nid){fKeepNegID=nid;}//set it to kTRUE only if you know what you are doing
44   void SetTrackCuts(AliESDtrackCuts* cuts){
45     if(fTrackCutsAll) delete fTrackCutsAll;
46     fTrackCutsAll=new AliESDtrackCuts(*cuts);
47   }
48   void SetPionTrackCuts(AliESDtrackCuts* cuts){
49     if(fTrackCutsPion) delete fTrackCutsPion;
50     fTrackCutsPion=new AliESDtrackCuts(*cuts);
51   }
52   void SetKaonTrackCuts(AliESDtrackCuts* cuts){
53     if(fTrackCutsKaon) delete fTrackCutsKaon;
54     fTrackCutsKaon=new AliESDtrackCuts(*cuts);
55   }
56   void SetPIDHF(AliAODPidHF* pid){
57     if(fPidHF) delete fPidHF;
58     fPidHF=new AliAODPidHF(*pid);
59   }
60   void SetRDHFCuts(AliRDHFCuts* cuts){
61     fAnalysisCuts=cuts;
62   }
63   void SetFilterMask(UInt_t mask=16){fFilterMask=mask;} 
64   void SetAnalysisLevel(Int_t level){fFullAnalysis=level;}
65   void ConfigureRotation(Int_t n, Double_t phimin, Double_t phimax){
66     fNRotations=n;
67     fMinAngleForRot=phimin;
68     fMaxAngleForRot=phimax;
69   }
70   Bool_t IsTrackSelected(AliAODTrack* track);
71   Bool_t IsKaon(AliAODTrack* track);
72   Bool_t IsPion(AliAODTrack* track);
73   Bool_t SelectAODTrack(AliAODTrack *track, AliESDtrackCuts *cuts);
74   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);
75   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);
76   void FillGenHistos(TClonesArray* arrayMC);
77   Bool_t CheckAcceptance(TClonesArray* arrayMC, Int_t nProng, Int_t *labDau);
78   enum EMesonSpecies {kDzero, kDplus, kDstar, kDs};
79   enum EPrompFd {kNone,kPrompt,kFeeddown,kBoth};
80   void SetMaxPforIDPion(Double_t maxpIdPion){fmaxPforIDPion=maxpIdPion;}
81   void SetMaxPforIDKaon(Double_t maxpIdKaon){fmaxPforIDKaon=maxpIdKaon;}
82
83  private:
84
85   AliAnalysisTaskCombinHF(const AliAnalysisTaskCombinHF &source);
86   AliAnalysisTaskCombinHF& operator=(const AliAnalysisTaskCombinHF& source); 
87
88   TList   *fOutput; //! list send on output slot 0
89   TH1F *fHistNEvents;         //!hist. for No. of events
90   TH1F *fHistTrackStatus;     //!hist. of status of tracks
91   TH1F *fHistCheckOrigin;     //!hist. of origin (c/b) of D meson
92   TH1F *fHistCheckOriginSel;  //!hist. of origin (c/b) of D meson
93   TH1F *fHistCheckDecChan;    //!hist. of decay channel of D meson
94   TH1F *fHistCheckDecChanAcc; //!hist. of decay channel of D meson in acc.
95   TH2F *fPtVsYGen;        //! hist. of Y vs. Pt generated (all D)
96   TH2F *fPtVsYGenLimAcc;  //! hist. of Y vs. Pt generated (|y|<0.5)
97   TH2F *fPtVsYGenAcc;     //! hist. of Y vs. Pt generated (D in acc)
98   TH2F *fPtVsYReco;       //! hist. of Y vs. Pt generated (Reco D)
99   TH3F *fMassVsPtVsY;     //! hist. of Y vs. Pt vs. Mass (all cand)
100   TH3F *fMassVsPtVsYRot;   //! hist. of Y vs. Pt vs. Mass (rotations)
101   TH3F *fMassVsPtVsYLSpp;  //! hist. of Y vs. Pt vs. Mass (like sign ++)
102   TH3F *fMassVsPtVsYLSmm;  //! hist. of Y vs. Pt vs. Mass (like sign --)
103   TH3F *fMassVsPtVsYSig;   //! hist. of Y vs. Pt vs. Mass (signal)
104   TH3F *fMassVsPtVsYRefl;  //! hist. of Y vs. Pt vs. Mass (reflections)
105   TH3F *fMassVsPtVsYBkg;   //! hist. of Y vs. Pt vs. Mass (background)
106   TH1F *fNSelected;        //! hist. of n. of selected D+
107   TH1F *fNormRotated;      //! hist. rotated/selected D+
108   TH1F *fDeltaMass;        //! hist. mass difference after rotations
109   THnSparse *fDeltaMassFullAnalysis; //! hist. mass difference after rotations with more details
110
111   UInt_t fFilterMask; // FilterMask
112   AliESDtrackCuts* fTrackCutsAll; // track selection
113   AliESDtrackCuts* fTrackCutsPion; // pion track selection
114   AliESDtrackCuts* fTrackCutsKaon; // kaon track selection
115   AliAODPidHF* fPidHF; // PID configuration
116   AliRDHFCuts *fAnalysisCuts; // Cuts for candidates
117   Double_t fMinMass; // minimum value of invariant mass
118   Double_t fMaxMass; // maximum value of invariant mass
119
120   Double_t fEtaAccCut; // eta limits for acceptance step
121   Double_t fPtAccCut; // pt limits for acceptance step
122
123
124   Int_t fNRotations; // number of rotations
125   Double_t fMinAngleForRot; // minimum angle for track rotation
126   Double_t fMaxAngleForRot; // maximum angle for track rotation
127   Double_t fMinAngleForRot3; // minimum angle for track rotation (3rd prong)
128   Double_t fMaxAngleForRot3; // maximum angle for track rotation (3rd prong)
129
130   AliNormalizationCounter *fCounter;//!Counter for normalization
131
132   Int_t fMeson;          // mesonSpecies (see enum)
133   Bool_t  fReadMC;       //  flag for access to MC
134   Int_t fPromptFeeddown; // flag to select prompt (1), feeddown (2) or all (3) 
135   Bool_t fGoUpToQuark;   // flag for definition of c,b origin
136   Int_t fFullAnalysis;   // flag to set analysis level (0 is the fastest)
137
138   Double_t fmaxPforIDPion; // flag for upper p limit for id band for pion
139   Double_t fmaxPforIDKaon; // flag for upper p limit for id band for kaon
140   Bool_t   fKeepNegID;    // flag to keep also track with negative ID (default kFALSE, change it only if you know what you are doing)
141   ClassDef(AliAnalysisTaskCombinHF,3); // D+ task from AOD tracks
142 };
143
144 #endif