]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/UserTasks/AliAnalysisTaskPID.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGJE / UserTasks / AliAnalysisTaskPID.h
1 #ifndef ALI_ANALYSIS_TASK_PID_H
2 #define ALI_ANALYSIS_TASK_PID_H
3
4 /*
5 This task collects PID output from different detectors.
6 Only tracks fulfilling some standard quality cuts are taken into account.
7 At the moment, only data from TPC and TOF is collected. But in future,
8 data from e.g. HMPID is also foreseen.
9
10 Class written by Benjamin Hess.
11 Contact: bhess@cern.ch
12 */
13
14 class TF1;
15 class TRandom3;
16 class AliAnalysisFilter;
17 class AliCFContainer;
18 class AliESDEvent;
19 class AliMCEvent;
20 class AliMCParticle;
21 class AliPID;
22 class AliPIDCombined;
23 class AliPIDResponse;
24 class AliTOFPIDResponse;
25 class AliVEvent;
26 class AliVTrack;
27
28 #include "TAxis.h"
29 #include "TH1D.h"
30 #include "TH2D.h"
31 #include "TH3D.h"
32 #include "THnSparse.h"
33 #include "TProfile.h"
34 #include "TString.h"
35
36 #include "AliCentrality.h"
37 #include "AliCFContainer.h"
38
39 #include "AliPID.h"
40 #include "AliAnalysisTaskPIDV0base.h"
41
42 class AliAnalysisTaskPID : public AliAnalysisTaskPIDV0base {
43  public:
44   AliAnalysisTaskPID();
45   AliAnalysisTaskPID(const char *name);
46   virtual ~AliAnalysisTaskPID();
47   
48   virtual void   UserCreateOutputObjects();
49   virtual void   UserExec(Option_t *option);
50   virtual void   Terminate(const Option_t*);
51   
52   enum ErrorCode { kNoErrors = 1, kWarning = 0, kError = -1};
53   
54   enum dataAxes { kDataMCID = 0, kDataSelectSpecies = 1, kDataPt = 2, kDataDeltaPrimeSpecies = 3, kDataCentrality = 4,
55                   kDataJetPt = 5, kDataZ = 6, kDataXi = 7, kDataCharge = 8, kDataTOFpidInfo = 9, kDataNumAxes = 10 };
56
57   enum genAxes { kGenMCID = 0, kGenSelectSpecies = 1, kGenPt = 2, kGenDeltaPrimeSpecies = 3, kGenCentrality = 4,
58                  kGenJetPt = 5, kGenZ = 6, kGenXi = 7, kGenCharge = 8, kGenTOFpidInfo = 9, kGenNumAxes = 10 };
59   
60   enum genYieldAxes { kGenYieldMCID = 0, kGenYieldPt = 1, kGenYieldCentrality = 2, kGenYieldJetPt = 3, kGenYieldZ = 4, kGenYieldXi = 5,
61                       kGenYieldCharge = 6, kGenYieldNumAxes = 7 };
62   
63   enum ptResolutionAxes { kPtResJetPt = 0, kPtResGenPt = 1, kPtResRecPt = 2, kPtResCharge = 3, kPtResCentrality = 4, kPtResNumAxes = 5 };
64   
65   enum dEdxCheckAxes { kDeDxCheckPID = 0, kDeDxCheckP = 1, kDeDxCheckJetPt = 2, kDeDxCheckEtaAbs = 3 , kDeDxCheckDeDx = 4,
66                        kDeDxCheckNumAxes = 5 };
67   
68   enum efficiencyAxes { kEffMCID = 0, kEffTrackPt = 1, kEffTrackEta = 2, kEffTrackCharge = 3, kEffCentrality = 4, kEffJetPt = 5,
69                         kEffZ = 6, kEffXi = 7, kEffNumAxes = 8 };
70   
71   enum EffSteps { kStepGenWithGenCuts = 0, kStepRecWithGenCuts = 1, kStepRecWithGenCutsMeasuredObs = 2,
72                   kStepRecWithRecCutsMeasuredObs = 3, kStepRecWithRecCutsMeasuredObsPrimaries = 4,
73                   kStepRecWithRecCutsMeasuredObsStrangenessScaled = 5, kStepRecWithRecCutsPrimaries = 6, kNumSteps = 7};
74   
75   enum TOFpidInfo { kNoTOFinfo = -2, kNoTOFpid = -1, kTOFpion = 0, kTOFkaon = 1, kTOFproton = 2, kNumTOFspecies = 3,
76                     kNumTOFpidInfoBins = 5 };
77   
78   enum EventCounterType { kTriggerSel = 0, kTriggerSelAndVtxCut = 1, kTriggerSelAndVtxCutAndZvtxCut = 2 };
79   
80   static Int_t PDGtoMCID(Int_t pdg);
81   
82   static void GetJetTrackObservables(Double_t trackPt, Double_t jetPt, Double_t& z, Double_t& xi);
83   
84   static Double_t GetMCStrangenessFactorCMS(Int_t motherPDG, Double_t motherGenPt);
85   static Double_t GetMCStrangenessFactorCMS(AliMCEvent* mcEvent, AliMCParticle* daughter);
86   
87   static Bool_t IsSecondaryWithStrangeMotherMC(AliMCEvent* mcEvent, Int_t partLabel);
88   
89   Int_t GetIndexOfChargeAxisData() const
90     { return fStoreAdditionalJetInformation ? kDataCharge : kDataCharge - fgkNumJetAxes; };
91   Int_t GetIndexOfChargeAxisGen() const
92     { return fStoreAdditionalJetInformation ? kGenCharge : kGenCharge - fgkNumJetAxes; };
93   Int_t GetIndexOfChargeAxisGenYield() const
94     { return fStoreAdditionalJetInformation ? kGenYieldCharge : kGenYieldCharge - fgkNumJetAxes; };
95   
96   Int_t GetIndexOfTOFpidInfoAxisData() const
97     { return fStoreAdditionalJetInformation ? kDataTOFpidInfo : kDataTOFpidInfo - fgkNumJetAxes; };
98   Int_t GetIndexOfTOFpidInfoAxisGen() const
99     { return fStoreAdditionalJetInformation ? kGenTOFpidInfo : kGenTOFpidInfo - fgkNumJetAxes; };
100   
101   Bool_t FillXsec(Double_t xsection)
102     { if (!fh1Xsec) return kFALSE; fh1Xsec->Fill("<#sigma>", xsection); return kTRUE; };
103   Bool_t FillPythiaTrials(Double_t avgTrials)
104     { if (!fh1Trials) return kFALSE; fh1Trials->Fill("#sum{ntrials}", avgTrials); return kTRUE; };
105     
106   Bool_t FillEfficiencyContainer(const Double_t* values, EffSteps step, Double_t weight = 1.0);
107   
108   Bool_t FillGeneratedYield(const Double_t* values, Double_t weight = 1.0);
109   Bool_t FillPtResolution(Int_t mcID, const Double_t* values);
110   Bool_t FillGenJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm = -1.);
111   Bool_t FillRecJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm = -1.);
112   
113   Bool_t IncrementEventCounter(Double_t centralityPercentile, EventCounterType type);
114   
115   void PostOutputData();
116   
117   void PrintSettings(Bool_t printSystematicsSettings = kFALSE) const;
118   
119   void PrintSystematicsSettings() const;
120   
121   Bool_t ProcessTrack(const AliVTrack* track, Int_t particlePDGcode, Double_t centralityPercentile, Double_t jetPt) ;
122   
123   ErrorCode GenerateDetectorResponse(ErrorCode errCode, Double_t mean, Double_t sigma, Double_t* responses,
124                                      Int_t nResponses,
125                                      Bool_t usePureGaus = kFALSE);
126   ErrorCode SetParamsForConvolutedGaus(Double_t gausMean, Double_t gausSigma);
127   
128   const TString GetCentralityEstimator() const { return fCentralityEstimator; };
129   void SetCentralityEstimator(TString estimator) { fCentralityEstimator = estimator; };
130   
131   Double_t GetCentralityPercentile(AliVEvent* evt) const;
132   
133   inline Double_t GetConvolutedGaussTransitionPar(Int_t index) const;
134   
135   Bool_t SetConvolutedGaussLambdaParameter(Double_t lambda);
136
137   Bool_t GetInputFromOtherTask() const { return fInputFromOtherTask; };
138   void SetInputFromOtherTask(Bool_t flag) { fInputFromOtherTask = flag; };
139   
140   Bool_t GetDoPID() const { return fDoPID; };
141   void SetDoPID(Bool_t flag) { fDoPID = flag; };
142   
143   Bool_t GetDoEfficiency() const { return fDoEfficiency; };
144   void SetDoEfficiency(Bool_t flag) { fDoEfficiency = flag; };
145   
146   Bool_t GetDoPtResolution() const { return fDoPtResolution; };
147   void SetDoPtResolution(Bool_t flag) { fDoPtResolution = flag; };
148   
149   Bool_t GetDoDeDxCheck() const { return fDoDeDxCheck; };
150   void SetDoDeDxCheck(Bool_t flag) { fDoDeDxCheck = flag; };
151   
152   Bool_t GetStoreCentralityPercentile() const { return fStoreCentralityPercentile; };
153   void SetStoreCentralityPercentile(Bool_t flag) { fStoreCentralityPercentile = flag; };
154   
155   Bool_t GetStoreAdditionalJetInformation() const { return fStoreAdditionalJetInformation; };
156   void SetStoreAdditionalJetInformation(Bool_t flag) { fStoreAdditionalJetInformation = flag; };
157   
158   Bool_t GetUseMCidForGeneration() const { return fUseMCidForGeneration; };
159   void SetUseMCidForGeneration(Bool_t flag) { fUseMCidForGeneration = flag; };
160   
161   Bool_t GetUseConvolutedGaus() const { return fUseConvolutedGaus; };
162   void SetUseConvolutedGaus(Bool_t flag) { fUseConvolutedGaus = flag; };
163   
164   Double_t GetAccuracyNonGaussianTail() const { return fAccuracyNonGaussianTail; };
165   void SetAccuracyNonGaussianTail(Double_t value) { fAccuracyNonGaussianTail = value; };
166   
167   Bool_t GetTakeIntoAccountMuons() const { return fTakeIntoAccountMuons; };
168   void SetTakeIntoAccountMuons(Bool_t flag) { fTakeIntoAccountMuons = flag; };
169   
170   Int_t GetTOFmode() const { return fTOFmode; };
171   void SetTOFmode(Int_t tofMode) { fTOFmode = tofMode; };
172   
173   Bool_t GetUseTPCDefaultPriors() const { return fTPCDefaultPriors; };
174   void SetUseTPCDefaultPriors(Bool_t flag) { fTPCDefaultPriors = flag; };
175   
176   Bool_t GetUsePriors() const { return fUsePriors; };
177   void SetUsePriors(Bool_t flag) { fUsePriors = flag; };
178   
179   Bool_t GetUseITS() const { return fUseITS; };
180   void SetUseITS(Bool_t flag) { fUseITS = flag; };
181   
182   Bool_t GetUseTOF() const { return fUseTOF; };
183   void SetUseTOF(Bool_t flag) { fUseTOF = flag; };
184   
185   Double_t GetEtaAbsCutLow() const { return fEtaAbsCutLow; };
186   Double_t GetEtaAbsCutUp() const { return fEtaAbsCutUp; };
187   Bool_t SetEtaAbsCutRange(Double_t lowerLimit, Double_t upperLimit);
188   
189   Bool_t IsInAcceptedEtaRange(Double_t etaAbs) const { return (etaAbs >= fEtaAbsCutLow && etaAbs <= fEtaAbsCutUp); };
190   
191   Double_t GetSystematicScalingSplinesThreshold() const { return fSystematicScalingSplinesThreshold; };
192   void SetSystematicScalingSplinesThreshold(Double_t threshold) { fSystematicScalingSplinesThreshold = threshold; };
193   
194   Double_t GetSystematicScalingSplinesBelowThreshold() const { return fSystematicScalingSplinesBelowThreshold; };
195   void SetSystematicScalingSplinesBelowThreshold(Double_t scaleFactor) 
196     { fSystematicScalingSplinesBelowThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
197     
198   Double_t GetSystematicScalingSplinesAboveThreshold() const { return fSystematicScalingSplinesAboveThreshold; };
199   void SetSystematicScalingSplinesAboveThreshold(Double_t scaleFactor) 
200     { fSystematicScalingSplinesAboveThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
201   
202   Double_t GetSystematicScalingEtaCorrectionMomentumThr() const { return fSystematicScalingEtaCorrectionMomentumThr; };
203   void SetSystematicScalingEtaCorrectionMomentumThr(Double_t threshold) { fSystematicScalingEtaCorrectionMomentumThr = threshold; };
204   
205   Double_t GetSystematicScalingEtaCorrectionLowMomenta() const { return fSystematicScalingEtaCorrectionLowMomenta; };
206   void SetSystematicScalingEtaCorrectionLowMomenta(Double_t scaleFactor) 
207     { fSystematicScalingEtaCorrectionLowMomenta = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
208   
209   Double_t GetSystematicScalingEtaCorrectionHighMomenta() const { return fSystematicScalingEtaCorrectionHighMomenta; };
210   void SetSystematicScalingEtaCorrectionHighMomenta(Double_t scaleFactor) 
211     { fSystematicScalingEtaCorrectionHighMomenta = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
212   
213   Double_t GetSystematicScalingEtaSigmaParaThreshold() const { return fSystematicScalingEtaSigmaParaThreshold; };
214   void SetSystematicScalingEtaSigmaParaThreshold(Double_t threshold) { fSystematicScalingEtaSigmaParaThreshold = threshold; };
215   
216   Double_t GetSystematicScalingEtaSigmaParaBelowThreshold() const { return fSystematicScalingEtaSigmaParaBelowThreshold; };
217   void SetSystematicScalingEtaSigmaParaBelowThreshold(Double_t scaleFactor)
218     { fSystematicScalingEtaSigmaParaBelowThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
219   
220   Double_t GetSystematicScalingEtaSigmaParaAboveThreshold() const { return fSystematicScalingEtaSigmaParaAboveThreshold; };
221   void SetSystematicScalingEtaSigmaParaAboveThreshold(Double_t scaleFactor)
222     { fSystematicScalingEtaSigmaParaAboveThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
223   
224   Double_t GetSystematicScalingMultCorrection() const { return fSystematicScalingMultCorrection; };
225   void SetSystematicScalingMultCorrection(Double_t scaleFactor) 
226     { fSystematicScalingMultCorrection = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
227   
228   Double_t GetMaxEtaVariation(Double_t dEdxSplines);
229   Bool_t CalculateMaxEtaVariationMapFromPIDResponse();
230   
231   void CleanupParticleFractionHistos();
232   Bool_t GetParticleFraction(Double_t trackPt, Double_t jetPt, Double_t multiplicity,
233                              AliPID::EParticleType species, Double_t& fraction, Double_t& fractionErrorStat,
234                              Double_t& fractionErrorSys) const;
235   Bool_t GetParticleFractions(Double_t trackPt, Double_t jetPt, Double_t centralityPercentile,
236                               Double_t* prob, Int_t smearSpeciesByError, Int_t takeIntoAccountSpeciesSysError,
237                               Bool_t uniformSystematicError = kFALSE) const;
238   const TH3D* GetParticleFractionHisto(Int_t species, Bool_t sysError = kFALSE) const;
239   Bool_t SetParticleFractionHisto(const TH3D* hist, Int_t species, Bool_t sysError = kFALSE);
240   Int_t GetParticleFractionHistoNbinsTrackPt() const;
241   Int_t GetParticleFractionHistoNbinsJetPt() const;
242   Int_t GetParticleFractionHistoNbinsCentrality() const;
243   Bool_t SetParticleFractionHistosFromFile(const TString filePathName, Bool_t sysError = kFALSE);
244   Int_t GetRandomParticleTypeAccordingToParticleFractions(Double_t trackPt, Double_t jetPt, 
245                                                           Double_t centralityPercentile,
246                                                           Bool_t smearByError,
247                                                           Bool_t takeIntoAccountSysError = kFALSE) const;
248   
249   TOFpidInfo GetTOFType(const AliVTrack* track, Int_t tofMode) const;
250   
251  protected:
252   void CheckDoAnyStematicStudiesOnTheExpectedSignal();
253   Double_t ConvolutedGaus(const Double_t* xx, const Double_t* par) const;
254   inline Double_t FastGaus(Double_t x, Double_t mean, Double_t sigma) const;
255   inline Double_t FastNormalisedGaus(Double_t x, Double_t mean, Double_t sigma) const;
256   Int_t FindBinWithinRange(TAxis* axis, Double_t value) const;
257   Int_t FindFirstBinAboveIn3dSubset(const TH3* hist, Double_t threshold, Int_t yValue, Int_t zValue) const;
258   Int_t FindLastBinAboveIn3dSubset(const TH3* hist, Double_t threshold, Int_t yValue, Int_t zValue) const;
259   virtual void SetUpGenHist(THnSparse* hist, Double_t* binsPt, Double_t* binsDeltaPrime, Double_t* binsCent, Double_t* binsJetPt) const;
260   virtual void SetUpGenYieldHist(THnSparse* hist, Double_t* binsPt, Double_t* binsCent, Double_t* binsJetPt) const;
261   virtual void SetUpHist(THnSparse* hist, Double_t* binsPt, Double_t* binsDeltaPrime, Double_t* binsCent, Double_t* binsJetPt) const;
262   virtual void SetUpPtResHist(THnSparse* hist, Double_t* binsPt, Double_t* binsJetPt, Double_t* binsCent) const;
263   virtual void SetUpDeDxCheckHist(THnSparse* hist, const Double_t* binsPt, const Double_t* binsJetPt, const Double_t* binsEtaAbs) const;
264   virtual void SetUpPIDcombined();
265   
266   static const Int_t fgkNumJetAxes; // Number of additional axes for jets
267   static const Double_t fgkEpsilon; // Double_t threshold above zero
268   static const Int_t fgkMaxNumGenEntries; // Maximum number of generated detector responses per track and delta(Prime) and associated species
269
270  private:
271   static const Double_t fgkOneOverSqrt2; // = 1. / TMath::Sqrt2();
272   
273   Int_t fRun; // Current run number
274   AliPIDCombined* fPIDcombined; //! PID combined object
275   
276   Bool_t fInputFromOtherTask; // If set to kTRUE, no events are processed and the input must be fed in from another task. If set to kFALSE, normal event processing
277   
278   Bool_t fDoPID; // Only do PID processing (and post the output), if flag is set to kTRUE
279   Bool_t fDoEfficiency; // Only do efficiency processing (and post the output), if flag is set to kTRUE
280   Bool_t fDoPtResolution; // Only do pT resolution processing (and post the output), if flag is set to kTRUE
281   Bool_t fDoDeDxCheck; // Only check dEdx, if flag set to kTRUE
282   
283   Bool_t fStoreCentralityPercentile; // If set to kTRUE, store centrality percentile for each event. In case of kFALSE (appropriate for pp), centrality percentile will be set to -1 for every event
284   Bool_t fStoreAdditionalJetInformation; // If set to kTRUE, additional jet information like jetPt, z, xi will be stored in the THnSparses
285
286   Bool_t fTakeIntoAccountMuons; // Also take into account muons for the generation of the expected response and the most probable PID
287   Bool_t fUseITS; // Use ITS for PID combined probabilities
288   Bool_t fUseTOF; // Use TOF for PID combined probabilities
289   Bool_t fUsePriors; // Use priors for PID combined probabilities
290   Bool_t fTPCDefaultPriors; // Use TPC default priors for PID combined probabilities, if priors are enabled
291     
292   Bool_t fUseMCidForGeneration; // If MC, use MCid instead of PIDcombined to generate the signals
293   
294   Bool_t fUseConvolutedGaus; // Use convoluted gaus to generate detector response instead of pure gaus  
295   const Int_t fkConvolutedGausNPar; // Number of parameters for convoluted gaus
296   Double_t fAccuracyNonGaussianTail; // Accuracy of the non-gaussian tail (fraction of the maximum)
297   const Double_t fkDeltaPrimeLowLimit; // Lower deltaPrime limit
298   const Double_t fkDeltaPrimeUpLimit; // Upper deltaPrime limit
299   TF1* fConvolutedGausDeltaPrime; // Gaus convoluted with exponential tail to generate detector response (deltaPrime)
300   
301   Int_t fTOFmode; // TOF mode used for TOF PID info (affects num sigma inclusion/exclusion)
302   Double_t fConvolutedGaussTransitionPars[3]; // Parameter for transition from gaussian parameters to asymmetric shape
303   static const Double_t fgkSigmaReferenceForTransitionPars; // Reference sigma chosen to calculate transition parameters
304   
305   Double_t fEtaAbsCutLow; // Lower cut value on |eta|
306   Double_t fEtaAbsCutUp;  // Upper cut value on |eta|
307   
308   // For systematic studies
309   Bool_t   fDoAnySystematicStudiesOnTheExpectedSignal; // Internal flag indicating whether any systematic studies are going to be performed
310   Double_t fSystematicScalingSplinesThreshold;         // beta-gamma threshold for the systematic spline scale factor
311   Double_t fSystematicScalingSplinesBelowThreshold;        // Systematic scale factor for the splines (1. = no systematics) below threshold
312   Double_t fSystematicScalingSplinesAboveThreshold;        // Systematic scale factor for the splines (1. = no systematics) above threshold
313   Double_t fSystematicScalingEtaCorrectionMomentumThr;  // Momentum threshold for the systematic scale factor for the eta correction (separates low-p from high-p
314   Double_t fSystematicScalingEtaCorrectionLowMomenta;   // Systematic scale factor for the eta correction (1. = no systematics) at low momenta
315   Double_t fSystematicScalingEtaCorrectionHighMomenta;  // Systematic scale factor for the eta correction (1. = no systematics) at high momenta
316   
317   Double_t fSystematicScalingEtaSigmaParaThreshold; // dEdx threshold for the systematic scale factor for the parametrisation of the relative signal width
318   Double_t fSystematicScalingEtaSigmaParaBelowThreshold; // Systematic scale factor for the parametrisation of the relative signal width (1. = no systematics) below threshold
319   Double_t fSystematicScalingEtaSigmaParaAboveThreshold; // Systematic scale factor for the parametrisation of the relative signal width (1. = no systematics) above threshold 
320   Double_t fSystematicScalingMultCorrection; // Systematic scale factor for the multiplicity correction (1. = no systematics) 
321   
322   TH3D* fFractionHists[AliPID::kSPECIES]; // 3D histos of particle fraction as function  with trackPt, jetPt (-1 for inclusive spectra), centralityPercentile (-1 for pp)
323   TH3D* fFractionSysErrorHists[AliPID::kSPECIES]; // 3D histos of sys. error of particle fraction as function  with trackPt, jetPt (-1 for inclusive spectra), centralityPercentile (-1 for pp)
324   
325   TString fCentralityEstimator; // Estimator for the centrality (e.g. V0A, V0M)
326   
327   THnSparseD* fhPIDdataAll; //! Data histo
328   
329   // Generated response information
330   THnSparseD* fhGenEl; //! Generated response for el
331   THnSparseD* fhGenKa; //! Generated response for ka
332   THnSparseD* fhGenPi; //! Generated response for pi
333   THnSparseD* fhGenMu; //! Generated response for mu
334   THnSparseD* fhGenPr; //! Generated response for pr
335   
336   // Generated responses for a single track
337   Double_t* fGenRespElDeltaPrimeEl; //! Generated responses for a single track
338   Double_t* fGenRespElDeltaPrimeKa; //! Generated responses for a single track
339   Double_t* fGenRespElDeltaPrimePi; //! Generated responses for a single track
340   Double_t* fGenRespElDeltaPrimePr; //! Generated responses for a single track
341   Double_t* fGenRespKaDeltaPrimeEl; //! Generated responses for a single track
342   Double_t* fGenRespKaDeltaPrimeKa; //! Generated responses for a single track
343   Double_t* fGenRespKaDeltaPrimePi; //! Generated responses for a single track
344   Double_t* fGenRespKaDeltaPrimePr; //! Generated responses for a single track
345   Double_t* fGenRespPiDeltaPrimeEl; //! Generated responses for a single track
346   Double_t* fGenRespPiDeltaPrimeKa; //! Generated responses for a single track
347   Double_t* fGenRespPiDeltaPrimePi; //! Generated responses for a single track
348   Double_t* fGenRespPiDeltaPrimePr; //! Generated responses for a single track
349   Double_t* fGenRespMuDeltaPrimeEl; //! Generated responses for a single track
350   Double_t* fGenRespMuDeltaPrimeKa; //! Generated responses for a single track
351   Double_t* fGenRespMuDeltaPrimePi; //! Generated responses for a single track
352   Double_t* fGenRespMuDeltaPrimePr; //! Generated responses for a single track
353   Double_t* fGenRespPrDeltaPrimeEl; //! Generated responses for a single track
354   Double_t* fGenRespPrDeltaPrimeKa; //! Generated responses for a single track
355   Double_t* fGenRespPrDeltaPrimePi; //! Generated responses for a single track
356   Double_t* fGenRespPrDeltaPrimePr; //! Generated responses for a single track
357   /*
358   Double_t* fGenRespElDeltaEl; //! Generated responses for a single track
359   Double_t* fGenRespElDeltaKa; //! Generated responses for a single track
360   Double_t* fGenRespElDeltaPi; //! Generated responses for a single track
361   Double_t* fGenRespElDeltaPr; //! Generated responses for a single track
362   Double_t* fGenRespKaDeltaEl; //! Generated responses for a single track
363   Double_t* fGenRespKaDeltaKa; //! Generated responses for a single track
364   Double_t* fGenRespKaDeltaPi; //! Generated responses for a single track
365   Double_t* fGenRespKaDeltaPr; //! Generated responses for a single track
366   Double_t* fGenRespPiDeltaEl; //! Generated responses for a single track
367   Double_t* fGenRespPiDeltaKa; //! Generated responses for a single track
368   Double_t* fGenRespPiDeltaPi; //! Generated responses for a single track
369   Double_t* fGenRespPiDeltaPr; //! Generated responses for a single track
370   Double_t* fGenRespMuDeltaEl; //! Generated responses for a single track
371   Double_t* fGenRespMuDeltaKa; //! Generated responses for a single track
372   Double_t* fGenRespMuDeltaPi; //! Generated responses for a single track
373   Double_t* fGenRespMuDeltaPr; //! Generated responses for a single track
374   Double_t* fGenRespPrDeltaEl; //! Generated responses for a single track
375   Double_t* fGenRespPrDeltaKa; //! Generated responses for a single track
376   Double_t* fGenRespPrDeltaPi; //! Generated responses for a single track
377   Double_t* fGenRespPrDeltaPr; //! Generated responses for a single track
378   */
379   
380   TAxis* fDeltaPrimeAxis; //! Axis holding the deltaPrime binning
381   TH1D* fhMaxEtaVariation; //! Histo holding the maximum deviation of the eta correction factor from unity vs. 1/dEdx(splines)
382   
383   TH1D* fhEventsProcessed; //! Histo holding the number of processed events (i.e. passing trigger selection, vtx and zvtx cuts
384   TH1D* fhEventsTriggerSel; //! Histo holding the number of events passing trigger selection
385   TH1D* fhEventsTriggerSelVtxCut; //! Histo holding the number of events passing trigger selection and vtx cut
386   
387   TH2D* fhSkippedTracksForSignalGeneration; //! Number of tracks that have been skipped for the signal generation
388   THnSparseD* fhMCgeneratedYieldsPrimaries; //! Histo holding the generated (no reco, no cuts) primary particle yields in considered eta range
389   
390   TH2D* fh2FFJetPtRec;            //! Number of reconstructed jets vs. jetPt and centrality
391   TH2D* fh2FFJetPtGen;            //! Number of generated jets vs. jetPt and centrality
392   
393   TProfile* fh1Xsec;              //! pythia cross section and trials
394   TH1D*     fh1Trials;            //! sum of trials
395   
396   AliCFContainer* fContainerEff; //! Container for efficiency determination
397   
398   THnSparseD* fPtResolution[AliPID::kSPECIES]; //! Pt Resolution for the individual species
399   
400   THnSparseD* fDeDxCheck; //! dEdx check
401   
402   TObjArray* fOutputContainer;  //! output data container
403   
404   TObjArray* fQAContainer; //! output data container for QA
405   
406   AliAnalysisTaskPID(const AliAnalysisTaskPID&); // not implemented
407   AliAnalysisTaskPID& operator=(const AliAnalysisTaskPID&); // not implemented
408   
409   ClassDef(AliAnalysisTaskPID, 19);
410 };
411
412
413 //_____________________________________________________________________________
414 inline Bool_t AliAnalysisTaskPID::FillEfficiencyContainer(const Double_t* values, AliAnalysisTaskPID::EffSteps step,
415                                                           Double_t weight) 
416 {
417   // Fill efficiency container at step "step" with the values
418   
419   if (!fDoEfficiency)
420     return kFALSE;
421   
422   if (!fContainerEff) {
423     AliError("Efficiency container not initialised -> cannot be filled!");
424     return kFALSE;
425   }
426   
427   fContainerEff->Fill(values, step, weight);    
428   
429   return kTRUE;
430 }
431
432
433 //_____________________________________________________________________________
434 inline Bool_t AliAnalysisTaskPID::FillGeneratedYield(const Double_t* values, Double_t weight)
435 {
436   // Fill histos with generated primary yields with provided values
437   
438   if (!fDoPID)
439     return kFALSE;
440   
441   if (!fhMCgeneratedYieldsPrimaries) {
442     AliError("Histo for generated primary yield not initialised -> cannot be filled!");
443     return kFALSE;
444   }
445   
446   fhMCgeneratedYieldsPrimaries->Fill(values, weight);
447     
448   return kTRUE;
449 }
450
451
452 //_____________________________________________________________________________
453 inline Bool_t AliAnalysisTaskPID::FillGenJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm)
454 {
455   if (!fDoPID && !fDoEfficiency)
456     return kFALSE;
457   
458   if (!fh2FFJetPtGen)
459     return kFALSE;
460   
461   if (norm > 0.)
462     fh2FFJetPtGen->Fill(centralityPercentile, jetPt, 1. / norm);
463   else
464     fh2FFJetPtGen->Fill(centralityPercentile, jetPt);
465   
466   return kTRUE;
467 }
468
469
470 //_____________________________________________________________________________
471 inline Bool_t AliAnalysisTaskPID::FillRecJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm)
472 {
473   if (!fDoPID && !fDoEfficiency)
474     return kFALSE;
475   
476   if (!fh2FFJetPtRec)
477     return kFALSE;
478   
479   if (norm > 0.)
480     fh2FFJetPtRec->Fill(centralityPercentile, jetPt, 1. / norm);
481   else
482     fh2FFJetPtRec->Fill(centralityPercentile, jetPt);
483   
484   return kTRUE;
485 }
486
487
488 //_____________________________________________________________________________
489 inline Bool_t AliAnalysisTaskPID::FillPtResolution(Int_t mcID, const Double_t* values)
490 {
491   // Fill histos with pT resolution with provided values
492   
493   if (!fDoPtResolution || mcID < 0 || mcID >= AliPID::kSPECIES)
494     return kFALSE;
495   
496   if (!fPtResolution[mcID]) {
497     AliError(Form("Histo for pT resolution (species: %s) not initialised -> cannot be filled!", AliPID::ParticleName(mcID)));
498     return kFALSE;
499   }
500   
501   fPtResolution[mcID]->Fill(values);
502     
503   return kTRUE;
504 }
505  
506
507 //_____________________________________________________________________________
508 inline Bool_t AliAnalysisTaskPID::IncrementEventCounter(Double_t centralityPercentile, AliAnalysisTaskPID::EventCounterType type)
509 {
510   // Increment the number of events for the given centrality percentile and the corresponding counter type
511   
512   if (type == kTriggerSel) {
513     if (!fhEventsTriggerSel) {
514       AliError("Histogram for number of events (kTriggerSel) not initialised -> cannot be incremented!");
515       return kFALSE;
516     }
517     
518     fhEventsTriggerSel->Fill(centralityPercentile);
519   }
520   else if (type == kTriggerSelAndVtxCut) {
521     if (!fhEventsTriggerSelVtxCut) {
522       AliError("Histogram for number of events (kTriggerSelAndVtxCut) not initialised -> cannot be incremented!");
523       return kFALSE;
524     }
525     
526     fhEventsTriggerSelVtxCut->Fill(centralityPercentile);
527   }
528   else if (type == kTriggerSelAndVtxCutAndZvtxCut) {
529     if (!fhEventsProcessed) {
530       AliError("Histogram for number of events (kTriggerSelAndVtxCutAndZvtxCut) not initialised -> cannot be incremented!");
531       return kFALSE;
532     }
533     
534     fhEventsProcessed->Fill(centralityPercentile);
535   }
536   
537   
538   return kTRUE;
539 };
540
541
542 //_____________________________________________________________________________
543 inline Bool_t AliAnalysisTaskPID::SetEtaAbsCutRange(Double_t lowerLimit, Double_t upperLimit)
544 {
545   if (lowerLimit >= upperLimit) {
546     AliError(Form("Requested lower |eta| cut limit >= upper |eta| cut limit. Old eta cut range will be used (low %f, high %f).",
547                   fEtaAbsCutLow, fEtaAbsCutUp));
548     return kFALSE;
549   }
550   
551   fEtaAbsCutLow = lowerLimit;
552   fEtaAbsCutUp = upperLimit;
553   
554   return kTRUE;
555 };
556
557
558 //_____________________________________________________________________________
559 inline Double_t AliAnalysisTaskPID::GetConvolutedGaussTransitionPar(Int_t index) const
560 {
561   if (index < 0 || index >= 3) {
562     printf("Invalid index %d!\n", index);
563     return -1;
564   }
565   return fConvolutedGaussTransitionPars[index];
566 }
567
568
569 //_____________________________________________________________________________
570 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsTrackPt() const
571 {
572   if (!fFractionHists[AliPID::kPion])
573     return -1;
574   
575   return fFractionHists[AliPID::kPion]->GetNbinsX();
576 }
577
578
579 //_____________________________________________________________________________
580 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsJetPt() const
581 {
582   if (!fFractionHists[AliPID::kPion])
583     return -1;
584   
585   return fFractionHists[AliPID::kPion]->GetNbinsY();
586 }
587
588
589 //_____________________________________________________________________________
590 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsCentrality() const
591 {
592   if (!fFractionHists[AliPID::kPion])
593     return -1;
594   
595   return fFractionHists[AliPID::kPion]->GetNbinsZ();
596 }
597
598
599 //_____________________________________________________________________________
600 inline Double_t AliAnalysisTaskPID::GetCentralityPercentile(AliVEvent* evt) const
601 {
602   if (!evt)
603     return -1;
604   
605   Double_t centralityPercentile = -1.;
606   if (fStoreCentralityPercentile)
607     centralityPercentile = evt->GetCentrality()->GetCentralityPercentile(fCentralityEstimator.Data());
608   
609   return centralityPercentile;
610 }
611
612
613 //_____________________________________________________________________________
614 inline void AliAnalysisTaskPID::PostOutputData()
615 {
616   PostData(1, fOutputContainer);
617   
618   if (fDoEfficiency)
619     PostData(2, fContainerEff);
620   
621   if (fDoPtResolution || fDoDeDxCheck)
622     PostData(3, fQAContainer);
623 }
624
625 #endif