]>
Commit | Line | Data |
---|---|---|
d486095a | 1 | #ifndef ALIANALYSISTASKDPLUS_H |
2 | #define ALIANALYSISTASKDPLUS_H | |
3 | ||
4 | /* Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. * | |
5 | * See cxx source for full Copyright notice */ | |
6 | ||
7 | //************************************************************************* | |
8 | // Class AliAnalysisTaskSEDplus | |
4afc48a2 | 9 | // AliAnalysisTaskSE for the D+ candidates Invariant Mass Histogram and |
10 | //comparison of heavy-flavour decay candidates | |
d486095a | 11 | // to MC truth (kinematics stored in the AOD) |
4afc48a2 | 12 | // Renu Bala, bala@to.infn.it |
993bfba1 | 13 | // F. Prino, prino@to.infn.it |
14 | // G. Ortona, ortona@to.infn.it | |
d486095a | 15 | //************************************************************************* |
16 | ||
17 | #include <TROOT.h> | |
18 | #include <TSystem.h> | |
19 | #include <TNtuple.h> | |
20 | #include <TH1F.h> | |
595cc7e2 | 21 | #include <TH2F.h> |
82bb8d4b | 22 | #include <TArrayD.h> |
d486095a | 23 | |
4c230f6d | 24 | #include "AliRDHFCutsDplustoKpipi.h" |
d486095a | 25 | #include "AliAnalysisTaskSE.h" |
26 | #include "AliAnalysisVertexingHF.h" | |
a96083b9 | 27 | #include "AliNormalizationCounter.h" |
d486095a | 28 | |
29 | class AliAnalysisTaskSEDplus : public AliAnalysisTaskSE | |
30 | { | |
31 | public: | |
32 | ||
33 | AliAnalysisTaskSEDplus(); | |
4c230f6d | 34 | AliAnalysisTaskSEDplus(const char *name, AliRDHFCutsDplustoKpipi* analysiscuts,AliRDHFCutsDplustoKpipi* productioncuts,Bool_t fillNtuple=kFALSE); |
d486095a | 35 | virtual ~AliAnalysisTaskSEDplus(); |
82bb8d4b | 36 | |
4afc48a2 | 37 | void SetReadMC(Bool_t readMC=kTRUE){fReadMC=readMC;} |
38 | void SetDoLikeSign(Bool_t dols=kTRUE){fDoLS=dols;} | |
a96083b9 | 39 | void SetUseStrangeness(Bool_t uses=kTRUE){fUseStrangeness=uses;} |
82bb8d4b | 40 | void SetMassLimits(Float_t range); |
41 | void SetMassLimits(Float_t lowlimit, Float_t uplimit); | |
4c230f6d | 42 | void SetPtBinLimit(Int_t n, Float_t *limitarray); |
993bfba1 | 43 | void SetBinWidth(Float_t w); |
82bb8d4b | 44 | |
45 | Float_t GetUpperMassLimit(){return fUpmasslimit;} | |
46 | Float_t GetLowerMassLimit(){return fLowmasslimit;} | |
47 | Int_t GetNBinsPt(){return fNPtBins;} | |
48 | Double_t GetPtBinLimit(Int_t ibin); | |
993bfba1 | 49 | Float_t GetBinWidth(){return fBinWidth;} |
50 | Int_t GetNBinsHistos(); | |
82bb8d4b | 51 | |
52 | void LSAnalysis(TClonesArray *arrayOppositeSign,TClonesArray *arrayLikeSign,AliAODEvent *aod,AliAODVertex *vtx1, Int_t nDplusOS); | |
53 | ||
d486095a | 54 | // Implementation of interface methods |
55 | virtual void UserCreateOutputObjects(); | |
56 | virtual void Init(); | |
57 | virtual void LocalInit() {Init();} | |
58 | virtual void UserExec(Option_t *option); | |
59 | virtual void Terminate(Option_t *option); | |
4afc48a2 | 60 | |
d486095a | 61 | private: |
62 | ||
63 | AliAnalysisTaskSEDplus(const AliAnalysisTaskSEDplus &source); | |
64 | AliAnalysisTaskSEDplus& operator=(const AliAnalysisTaskSEDplus& source); | |
82bb8d4b | 65 | Int_t GetHistoIndex(Int_t iPtBin) const { return iPtBin*3;} |
66 | Int_t GetSignalHistoIndex(Int_t iPtBin) const { return iPtBin*3+1;} | |
67 | Int_t GetBackgroundHistoIndex(Int_t iPtBin) const { return iPtBin*3+2;} | |
68 | Int_t GetLSHistoIndex(Int_t iPtBin)const { return iPtBin*5;} | |
69 | ||
73173a6a | 70 | enum {kMaxPtBins=20}; |
82bb8d4b | 71 | |
d486095a | 72 | TList *fOutput; //! list send on output slot 0 |
4afc48a2 | 73 | TH1F *fHistNEvents; //!hist. for No. of events |
82bb8d4b | 74 | TH1F *fMassHist[3*kMaxPtBins]; //!hist. for inv mass (LC) |
ba9ae5b2 | 75 | TH1F* fCosPHist[3*kMaxPtBins]; //!hist. for PointingAngle (LC) |
76 | TH1F* fDLenHist[3*kMaxPtBins]; //!hist. for Dec Length (LC) | |
77 | TH1F* fSumd02Hist[3*kMaxPtBins]; //!hist. for sum d02 (LC) | |
78 | TH1F* fSigVertHist[3*kMaxPtBins]; //!hist. for sigVert (LC) | |
79 | TH1F* fPtMaxHist[3*kMaxPtBins]; //!hist. for Pt Max (LC) | |
80 | TH1F* fDCAHist[3*kMaxPtBins]; //!hist. for DCA (LC) | |
82bb8d4b | 81 | TH1F *fMassHistTC[3*kMaxPtBins]; //!hist. for inv mass (TC) |
bb69f127 | 82 | TH1F *fMassHistTCPlus[3*kMaxPtBins]; //!hist. for D+ inv mass (TC) |
83 | TH1F *fMassHistTCMinus[3*kMaxPtBins]; //!hist. for D- inv mass (TC) | |
82bb8d4b | 84 | TH1F *fMassHistLS[5*kMaxPtBins];//!hist. for LS inv mass (LC) |
ba9ae5b2 | 85 | TH1F *fCosPHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 1 (LC) |
86 | TH1F *fDLenHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 2 (LC) | |
87 | TH1F *fSumd02HistLS[3*kMaxPtBins];//!hist. for LS cuts variable 3 (LC) | |
88 | TH1F *fSigVertHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 4 (LC) | |
89 | TH1F *fPtMaxHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 5 (LC) | |
90 | TH1F *fDCAHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 6 (LC) | |
82bb8d4b | 91 | TH1F *fMassHistLSTC[5*kMaxPtBins];//!hist. for LS inv mass (TC) |
595cc7e2 | 92 | TH2F *fPtVsMass; //! hist. of pt vs. mass (prod. cuts) |
93 | TH2F *fPtVsMassTC; //! hist. of pt vs. mass (analysis cuts) | |
94 | TH2F *fYVsPt; //! hist. of Y vs. Pt (prod. cuts) | |
95 | TH2F *fYVsPtTC; //! hist. of Y vs. Pt (analysis cuts) | |
96 | TH2F *fYVsPtSig; //! hist. of Y vs. Pt (MC, only sig, prod. cuts) | |
97 | TH2F *fYVsPtSigTC; //! hist. of Y vs. Pt (MC, only sig, analysis cuts) | |
1f4e9722 | 98 | TNtuple *fNtupleDplus; //! output ntuple |
82bb8d4b | 99 | Float_t fUpmasslimit; //upper inv mass limit for histos |
100 | Float_t fLowmasslimit; //lower inv mass limit for histos | |
4c230f6d | 101 | Int_t fNPtBins; //Number of Pt Bins |
993bfba1 | 102 | Float_t fBinWidth;//width of one bin in output histos |
4c230f6d | 103 | TList *fListCuts; //list of cuts |
104 | AliRDHFCutsDplustoKpipi *fRDCutsProduction; //Production D+ Cuts | |
105 | AliRDHFCutsDplustoKpipi *fRDCutsAnalysis; //Cuts for Analysis | |
a96083b9 | 106 | AliNormalizationCounter *fCounter;//!Counter for normalization |
82bb8d4b | 107 | Double_t fArrayBinLimits[kMaxPtBins+1]; //limits for the Pt bins |
1f4e9722 | 108 | Bool_t fFillNtuple; // flag for filling ntuple |
82bb8d4b | 109 | Bool_t fReadMC; //flag for access to MC |
a96083b9 | 110 | Bool_t fUseStrangeness;//flag to enhance strangeness in MC to fit to data |
82bb8d4b | 111 | Bool_t fDoLS; //flag to do LS analysis |
d486095a | 112 | |
a96083b9 | 113 | ClassDef(AliAnalysisTaskSEDplus,9); // AliAnalysisTaskSE for the MC association of heavy-flavour decay candidates |
d486095a | 114 | }; |
115 | ||
116 | #endif |