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