AliPhysics  8417398 (8417398)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliAnalysisTaskSED0Mass.h
Go to the documentation of this file.
1 #ifndef ALIANALYSISTASKSED0MASS_H
2 #define ALIANALYSISTASKSED0MASS_H
3 
4 /* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice */
6 
7 /* $Id$ */
8 
17 
18 #include <TROOT.h>
19 #include <TSystem.h>
20 #include <TNtuple.h>
21 #include <TTree.h>
22 #include <TH1F.h>
23 #include <THnSparse.h>
24 
25 #include "AliAnalysisTaskSE.h"
26 #include "AliRDHFCutsD0toKpi.h"
28 
29 class AliAODEvent;
30 
31 class AliAnalysisTaskSED0Mass : public AliAnalysisTaskSE
32 {
33  public:
34 
36  AliAnalysisTaskSED0Mass(const char *name,AliRDHFCutsD0toKpi* cuts);
37  virtual ~AliAnalysisTaskSED0Mass();
38 
39 
41  virtual void UserCreateOutputObjects();
42  virtual void Init();
43  virtual void LocalInit() {Init();}
44  virtual void UserExec(Option_t *option);
45  virtual void Terminate(Option_t *option);
46 
48  Bool_t CheckAcc(TClonesArray* arrayMC,Int_t nProng, Int_t *labDau);
49  void FillMCAcceptanceHistos(TClonesArray *arrayMC, AliAODMCHeader *mcHeader);
50 
51  void NormIPvar(AliAODEvent *aod, AliAODRecoDecayHF2Prong *part,TClonesArray *arrMC);
53  enum{kD0,kLS};
54 
55  void SetReadMC(Bool_t readMC=kFALSE){fReadMC=readMC;}
56  void SetDoMCAcceptanceHistos(Bool_t doMCAcc=kTRUE){fStepMCAcc=doMCAcc;}
57  void SetCutOnDistr(Bool_t cutondistr=kFALSE){fCutOnDistr=cutondistr;}
58  void SetUsePid4Distr(Bool_t usepid=kTRUE){fUsePid4Distr=usepid;}
59  void SetFillOnlyD0D0bar(Int_t flagfill){fFillOnlyD0D0bar=flagfill;}
60  void SetFillVarHists(Bool_t flag) {fFillVarHists=flag;}
61  void SetFillPtHistos(Bool_t flag) {fFillPtHist=flag;}
62  void SetFillYHistos(Bool_t flag) {fFillYHist=flag;}
64  void SetSystem(Int_t sys){fSys=sys; if(fSys==1) SetFillVarHists(kFALSE);}
65  void SetRejectSDDClusters(Bool_t flag) { fIsRejectSDDClusters=flag; }
66  void SetUseSelectionBit(Bool_t flag) { fUseSelectionBit=flag; }
67  void SetWriteVariableTree(Bool_t flag) { fWriteVariableTree=flag; }
68  void SetDrawDetSignal(Bool_t flag) { fDrawDetSignal=flag; }
69  void SetPIDCheck(Bool_t flag) { fPIDCheck=flag; }
70  void SetUseQuarkLevelTag(Bool_t opt){fUseQuarkTagInKine=opt;}
71 
72 
73  Bool_t GetCutOnDistr() const {return fCutOnDistr;}
74  Bool_t GetUsePid4Distr() const {return fUsePid4Distr;}
75  Int_t GetFillOnlyD0D0bar() const {return fFillOnlyD0D0bar;}
76  Bool_t GetFillVarHists() const {return fFillVarHists;}
77  Bool_t GetFillPtHistos() const {return fFillPtHist;}
78  Bool_t GetFillYHistos() const {return fFillYHist;}
80  Int_t GetSystem() const {return fSys;}
81  Bool_t GetRejectSDDClusters() const { return fIsRejectSDDClusters; }
82  Bool_t GetUseSelectionBit() const { return fUseSelectionBit; }
83  Bool_t GetWriteVariableTree() const {return fWriteVariableTree;}
84  Bool_t GetDrawDetSignal() const {return fDrawDetSignal;}
85  Bool_t GetPIDCheck() const {return fPIDCheck;}
86 
87  private:
88 
91  void DrawDetSignal(AliAODRecoDecayHF2Prong *part, TList *ListDetSignal);
92 
93  void FillMassHists(AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliAODMCHeader *mcHeader, AliRDHFCutsD0toKpi *cuts, TList *listout);
94  void FillVarHists(AliAODEvent *aodev,AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliRDHFCutsD0toKpi *cuts, TList *listout);
95  AliAODVertex* GetPrimaryVtxSkipped(AliAODEvent *aodev);
97  Int_t CheckOrigin(TClonesArray* arrayMC, AliAODMCParticle *mcPartCandidate) const;
98  Float_t GetTrueImpactParameter(AliAODMCHeader *mcHeader, TClonesArray* arrayMC, AliAODMCParticle *partD0) const ;
99 
100  TList *fOutputMass;
101  TList *fOutputMassPt;
102  TList *fOutputMassY;
103  TList *fDistr;
104  TH1F *fNentries;
105  THnSparseF *fMCAccPrompt;
106  THnSparseF *fMCAccBFeed;
107  Bool_t fStepMCAcc; // flag to activate histos for StepMCAcc
108  AliRDHFCutsD0toKpi *fCuts; // Cuts - sent to output slot 4
109  THnSparseF *fHistMassPtImpParTC[5];
110  Int_t fArray;
111  Bool_t fReadMC;
112  Bool_t fCutOnDistr;
113  Bool_t fUsePid4Distr; // flag to use the particle identification to fill the signal histograms of distributions. It has effect only with fReadMC=kFALSE
115  Int_t fNPtBins;
116  Double_t fLsNormalization;
118  TObjArray fDaughterTracks;
120  Bool_t fFillVarHists;
121  Int_t fSys;
123  Bool_t fFillPtHist;
124  Bool_t fFillYHist;
127 
129  TTree *fVariablesTree;
131  Bool_t fPIDCheck;
132  Bool_t fDrawDetSignal;
133  Bool_t fUseQuarkTagInKine; // flag for quark/hadron level identification of prompt and feeddown
137  TList *fDetSignal;
138 
140  ClassDef(AliAnalysisTaskSED0Mass,20);
141 };
143 
144 #endif
145 
THnSparseF * fMCAccBFeed
!histo for StepMCAcc for D0 FD (pt,y,ptB)
TList * fOutputMassPt
! list send on output slot 6
Bool_t fDrawDetSignal
flag to decide whether to fill "PID = x" bins in fNentrie
void SetUseSelectionBit(Bool_t flag)
Bool_t fFillYHist
flag to fill Pt and Impact Parameter Histograms
void SetDoMCAcceptanceHistos(Bool_t doMCAcc=kTRUE)
void FillVarHists(AliAODEvent *aodev, AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliRDHFCutsD0toKpi *cuts, TList *listout)
Double_t * fCandidateVariables
! variables to be written to the tree
Bool_t fWriteVariableTree
flag to check or not the selection bit
Bool_t GetFillImpactParameterHistos() const
void SetArray(Int_t type=AliAnalysisTaskSED0Mass::kD0)
TList * fDistr
! list send on output slot 2
void NormIPvar(AliAODEvent *aod, AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC)
TList * fDetSignal
!Detector signal histograms (on output slot 8)
Bool_t fUseSelectionBit
flag to fill Pt and Impact Parameter Histograms
void SetFillOnlyD0D0bar(Int_t flagfill)
void SetUseQuarkLevelTag(Bool_t opt)
void FillMCAcceptanceHistos(TClonesArray *arrayMC, AliAODMCHeader *mcHeader)
Int_t fSys
flag to enable filling variable histos
void SetWriteVariableTree(Bool_t flag)
TList * fOutputMass
! list send on output slot 1
Bool_t CheckAcc(TClonesArray *arrayMC, Int_t nProng, Int_t *labDau)
Bool_t fReadMC
can be D0 or Like Sign candidates
void SetReadMC(Bool_t readMC=kFALSE)
THnSparseF * fHistMassPtImpParTC[5]
! histograms for impact paramter studies
AliAnalysisTaskSED0Mass & operator=(const AliAnalysisTaskSED0Mass &source)
TH1F * fNentries
! histogram with number of events on output slot 3
virtual void UserCreateOutputObjects()
Implementation of interface methods.
Bool_t fFillPtHist
flag to reject events with SDD clusters
Bool_t fFillImpParHist
flag to fill Y Histograms
Bool_t fUseQuarkTagInKine
flag to decide whether to draw the TPC dE/dx and TOF signal before/after PID
void SetRejectSDDClusters(Bool_t flag)
Bool_t fIsRejectSDDClusters
fSys=0 -> p-p; fSys=1 ->PbPb (in this case fFillVarHists=kFALSE by default: set it to kTRUE after if ...
void SetFillImpactParameterHistos(Bool_t flag)
Int_t CheckOrigin(TClonesArray *arrayMC, AliAODMCParticle *mcPartCandidate) const
virtual void Terminate(Option_t *option)
Float_t GetTrueImpactParameter(AliAODMCHeader *mcHeader, TClonesArray *arrayMC, AliAODMCParticle *partD0) const
void DrawDetSignal(AliAODRecoDecayHF2Prong *part, TList *ListDetSignal)
void SetDrawDetSignal(Bool_t flag)
Double_t fLsNormalization
number of pt bins
Bool_t fUsePid4Distr
flag to decide if apply cut also on distributions: 0 no cuts, 1 looser cuts, 2 tighter/ cuts ...
THnSparseF * fhStudyImpParSingleTrackFd
! sparse with imp par residual cuts for MC
Bool_t fCutOnDistr
flag for MC array: kTRUE = read it, kFALSE = do not read it
AliAODVertex * GetPrimaryVtxSkipped(AliAODEvent *aodev)
THnSparseF * fhStudyImpParSingleTrackSign
! sparse with imp par residual cuts for MC
void FillMassHists(AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliAODMCHeader *mcHeader, AliRDHFCutsD0toKpi *cuts, TList *listout)
Int_t fIsSelectedCandidate
keeps the daughter tracks
THnSparseF * fhStudyImpParSingleTrackCand
! sparse with imp par residual cuts for Data
TList * fOutputMassY
! list send on output slot 9
TTree * fVariablesTree
flag to decide whether to write the candidate variables on a tree variables
virtual void UserExec(Option_t *option)
TObjArray fDaughterTracks
flag to fill mass histogram with D0/D0bar only (0 = fill with both, 1 = fill with D0 only...
void SetCutOnDistr(Bool_t cutondistr=kFALSE)
void SetUsePid4Distr(Bool_t usepid=kTRUE)
Int_t fFillOnlyD0D0bar
normalization
THnSparseF * fMCAccPrompt
!histo for StepMCAcc for D0 prompt (pt,y,ptB)
Bool_t fFillVarHists
selection outcome
AliNormalizationCounter * fCounter
! AliNormalizationCounter on output slot 5