]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/AliAnalysisTaskSECharmFraction.h
Coveriry
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliAnalysisTaskSECharmFraction.h
1 #ifndef ALIANALYSISTASKSECHARMFRACTION_H
2 #define ALIANALYSISTASKSECHARMFRACTION_H
3
4 /* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 //*************************************************************************
8 // Class AliAnalysisTaskSECharmFraction
9 // AliAnalysisTask for the extraction of the fraction of prompt charm
10 // using the charm hadron impact parameter to the primary vertex
11 //
12 //
13 // Author: Andrea Rossi andrea.rossi@pd.infn.it
14 //*************************************************************************
15
16 class TH1F;
17 class TH2F;
18 class AliAODDEvent;
19 class AliAODMCHeader;
20 class AliAODRecoDecayHF2Prong;
21 class AliAODRecoDecayHF;
22 class AliAODMCParticle;
23 class AliAnalysisVertexingHF;
24 class AliRDHFCutsD0toKpi;
25 class AliNormalizationCounter;
26
27 #include "AliAnalysisTaskSE.h"
28
29 class AliAnalysisTaskSECharmFraction : public AliAnalysisTaskSE {
30  public:
31   AliAnalysisTaskSECharmFraction();
32   AliAnalysisTaskSECharmFraction(const char *name);
33   AliAnalysisTaskSECharmFraction(const char *name,AliRDHFCutsD0toKpi *cutsA,AliRDHFCutsD0toKpi *cutsB);
34
35   virtual ~AliAnalysisTaskSECharmFraction(); 
36
37   // Implementation of interface methods
38   virtual void UserCreateOutputObjects();
39   virtual void Init();
40   virtual void LocalInit() {Init();}
41   virtual void UserExec(Option_t *option);
42   virtual void Terminate(Option_t *option);  
43   void SetReadMC(Bool_t readMC=kTRUE){fReadMC=readMC;}
44   void SetSplitMassD0D0bar(Bool_t splitD0D0bar=kTRUE){fsplitMassD0D0bar=splitD0D0bar;}
45   Bool_t GetIsSplitMassD0D0bar(){return fsplitMassD0D0bar;}
46   void SetUsePID(Bool_t pid){fusePID=pid;}
47   void SetAnalyzeLikeSign(Bool_t likesign=kFALSE){fLikeSign=likesign;}
48   void SetNMaxTrForVtx(const Int_t ntrMaxforVtx){fNtrMaxforVtx=ntrMaxforVtx;}
49   Int_t GetNMaxTrForVtx(){return fNtrMaxforVtx;}
50   void SetPtBins(Int_t nbins,const Float_t *ptbins);
51   void SetSignalInvMassCut(const Double_t signalInvMassCut=0.027){fsignalInvMassCut=signalInvMassCut;}
52   void SetLargeInvMassCut(const Double_t largeInvMassCut=2.){flargeInvMassCut=largeInvMassCut;}
53   void SetSideBandInvMassCut(const Double_t sidebandInvMassCut=0.054){// default value ~ 2x3 times inv mass resol.: a factor 2 is applied w.r.t. 3sigma, should be safe enough to exclude most of the reflections 
54     fsidebandInvMassCut=sidebandInvMassCut;  
55   }
56   void SetSideBandInvMassWindow(const Double_t sidebandInvMassWindow=0.108){//~ 6 times inv. mass resol.
57     fsidebandInvMassWindow=sidebandInvMassWindow;
58   }  
59   void SetAcceptanceCut(const Double_t eta=0.8,const Double_t nITSpoints=5.,const Double_t nSPDpoints=2.){fAcceptanceCuts[0]=eta;fAcceptanceCuts[1]=nITSpoints;fAcceptanceCuts[2]=nSPDpoints;}
60   void SetStandardMassSelection();
61   Int_t SetStandardCuts(Double_t pt,Double_t invMassCut);
62   Int_t SetStandardCuts(Float_t *&ptbinlimits);
63   void CheckInvMassD0(AliAODRecoDecayHF2Prong *d,Double_t &invMassD0,Double_t &invMassD0bar,Bool_t &isPeakD0,Bool_t &isPeakD0bar,Bool_t &isSideBandD0,Bool_t &isSideBandD0bar);
64   void SetAnalysisLevel(Int_t level){fFastAnalysis=level;}
65   Int_t GetAnalysisLevel(){return fFastAnalysis;}
66   Int_t CheckOrigin(const TClonesArray* arrayMC, const AliAODMCParticle *mcPartCandidate)const;
67   AliAODRecoDecayHF *GetD0toKPiSignalType(const AliAODRecoDecayHF2Prong *d,TClonesArray *arrayMC,Int_t &signaltype,Double_t &massMumTrue,Double_t *primaryVtx);
68   AliAODRecoDecayHF *GetD0toKPiSignalTypeObsolete(const AliAODRecoDecayHF2Prong *d,TClonesArray *arrayMC,Int_t &signaltype,Double_t &massMumTrue,Double_t *primaryVtx);
69   AliAODRecoDecayHF* ConstructFakeTrueSecVtx(const AliAODMCParticle *b1,const AliAODMCParticle *b2,const AliAODMCParticle *mum,Double_t *primaryVtxTrue);
70   void SetUseMC(Bool_t useMC){fUseMC=useMC;}
71   Bool_t SpecialSelD0(AliAODRecoDecayHF2Prong *d,Int_t &nusedforVtx);
72   Bool_t FillAziList(AliAODEvent *aod,Double_t azilist[30000],Int_t trkIDlist[30000],Int_t &nprim)const;
73   void FillAziHistos(AliAODRecoDecayHF2Prong *d,TList *&list,Int_t ptbin,Double_t azilist[30000],Int_t trkIDlist[30000],Int_t nprim,Int_t okD0,Int_t okD0bar,Bool_t isPeakD0,Bool_t isPeakD0bar,Bool_t isSideBandD0,Bool_t isSideBandD0bar)const;
74
75   AliAODVertex* GetPrimaryVtxSkipped(AliAODEvent *aodev,AliAODRecoDecayHF2Prong *d);
76  
77   /* ######### THE FOLLOWING IS FOR FURTHER IMPLEMENATION ############
78      Int_t GetPtBin(Double_t pt)const;
79      void SetD0Cuts(Int_t ptbin,Double_t &*d0cutsLoose,Double_t &*d0cutsTight);
80      
81      //  void InvMassSelection();
82      
83      void SetCheckMC(Bool_t checkMC){fcheckMC=checkMC;}
84      void SetCheckMC_D0(Bool_t check_D0){fcheckMCD0=check_D0;}
85      void SetCheckMC_2prongs(Bool_t check2prongs){fcheckMC2prongs=check2prongs;}
86      void SetCheckMC_prompt(Bool_t checkprompt){fcheckMCprompt=checkprompt;}
87      void SetCheckMC_fromB(Bool_t checkfromB){fcheckMCfromB=checkfromB;}
88      void SetCheckMC_fromDstar(Bool_t skipD0star){fSkipD0star=skipD0star;}
89      void SetUseCuts(Bool_t usecuts){fD0usecuts=usecuts;}
90      void SetSideBands(Double_t sideband){fSideBands=sideband;}
91      void SetStudyPureBackground(Bool_t back){fStudyPureBackground=back;}
92   */
93   AliRDHFCutsD0toKpi* GetLooseCut(){
94     return fCutsLoose;
95   }
96   AliRDHFCutsD0toKpi* GetTightCut(){
97     return fCutsTight;
98   }
99   /* void SetCutFunction(Int_t (*setcuts)(AliAnalysisTaskSECharmFraction*,Double_t,Double_t)){
100      fSetCuts=setcuts;
101      fStandCuts=kFALSE;
102      }
103   */
104   //  Int_t SetCuts(AliAnalysisTaskSECharmFraction *alchfr,Double_t pt,Double_t invMassCut);
105   
106  private:
107   Bool_t FillHistos(AliAODRecoDecayHF2Prong *d,TList *&list,Int_t ptbin,Int_t okD0,Int_t okD0bar,Double_t invMassD0,Double_t invMassD0bar,Bool_t isPeakD0,Bool_t isPeakD0bar,Bool_t isSideBandD0,Bool_t isSideBandD0bar,Double_t massmumtrue,AliAODRecoDecayHF *aodDMC,Double_t *vtxTrue);
108   void FillHistoMCproperties(TClonesArray *arrayMC);
109
110   AliRDHFCutsD0toKpi *fCutsLoose;        // Loose cuts object
111   AliRDHFCutsD0toKpi *fCutsTight;      // Vertexer heavy flavour
112   Int_t fFastAnalysis;                  // Level of analysis speed: default is 1, switch it to 2 to fill the THnSparse
113   Bool_t  fReadMC;                          // Flag To switch on/off access to MC 
114   Bool_t  fsplitMassD0D0bar;                // Flag to use two shistos for D0 and D0bar invariant masses
115   Bool_t  fLikeSign;                        // Flag to analyse Like Sign array
116   Bool_t  fusePID;                          // Flag to use PID
117   Double_t    fmD0PDG;                      //  MC D0 mass
118   Int_t        fnbins;                      // Number of pt bins
119   Float_t *fptbins;                        //[fnbins] ptbins 
120   Int_t fNtrMaxforVtx;                      // N Max acceptable tracks used for vertex (0,1,2)
121   Double_t fptAll;                          //!Sum of pt of the reco tracks
122   Double_t fptAllSq;                        //!Sum of the square of the pt of the reco tracks
123   Double_t fptMax[3];                       //!Three largest track pt in the event
124   Double_t fAcceptanceCuts[3];                // array with acceptance cuts
125   Double_t fsignalInvMassCut;               // invariant mass cut to define signal region
126   Double_t flargeInvMassCut;                // invariant mass cut to accept all inv mass window
127   Double_t fsidebandInvMassCut;             // invariant mass cut to define side band region lower limit
128   Double_t fsidebandInvMassWindow;          // invariant mass cut to define side band region width
129   Bool_t fUseMC;                            // flag to use or not MC info
130   Bool_t fCleanCandOwnVtx;                  // flag to switch on/off cleaning of the candidate own vtx
131   TH1F *fNentries;                          //!histo for #AOD analysed, container 1
132   TH1F *fSignalType;                        //!histo for the type of MC signal , container 2
133   TH1F *fSignalTypeLsCuts;                 //!histo for the type of MC signal with loose cuts , container 3
134   TH1F *fSignalTypeTghCuts;                //!histo for the type of MC signal with tight cuts, container 4
135   AliNormalizationCounter *fCounter;        //!counter for the normalization 
136   TList *flistMCproperties;               //!TLists for MC properties of D0 w.r.t. B mesons and c quarks cntainer 5
137   TList *flistNoCutsSignal;               //!TList for signal (D prompt) with nocuts, container 6
138   TList *flistNoCutsBack;               //!TList for background with nocuts, container 7
139   TList *flistNoCutsFromB;               //!TList for D from B or D from Dstar from Bwith nocuts, container 8
140   TList *flistNoCutsFromDstar;               //!TList for D from Dstar with nocuts, container 9
141   TList *flistNoCutsOther;               //!TList for others with nocuts, container 10
142   TList *flistLsCutsSignal;               //!TList for signal (D prompt) with loose cuts, container 11
143   TList *flistLsCutsBack;               //!TList for background with loose cuts, container 12
144   TList *flistLsCutsFromB;               //!TList for D from B or D from Dstar from B with loose cuts, container 13
145   TList *flistLsCutsFromDstar;               //!TList for D from Dstar with loose cuts, container 14
146   TList *flistLsCutsOther;               //!TList for others with loose cuts, container 15
147   TList *flistTghCutsSignal;               //!TList for signal (D prompt) with tight cuts, container 16
148   TList *flistTghCutsBack;               //!TList for backgrnd with tight cuts, container 17
149   TList *flistTghCutsFromB;               //!TList for D from B or D from Dstar from Bwith tight cuts, container 18
150   TList *flistTghCutsFromDstar;               //!TList for D from Dstar Dstar with tight cuts, container 19
151   TList *flistTghCutsOther;               //!TList for others with tight cuts, container 20
152   /*  Bool_t       fD0usecuts;            // Switch on the use of the cuts             TO BE IMPLEMENTED 
153       Bool_t       fcheckMC;              //  Switch on MC check: minimum check is same mother  TO BE IMPLEMENTED
154       Bool_t       fcheckMCD0;           //  check the mother is a D0                  TO BE IMPLEMENTED
155       Bool_t       fcheckMC2prongs;         //  check the decay is in two prongs       TO BE IMPLEMENTED  
156       Bool_t       fcheckMCprompt;       //  check the D0 comes from a c quark         TO BE IMPLEMENTED
157       Bool_t       fcheckMCfromB;        //  check the D0 comes from a b quark         TO BE IMPLEMENTED
158       Bool_t       fSkipD0star;           // skip if D0 comes from a D*                TO BE IMPLEMENTED
159       Bool_t  fStudyd0fromBTrue;         // Flag for analyze true impact par of D0 from B       TO BE IMPLEMENTED 
160       Bool_t  fStudyPureBackground;      // Flag to study the background (reverse the selection on the signal)     TO BE IMPLEMENTED 
161       Double_t  fSideBands;                //Side bands selection (see cxx)            TO BE IMPLEMENTED
162   */
163   AliAnalysisTaskSECharmFraction(const AliAnalysisTaskSECharmFraction&); // not implemented
164   AliAnalysisTaskSECharmFraction& operator=(const AliAnalysisTaskSECharmFraction&); // not implemented
165   
166   ClassDef(AliAnalysisTaskSECharmFraction,3); // analysis task for prompt charm fraction
167 };
168
169 #endif