]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/AliCFTaskVertexingHF.h
Possibility to 1) apply multiplicative weights, 2) consider the z-vertex correction...
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / AliCFTaskVertexingHF.h
1 #ifndef ALICFTASKVERTEXINGHF_H
2 #define ALICFTASKVERTEXINGHF_H
3 /**************************************************************************
4  * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 /* $Id$ */ 
19
20 //-----------------------------------------------------------------------
21 // Class for HF corrections as a function of many variables and step 
22 // Author : C. Zampolli, CERN
23 //                      D. Caffarri, Univ & INFN Padova caffarri@pd.infn.it
24 // Base class for HF Unfolding - agrelli@uu.nl
25 //-----------------------------------------------------------------------
26
27
28 #include "AliAnalysisTaskSE.h"
29 #include "AliCFVertexingHF2Prong.h"
30 #include "AliCFVertexingHF3Prong.h"
31 #include "AliCFVertexingHFLctoV0bachelor.h"
32 #include "AliCFVertexingHF.h"
33 #include <TH1F.h>
34 #include <TProfile.h>
35
36 class TH1I;
37 class TParticle ;
38 class TFile ;
39 class TClonesArray ;
40 class AliCFManager;
41 class AliAODRecoDecay;
42 class AliAODRecoDecayHF2Prong;
43 class AliAODMCParticle;
44 class THnSparse;
45 class TF1;
46 class AliRDHFCuts;
47 class AliCFVertexingHF2Prong;
48 class AliCFVertexingHF3Prong;
49
50 class AliCFTaskVertexingHF: public AliAnalysisTaskSE {
51 public:
52         
53         enum {
54                 kStepGeneratedLimAcc = 0,
55                 kStepGenerated       = 1,
56                 kStepAcceptance      = 2,
57                 kStepVertex          = 3,
58                 kStepRefit           = 4,
59                 kStepReconstructed   = 5,
60                 kStepRecoAcceptance  = 6,
61                 kStepRecoITSClusters = 7,
62                 kStepRecoPPR         = 8,
63                 kStepRecoPID         = 9
64         };
65
66         enum {
67                 kSnail = 0,    // slow configuration, all variables
68                 kCheetah = 1   // fast configuration, only a subset of variables
69         };
70
71         enum {
72           kAll = 0,   // all decays (resonant + non-resonant)
73           kNonResonant = 1, // only non resonant
74           kL1520 = 2,  // Lc --> L(1520) + p
75           kKstar = 3,  // Lc --> K* + pi
76           kDelta = 4   // Lc --> Delta + K
77         };
78
79         enum { kNtrk10=0, kNtrk10to16=1, kVZERO=2 }; // multiplicity estimators
80         
81         AliCFTaskVertexingHF();
82         AliCFTaskVertexingHF(const Char_t* name, AliRDHFCuts* cuts, TF1* func = 0x0);
83         AliCFTaskVertexingHF& operator= (const AliCFTaskVertexingHF& c);
84         AliCFTaskVertexingHF(const AliCFTaskVertexingHF& c);
85         virtual ~AliCFTaskVertexingHF();
86         
87         // ANALYSIS FRAMEWORK STUFF to loop on data and fill output objects
88         void     UserCreateOutputObjects();
89         void     UserExec(Option_t *option);
90         void     Init();
91         void     LocalInit() {Init();}
92         void     Terminate(Option_t *);
93         
94         // UNFOLDING
95         void     SetCorrelationMatrix(THnSparse* h) {fCorrelation=h;}
96         void     SetAcceptanceUnf(Bool_t AcceptanceUnf) {fAcceptanceUnf = AcceptanceUnf;}
97         Bool_t   GetAcceptanceUnf() const {return fAcceptanceUnf;}
98         
99         
100         // CORRECTION FRAMEWORK RELATED FUNCTIONS
101         void           SetCFManager(AliCFManager* io) {fCFManager = io;}   // global correction manager
102         AliCFManager * GetCFManager()                 {return fCFManager;} // get corr manager
103         
104         // Setters (and getters) for the config macro
105         void    SetFillFromGenerated(Bool_t flag) {fFillFromGenerated = flag;}
106         Bool_t  GetFillFromGenerated() const {return fFillFromGenerated;}
107         void    SetDecayChannel (Int_t decayChannel) {fDecayChannel = decayChannel;}
108         Int_t   GetDecayChannel () {return fDecayChannel;}
109         void     SetUseWeight(Bool_t useWeight){fUseWeight=useWeight;}
110         Bool_t   GetUseWeight() const {return fUseWeight;}
111         Double_t GetWeight(Float_t pt);
112         Double_t dNdptFit(Float_t pt, Double_t* par);
113
114         void SetUseFlatPtWeight(Bool_t useWeight){fUseFlatPtWeight=useWeight; fUseWeight=useWeight;}
115         Bool_t GetUseFlatPtWeight() const {return fUseFlatPtWeight;}
116         void SetUseZWeight(Bool_t useWeight){fUseZWeight=useWeight;}
117         Bool_t GetUseZWeight() const {return fUseZWeight;}
118         Double_t GetZWeight(Float_t z, Int_t runnumber);
119         Double_t DodzFit(Float_t z, Double_t* par);
120
121         void SetUseNchWeight(Bool_t useWeight){fUseNchWeight=useWeight;}
122         Bool_t GetUseNchWeight() const {return fUseNchWeight;}
123         void SetMCNchHisto(TH1F* h){
124           if(fHistoMCNch) delete fHistoMCNch;
125           fHistoMCNch=new TH1F(*h);
126         }
127         void CreateMeasuredNchHisto();
128         Double_t GetNchWeight(Int_t nch);
129         void SetMultiplicityEstimator(Int_t value){ fMultiplicityEstimator=value; }
130         Int_t GetMultiplicityEstimator(){ return fMultiplicityEstimator; }
131         void SetIsPPData(Bool_t flag){ fIsPPData = flag; }
132
133         void SetUseNchTrackletsWeight(Bool_t useWeight = kTRUE) { fUseNchWeight=useWeight; fUseTrackletsWeight=useWeight; }
134         Bool_t GetUseNchTrackletsWeight() const {return fUseTrackletsWeight;}
135
136         void SetUseZvtxCorrectedNtrkEstimator(Bool_t flag) { fZvtxCorrectedNtrkEstimator=flag; }
137         Bool_t GetUseZvtxCorrectedNtrkEstimator() { return fZvtxCorrectedNtrkEstimator; }
138         void SetMultiplVsZProfileLHC10b(TProfile* hprof){
139           if(fMultEstimatorAvg[0]) delete fMultEstimatorAvg[0];
140           fMultEstimatorAvg[0]=new TProfile(*hprof);
141         }
142         void SetMultiplVsZProfileLHC10c(TProfile* hprof){
143           if(fMultEstimatorAvg[1]) delete fMultEstimatorAvg[1];
144           fMultEstimatorAvg[1]=new TProfile(*hprof);
145         }
146         void SetMultiplVsZProfileLHC10d(TProfile* hprof){
147           if(fMultEstimatorAvg[2]) delete fMultEstimatorAvg[2];
148           fMultEstimatorAvg[2]=new TProfile(*hprof);
149         }
150         void SetMultiplVsZProfileLHC10e(TProfile* hprof){
151           if(fMultEstimatorAvg[3]) delete fMultEstimatorAvg[3];
152           fMultEstimatorAvg[3]=new TProfile(*hprof);
153         }
154         TProfile* GetEstimatorHistogram(const AliVEvent* event);
155         void SetReferenceMultiplcity(Double_t rmu){fRefMult=rmu;}
156         
157         void   SetDselection(UShort_t originDselection) {fOriginDselection=originDselection;}
158         UShort_t GetDselection (){return fOriginDselection;}
159         void SetSign(Char_t isSign) {fSign = isSign;}
160         Char_t GetSign() {return fSign;}
161          
162         void SetCentralitySelection(Bool_t centSelec = kTRUE) {fCentralitySelection = centSelec;}   
163         Bool_t GetCentralitySelection() {return fCentralitySelection;} 
164
165         void SetFakeSelection(Int_t fakeSel = 0) {fFakeSelection=fakeSel;}
166         Int_t GetFakeSelection(){return fFakeSelection;}
167
168         void SetRejectCandidateIfNotFromQuark(Bool_t opt){fRejectIfNoQuark=opt;}
169         Bool_t GetRejectCandidateIfNotFromQuark(){return fRejectIfNoQuark;}
170
171         void SetUseMCVertex(Bool_t opt){fUseMCVertex=opt;}
172         Bool_t GetUseMCVertex(){return fUseMCVertex;}
173         
174
175         void SetKeepDsViaPhi(){fDsOption=1;}
176         void SetKeepDsViaK0star(){fDsOption=2;}
177         void SetKeepAllDs(){fDsOption=3;}
178         void SetCountAllDs(){fGenDsOption=AliCFVertexingHF3Prong::kCountAllDsKKpi;}
179         void SetCountDsViaPhi(){fGenDsOption=AliCFVertexingHF3Prong::kCountPhipi;}
180         void SetCountDsViaK0star(){fGenDsOption=AliCFVertexingHF3Prong::kCountK0stK;}
181         void SetCountResonantDs(){fGenDsOption=AliCFVertexingHF3Prong::kCountResonant;}
182         void SetCountNonResonantDs(){fGenDsOption=AliCFVertexingHF3Prong::kCountNonResonant;}
183
184         Bool_t ProcessDs(Int_t returnCodeDs) const;
185
186         void SetConfiguration(Int_t configuration) {(configuration == kSnail) ? Printf("Slow configuration chosen, all variables will be used!") : Printf("Fast configuration chosen, all variablesOnly pt, y, phi, ct, fake, z_vtx, centrality and multiplicity will be used!"); fConfiguration = configuration;} 
187         Int_t GetConfiguration() const {return fConfiguration;} 
188         
189         void SetWeightFunction(TF1* func) {fFuncWeight = func;}
190         TF1* GetWeightFunction() const {return fFuncWeight;}
191         void SetPtWeightsFromFONLL276overLHC12a17a();
192         void SetPtWeightsFromDataPbPb276overLHC12a17a();
193         void SetPtWeightsFromFONLL276overLHC12a17b();
194         void SetPtWeightsFromFONLL276andBAMPSoverLHC12a17b();
195         void SetPtWeightsFromFONLL5overLHC13d3();
196
197         void SetResonantDecay(UInt_t resonantDecay) {fResonantDecay = resonantDecay;}
198         UInt_t GetResonantDecay() const {return fResonantDecay;}
199
200         void SetKeepLctoK0Sp() {fLctoV0bachelorOption=1;}
201         void SetKeepLctoLambdaBarpi() {fLctoV0bachelorOption=2;}
202         void SetKeepLctoLambdapi() {fLctoV0bachelorOption=4;}
203         void SetKeepLctoV0bachelor() {fLctoV0bachelorOption=7;}
204
205         void SetCountLctoK0Sp(){fGenLctoV0bachelorOption=AliCFVertexingHFLctoV0bachelor::kCountK0Sp;}
206         void SetCountLctoLambdapi(){fGenLctoV0bachelorOption=AliCFVertexingHFLctoV0bachelor::kCountLambdapi;}
207     
208         void SetUseSelectionBit(Bool_t flag) { fUseSelectionBit=flag; }
209         Bool_t GetUseSelectionBit() const { return fUseSelectionBit; }
210
211         Bool_t ProcessLctoV0Bachelor(Int_t returnCodeDs) const;
212
213 protected:
214         AliCFManager   *fCFManager;   //  pointer to the CF manager
215         TH1I *fHistEventsProcessed;   //! simple histo for monitoring the number of events processed
216         THnSparse* fCorrelation;      //  response matrix for unfolding
217         TList  *fListProfiles; //list of profile histos for z-vtx correction
218         Int_t fCountMC;               //  MC particle found
219         Int_t fCountAcc;              //  MC particle found that satisfy acceptance cuts
220         Int_t fCountVertex;       //  Reco particle found that satisfy vertex constrained
221         Int_t fCountRefit;        //  Reco particle found that satisfy kTPCrefit and kITSrefit
222         Int_t fCountReco;             //  Reco particle found that satisfy cuts
223         Int_t fCountRecoAcc;          //  Reco particle found that satisfy cuts in requested acceptance
224         Int_t fCountRecoITSClusters;  //  Reco particle found that satisfy cuts in n. of ITS clusters
225         Int_t fCountRecoPPR;          //  Reco particle found that satisfy cuts in PPR
226         Int_t fCountRecoPID;          //Reco PID step 
227         Int_t fEvents;                //  n. of events
228         Int_t fDecayChannel;          // decay channel to configure the task
229         Bool_t fFillFromGenerated;    //  flag to indicate whether data container should be filled with generated values also for reconstructed particles
230         UShort_t fOriginDselection;      // flag to select D0 origins. 0 Only from charm 1 only from beauty 2 both from charm and beauty
231         Bool_t fAcceptanceUnf;        //  flag for unfolding before or after cuts.
232         AliRDHFCuts* fCuts;            // cuts
233         Bool_t fUseWeight;             //flag to decide whether to use pt-weights != 1 when filling the container or not
234         Double_t fWeight;              //weight used to fill the container
235         Bool_t fUseFlatPtWeight;       // flag to decide to use a flat pt shape
236         Bool_t fUseZWeight;           // flag to decide whether to use z-vtx weights != 1 when filling the container or not
237         Bool_t fUseNchWeight;         // flag to decide whether to use Ncharged weights != 1 when filling the container or not
238         Bool_t fUseTrackletsWeight;   // flag to decide whether to use Ncharged weights != 1 when filling the container or not
239         Int_t fNvar;                  // number of variables for the container
240         TString fPartName;    // D meson name
241         TString fDauNames;    // daughter in fin state
242         Char_t fSign;                 // flag to decide wheter to keep D0 only (0), D0bar only (1), or both D0 and D0bar (2)
243         Bool_t fCentralitySelection;  //flag to switch off the centrality selection
244         Int_t  fFakeSelection;  //selection flag for fakes tracks 
245         Bool_t fRejectIfNoQuark;  // flag to remove events not geenrated with PYTHIA
246         Bool_t fUseMCVertex;  // flag to use MC vertex (useful when runnign in pp)
247         Int_t  fDsOption;     // Ds decay option (selection level)
248         Int_t  fGenDsOption;     // Ds decay option (generation level)
249         Int_t fConfiguration; // configuration (slow / fast) of the CF --> different variables will be allocated (all / reduced number)
250         TF1* fFuncWeight;     // user-defined function to be used to calculate weights
251         TH1F* fHistoMeasNch;  // histogram with measured Nch distribution (pp 7 TeV)
252         TH1F* fHistoMCNch;  // histogram with Nch distribution from MC production
253         UInt_t fResonantDecay;  // resonant deacy channel to be used if the CF should be run on resonant channels only
254         Int_t fLctoV0bachelorOption; // Lc->V0+bachelor decay option (selection level)
255         Int_t fGenLctoV0bachelorOption; // Lc->V0+bachelor decay option (generation level)
256         Bool_t fUseSelectionBit;     // flag to use selection bit
257         UInt_t fPDGcode; // PDG code
258
259         Int_t fMultiplicityEstimator; // Definition of the multiplicity estimator: kNtrk10=0, kNtrk10to16=1, kVZERO=2
260         TProfile* fMultEstimatorAvg[4]; // TProfile with mult vas. Z per period
261         Double_t fRefMult;   // refrence multiplcity (period b)
262         Bool_t fZvtxCorrectedNtrkEstimator; // flag to use the z-vtx corrected (if not use uncorrected) multiplicity estimator
263         Bool_t fIsPPData; // flag for pp data (not checking centrality)
264
265         ClassDef(AliCFTaskVertexingHF,18); // class for HF corrections as a function of many variables
266 };
267
268 #endif