1 #ifndef ALI_ANALYSIS_TASK_PID_H
2 #define ALI_ANALYSIS_TASK_PID_H
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.
10 Class written by Benjamin Hess.
11 Contact: bhess@cern.ch
16 class AliAnalysisFilter;
24 class AliTOFPIDResponse;
31 #include "THnSparse.h"
35 #include "AliCentrality.h"
36 #include "AliCFContainer.h"
39 #include "AliAnalysisTaskPIDV0base.h"
41 class AliAnalysisTaskPID : public AliAnalysisTaskPIDV0base {
44 AliAnalysisTaskPID(const char *name);
45 virtual ~AliAnalysisTaskPID();
47 virtual void UserCreateOutputObjects();
48 virtual void UserExec(Option_t *option);
49 virtual void Terminate(const Option_t*);
51 enum ErrorCode { kNoErrors = 1, kWarning = 0, kError = -1};
53 enum dataAxes { kDataMCID = 0, kDataSelectSpecies = 1, kDataPt = 2, kDataDeltaPrimeSpecies = 3, kDataCentrality = 4,
54 kDataJetPt = 5, kDataZ = 6, kDataXi = 7, kDataCharge = 8, kDataTOFpidInfo = 9, kDataNumAxes = 10 };
56 enum genAxes { kGenMCID = 0, kGenSelectSpecies = 1, kGenPt = 2, kGenDeltaPrimeSpecies = 3, kGenCentrality = 4,
57 kGenJetPt = 5, kGenZ = 6, kGenXi = 7, kGenCharge = 8, kGenTOFpidInfo = 9, kGenNumAxes = 10 };
59 enum genYieldAxes { kGenYieldMCID = 0, kGenYieldPt = 1, kGenYieldCentrality = 2, kGenYieldJetPt = 3, kGenYieldZ = 4, kGenYieldXi = 5,
60 kGenYieldCharge = 6, kGenYieldNumAxes = 7 };
62 enum ptResolutionAxes { kPtResJetPt = 0, kPtResGenPt = 1, kPtResRecPt = 2, kPtResCharge = 3, kPtResCentrality = 4, kPtResNumAxes = 5 };
64 enum efficiencyAxes { kEffMCID = 0, kEffTrackPt = 1, kEffTrackEta = 2, kEffTrackCharge = 3, kEffCentrality = 4, kEffJetPt = 5,
65 kEffZ = 6, kEffXi = 7, kEffNumAxes = 8 };
67 enum EffSteps { kStepGenWithGenCuts = 0, kStepRecWithGenCuts = 1, kStepRecWithGenCutsMeasuredObs = 2,
68 kStepRecWithRecCutsMeasuredObs = 3, kStepRecWithRecCutsMeasuredObsPrimaries = 4,
69 kStepRecWithRecCutsMeasuredObsStrangenessScaled = 5, kStepRecWithRecCutsPrimaries = 6, kNumSteps = 7};
71 enum TOFpidInfo { kNoTOFinfo = -2, kNoTOFpid = -1, kTOFpion = 0, kTOFkaon = 1, kTOFproton = 2, kNumTOFspecies = 3,
72 kNumTOFpidInfoBins = 5 };
74 enum EventCounterType { kTriggerSel = 0, kTriggerSelAndVtxCut = 1, kTriggerSelAndVtxCutAndZvtxCut = 2 };
76 static Int_t PDGtoMCID(Int_t pdg);
78 static void GetJetTrackObservables(Double_t trackPt, Double_t jetPt, Double_t& z, Double_t& xi);
80 static Double_t GetMCStrangenessFactorCMS(Int_t motherPDG, Double_t motherGenPt);
81 static Double_t GetMCStrangenessFactorCMS(AliMCEvent* mcEvent, AliMCParticle* daughter);
83 static Bool_t IsSecondaryWithStrangeMotherMC(AliMCEvent* mcEvent, Int_t partLabel);
85 Int_t GetIndexOfChargeAxisData() const
86 { return fStoreAdditionalJetInformation ? kDataCharge : kDataCharge - fgkNumJetAxes; };
87 Int_t GetIndexOfChargeAxisGen() const
88 { return fStoreAdditionalJetInformation ? kGenCharge : kGenCharge - fgkNumJetAxes; };
89 Int_t GetIndexOfChargeAxisGenYield() const
90 { return fStoreAdditionalJetInformation ? kGenYieldCharge : kGenYieldCharge - fgkNumJetAxes; };
92 Int_t GetIndexOfTOFpidInfoAxisData() const
93 { return fStoreAdditionalJetInformation ? kDataTOFpidInfo : kDataTOFpidInfo - fgkNumJetAxes; };
94 Int_t GetIndexOfTOFpidInfoAxisGen() const
95 { return fStoreAdditionalJetInformation ? kGenTOFpidInfo : kGenTOFpidInfo - fgkNumJetAxes; };
97 Bool_t FillXsec(Double_t xsection)
98 { if (!fh1Xsec) return kFALSE; fh1Xsec->Fill("<#sigma>", xsection); return kTRUE; };
99 Bool_t FillPythiaTrials(Double_t avgTrials)
100 { if (!fh1Trials) return kFALSE; fh1Trials->Fill("#sum{ntrials}", avgTrials); return kTRUE; };
102 Bool_t FillEfficiencyContainer(const Double_t* values, EffSteps step, Double_t weight = 1.0);
104 Bool_t FillGeneratedYield(const Double_t* values, Double_t weight = 1.0);
105 Bool_t FillPtResolution(Int_t mcID, const Double_t* values);
106 Bool_t FillGenJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm = -1.);
107 Bool_t FillRecJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm = -1.);
109 Bool_t IncrementEventCounter(Double_t centralityPercentile, EventCounterType type);
111 void PostOutputData();
113 void PrintSettings(Bool_t printSystematicsSettings = kFALSE) const;
115 void PrintSystematicsSettings() const;
117 Bool_t ProcessTrack(const AliVTrack* track, Int_t particlePDGcode, Double_t centralityPercentile, Double_t jetPt) ;
119 ErrorCode GenerateDetectorResponse(ErrorCode errCode, Double_t mean, Double_t sigma, Double_t* responses,
121 Bool_t usePureGaus = kFALSE);
122 ErrorCode SetParamsForConvolutedGaus(Double_t gausMean, Double_t gausSigma);
124 const TString GetCentralityEstimator() const { return fCentralityEstimator; };
125 void SetCentralityEstimator(TString estimator) { fCentralityEstimator = estimator; };
127 Double_t GetCentralityPercentile(AliVEvent* evt) const;
129 inline Double_t GetConvolutedGaussTransitionPar(Int_t index) const;
131 Bool_t SetConvolutedGaussLambdaParameter(Double_t lambda);
133 Bool_t GetInputFromOtherTask() const { return fInputFromOtherTask; };
134 void SetInputFromOtherTask(Bool_t flag) { fInputFromOtherTask = flag; };
136 Bool_t GetDoPID() const { return fDoPID; };
137 void SetDoPID(Bool_t flag) { fDoPID = flag; };
139 Bool_t GetDoEfficiency() const { return fDoEfficiency; };
140 void SetDoEfficiency(Bool_t flag) { fDoEfficiency = flag; };
142 Bool_t GetDoPtResolution() const { return fDoPtResolution; };
143 void SetDoPtResolution(Bool_t flag) { fDoPtResolution = flag; };
145 Bool_t GetStoreCentralityPercentile() const { return fStoreCentralityPercentile; };
146 void SetStoreCentralityPercentile(Bool_t flag) { fStoreCentralityPercentile = flag; };
148 Bool_t GetStoreAdditionalJetInformation() const { return fStoreAdditionalJetInformation; };
149 void SetStoreAdditionalJetInformation(Bool_t flag) { fStoreAdditionalJetInformation = flag; };
151 Bool_t GetUseMCidForGeneration() const { return fUseMCidForGeneration; };
152 void SetUseMCidForGeneration(Bool_t flag) { fUseMCidForGeneration = flag; };
154 Bool_t GetUseConvolutedGaus() const { return fUseConvolutedGaus; };
155 void SetUseConvolutedGaus(Bool_t flag) { fUseConvolutedGaus = flag; };
157 Double_t GetAccuracyNonGaussianTail() const { return fAccuracyNonGaussianTail; };
158 void SetAccuracyNonGaussianTail(Double_t value) { fAccuracyNonGaussianTail = value; };
160 Bool_t GetTakeIntoAccountMuons() const { return fTakeIntoAccountMuons; };
161 void SetTakeIntoAccountMuons(Bool_t flag) { fTakeIntoAccountMuons = flag; };
163 Int_t GetTOFmode() const { return fTOFmode; };
164 void SetTOFmode(Int_t tofMode) { fTOFmode = tofMode; };
166 Bool_t GetUseTPCDefaultPriors() const { return fTPCDefaultPriors; };
167 void SetUseTPCDefaultPriors(Bool_t flag) { fTPCDefaultPriors = flag; };
169 Bool_t GetUsePriors() const { return fUsePriors; };
170 void SetUsePriors(Bool_t flag) { fUsePriors = flag; };
172 Bool_t GetUseITS() const { return fUseITS; };
173 void SetUseITS(Bool_t flag) { fUseITS = flag; };
175 Bool_t GetUseTOF() const { return fUseTOF; };
176 void SetUseTOF(Bool_t flag) { fUseTOF = flag; };
178 Double_t GetEtaAbsCutLow() const { return fEtaAbsCutLow; };
179 Double_t GetEtaAbsCutUp() const { return fEtaAbsCutUp; };
180 Bool_t SetEtaAbsCutRange(Double_t lowerLimit, Double_t upperLimit);
182 Bool_t IsInAcceptedEtaRange(Double_t etaAbs) const { return (etaAbs >= fEtaAbsCutLow && etaAbs <= fEtaAbsCutUp); };
184 Double_t GetSystematicScalingSplinesThreshold() const { return fSystematicScalingSplinesThreshold; };
185 void SetSystematicScalingSplinesThreshold(Double_t threshold) { fSystematicScalingSplinesThreshold = threshold; };
187 Double_t GetSystematicScalingSplinesBelowThreshold() const { return fSystematicScalingSplinesBelowThreshold; };
188 void SetSystematicScalingSplinesBelowThreshold(Double_t scaleFactor)
189 { fSystematicScalingSplinesBelowThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
191 Double_t GetSystematicScalingSplinesAboveThreshold() const { return fSystematicScalingSplinesAboveThreshold; };
192 void SetSystematicScalingSplinesAboveThreshold(Double_t scaleFactor)
193 { fSystematicScalingSplinesAboveThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
195 Double_t GetSystematicScalingEtaCorrectionMomentumThr() const { return fSystematicScalingEtaCorrectionMomentumThr; };
196 void SetSystematicScalingEtaCorrectionMomentumThr(Double_t threshold) { fSystematicScalingEtaCorrectionMomentumThr = threshold; };
198 Double_t GetSystematicScalingEtaCorrectionLowMomenta() const { return fSystematicScalingEtaCorrectionLowMomenta; };
199 void SetSystematicScalingEtaCorrectionLowMomenta(Double_t scaleFactor)
200 { fSystematicScalingEtaCorrectionLowMomenta = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
202 Double_t GetSystematicScalingEtaCorrectionHighMomenta() const { return fSystematicScalingEtaCorrectionHighMomenta; };
203 void SetSystematicScalingEtaCorrectionHighMomenta(Double_t scaleFactor)
204 { fSystematicScalingEtaCorrectionHighMomenta = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
206 Double_t GetSystematicScalingEtaSigmaPara() const { return fSystematicScalingEtaSigmaPara; };
207 void SetSystematicScalingEtaSigmaPara(Double_t scaleFactor)
208 { fSystematicScalingEtaSigmaPara = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
210 Double_t GetSystematicScalingMultCorrection() const { return fSystematicScalingMultCorrection; };
211 void SetSystematicScalingMultCorrection(Double_t scaleFactor)
212 { fSystematicScalingMultCorrection = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
214 void CleanupParticleFractionHistos();
215 Bool_t GetParticleFraction(Double_t trackPt, Double_t jetPt, Double_t multiplicity,
216 AliPID::EParticleType species, Double_t& fraction, Double_t& fractionErrorStat,
217 Double_t& fractionErrorSys) const;
218 Bool_t GetParticleFractions(Double_t trackPt, Double_t jetPt, Double_t centralityPercentile,
219 Double_t* prob, Int_t smearSpeciesByError, Int_t takeIntoAccountSpeciesSysError,
220 Bool_t uniformSystematicError = kFALSE) const;
221 const TH3D* GetParticleFractionHisto(Int_t species, Bool_t sysError = kFALSE) const;
222 Bool_t SetParticleFractionHisto(const TH3D* hist, Int_t species, Bool_t sysError = kFALSE);
223 Int_t GetParticleFractionHistoNbinsTrackPt() const;
224 Int_t GetParticleFractionHistoNbinsJetPt() const;
225 Int_t GetParticleFractionHistoNbinsCentrality() const;
226 Bool_t SetParticleFractionHistosFromFile(const TString filePathName, Bool_t sysError = kFALSE);
227 Int_t GetRandomParticleTypeAccordingToParticleFractions(Double_t trackPt, Double_t jetPt,
228 Double_t centralityPercentile,
230 Bool_t takeIntoAccountSysError = kFALSE) const;
232 TOFpidInfo GetTOFType(const AliVTrack* track, Int_t tofMode) const;
235 void CheckDoAnyStematicStudiesOnTheExpectedSignal();
236 Double_t ConvolutedGaus(const Double_t* xx, const Double_t* par) const;
237 inline Double_t FastGaus(Double_t x, Double_t mean, Double_t sigma) const;
238 inline Double_t FastNormalisedGaus(Double_t x, Double_t mean, Double_t sigma) const;
239 Int_t FindBinWithinRange(TAxis* axis, Double_t value) const;
240 Int_t FindFirstBinAboveIn3dSubset(const TH3* hist, Double_t threshold, Int_t yValue, Int_t zValue) const;
241 Int_t FindLastBinAboveIn3dSubset(const TH3* hist, Double_t threshold, Int_t yValue, Int_t zValue) const;
242 virtual void SetUpGenHist(THnSparse* hist, Double_t* binsPt, Double_t* binsDeltaPrime, Double_t* binsCent, Double_t* binsJetPt) const;
243 virtual void SetUpGenYieldHist(THnSparse* hist, Double_t* binsPt, Double_t* binsCent, Double_t* binsJetPt) const;
244 virtual void SetUpHist(THnSparse* hist, Double_t* binsPt, Double_t* binsDeltaPrime, Double_t* binsCent, Double_t* binsJetPt) const;
245 virtual void SetUpPtResHist(THnSparse* hist, Double_t* binsPt, Double_t* binsJetPt, Double_t* binsCent) const;
246 virtual void SetUpPIDcombined();
248 static const Int_t fgkNumJetAxes; // Number of additional axes for jets
249 static const Double_t fgkEpsilon; // Double_t threshold above zero
250 static const Int_t fgkMaxNumGenEntries; // Maximum number of generated detector responses per track and delta(Prime) and associated species
253 static const Double_t fgkOneOverSqrt2; // = 1. / TMath::Sqrt2();
255 AliPIDCombined* fPIDcombined; //! PID combined object
257 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
259 Bool_t fDoPID; // Only do PID processing (and post the output), if flag is set to kTRUE
260 Bool_t fDoEfficiency; // Only do efficiency processing (and post the output), if flag is set to kTRUE
261 Bool_t fDoPtResolution; // Only do pT resolution processing (and post the output), if flag is set to kTRUE
263 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
264 Bool_t fStoreAdditionalJetInformation; // If set to kTRUE, additional jet information like jetPt, z, xi will be stored in the THnSparses
266 Bool_t fTakeIntoAccountMuons; // Also take into account muons for the generation of the expected response and the most probable PID
267 Bool_t fUseITS; // Use ITS for PID combined probabilities
268 Bool_t fUseTOF; // Use TOF for PID combined probabilities
269 Bool_t fUsePriors; // Use priors for PID combined probabilities
270 Bool_t fTPCDefaultPriors; // Use TPC default priors for PID combined probabilities, if priors are enabled
272 Bool_t fUseMCidForGeneration; // If MC, use MCid instead of PIDcombined to generate the signals
274 Bool_t fUseConvolutedGaus; // Use convoluted gaus to generate detector response instead of pure gaus
275 const Int_t fkConvolutedGausNPar; // Number of parameters for convoluted gaus
276 Double_t fAccuracyNonGaussianTail; // Accuracy of the non-gaussian tail (fraction of the maximum)
277 const Double_t fkDeltaPrimeLowLimit; // Lower deltaPrime limit
278 const Double_t fkDeltaPrimeUpLimit; // Upper deltaPrime limit
279 TF1* fConvolutedGausDeltaPrime; // Gaus convoluted with exponential tail to generate detector response (deltaPrime)
281 Int_t fTOFmode; // TOF mode used for TOF PID info (affects num sigma inclusion/exclusion)
282 Double_t fConvolutedGaussTransitionPars[3]; // Parameter for transition from gaussian parameters to asymmetric shape
283 static const Double_t fgkSigmaReferenceForTransitionPars; // Reference sigma chosen to calculate transition parameters
285 Double_t fEtaAbsCutLow; // Lower cut value on |eta|
286 Double_t fEtaAbsCutUp; // Upper cut value on |eta|
288 // For systematic studies
289 Bool_t fDoAnySystematicStudiesOnTheExpectedSignal; // Internal flag indicating whether any systematic studies are going to be performed
290 Double_t fSystematicScalingSplinesThreshold; // beta-gamma threshold for the systematic spline scale factor
291 Double_t fSystematicScalingSplinesBelowThreshold; // Systematic scale factor for the splines (1. = no systematics) below threshold
292 Double_t fSystematicScalingSplinesAboveThreshold; // Systematic scale factor for the splines (1. = no systematics) above threshold
293 Double_t fSystematicScalingEtaCorrectionMomentumThr; // Momentum threshold for the systematic scale factor for the eta correction (separates low-p from high-p
294 Double_t fSystematicScalingEtaCorrectionLowMomenta; // Systematic scale factor for the eta correction (1. = no systematics) at low momenta
295 Double_t fSystematicScalingEtaCorrectionHighMomenta; // Systematic scale factor for the eta correction (1. = no systematics) at high momenta
296 Double_t fSystematicScalingEtaSigmaPara; // Systematic scale factor for the parametrisation of the relative signal width (1. = no systematics)
297 Double_t fSystematicScalingMultCorrection; // Systematic scale factor for the multiplicity correction (1. = no systematics)
299 TH3D* fFractionHists[AliPID::kSPECIES]; // 3D histos of particle fraction as function with trackPt, jetPt (-1 for inclusive spectra), centralityPercentile (-1 for pp)
300 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)
302 TString fCentralityEstimator; // Estimator for the centrality (e.g. V0A, V0M)
304 THnSparseD* fhPIDdataAll; //! Data histo
306 // Generated response information
307 THnSparseD* fhGenEl; //! Generated response for el
308 THnSparseD* fhGenKa; //! Generated response for ka
309 THnSparseD* fhGenPi; //! Generated response for pi
310 THnSparseD* fhGenMu; //! Generated response for mu
311 THnSparseD* fhGenPr; //! Generated response for pr
313 // Generated responses for a single track
314 Double_t* fGenRespElDeltaPrimeEl; //! Generated responses for a single track
315 Double_t* fGenRespElDeltaPrimeKa; //! Generated responses for a single track
316 Double_t* fGenRespElDeltaPrimePi; //! Generated responses for a single track
317 Double_t* fGenRespElDeltaPrimePr; //! Generated responses for a single track
318 Double_t* fGenRespKaDeltaPrimeEl; //! Generated responses for a single track
319 Double_t* fGenRespKaDeltaPrimeKa; //! Generated responses for a single track
320 Double_t* fGenRespKaDeltaPrimePi; //! Generated responses for a single track
321 Double_t* fGenRespKaDeltaPrimePr; //! Generated responses for a single track
322 Double_t* fGenRespPiDeltaPrimeEl; //! Generated responses for a single track
323 Double_t* fGenRespPiDeltaPrimeKa; //! Generated responses for a single track
324 Double_t* fGenRespPiDeltaPrimePi; //! Generated responses for a single track
325 Double_t* fGenRespPiDeltaPrimePr; //! Generated responses for a single track
326 Double_t* fGenRespMuDeltaPrimeEl; //! Generated responses for a single track
327 Double_t* fGenRespMuDeltaPrimeKa; //! Generated responses for a single track
328 Double_t* fGenRespMuDeltaPrimePi; //! Generated responses for a single track
329 Double_t* fGenRespMuDeltaPrimePr; //! Generated responses for a single track
330 Double_t* fGenRespPrDeltaPrimeEl; //! Generated responses for a single track
331 Double_t* fGenRespPrDeltaPrimeKa; //! Generated responses for a single track
332 Double_t* fGenRespPrDeltaPrimePi; //! Generated responses for a single track
333 Double_t* fGenRespPrDeltaPrimePr; //! Generated responses for a single track
335 Double_t* fGenRespElDeltaEl; //! Generated responses for a single track
336 Double_t* fGenRespElDeltaKa; //! Generated responses for a single track
337 Double_t* fGenRespElDeltaPi; //! Generated responses for a single track
338 Double_t* fGenRespElDeltaPr; //! Generated responses for a single track
339 Double_t* fGenRespKaDeltaEl; //! Generated responses for a single track
340 Double_t* fGenRespKaDeltaKa; //! Generated responses for a single track
341 Double_t* fGenRespKaDeltaPi; //! Generated responses for a single track
342 Double_t* fGenRespKaDeltaPr; //! Generated responses for a single track
343 Double_t* fGenRespPiDeltaEl; //! Generated responses for a single track
344 Double_t* fGenRespPiDeltaKa; //! Generated responses for a single track
345 Double_t* fGenRespPiDeltaPi; //! Generated responses for a single track
346 Double_t* fGenRespPiDeltaPr; //! Generated responses for a single track
347 Double_t* fGenRespMuDeltaEl; //! Generated responses for a single track
348 Double_t* fGenRespMuDeltaKa; //! Generated responses for a single track
349 Double_t* fGenRespMuDeltaPi; //! Generated responses for a single track
350 Double_t* fGenRespMuDeltaPr; //! Generated responses for a single track
351 Double_t* fGenRespPrDeltaEl; //! Generated responses for a single track
352 Double_t* fGenRespPrDeltaKa; //! Generated responses for a single track
353 Double_t* fGenRespPrDeltaPi; //! Generated responses for a single track
354 Double_t* fGenRespPrDeltaPr; //! Generated responses for a single track
357 TH1D* fhEventsProcessed; //! Histo holding the number of processed events (i.e. passing trigger selection, vtx and zvtx cuts
358 TH1D* fhEventsTriggerSel; //! Histo holding the number of events passing trigger selection
359 TH1D* fhEventsTriggerSelVtxCut; //! Histo holding the number of events passing trigger selection and vtx cut
361 TH2D* fhSkippedTracksForSignalGeneration; //! Number of tracks that have been skipped for the signal generation
362 THnSparseD* fhMCgeneratedYieldsPrimaries; //! Histo holding the generated (no reco, no cuts) primary particle yields in considered eta range
364 TH2D* fh2FFJetPtRec; //! Number of reconstructed jets vs. jetPt and centrality
365 TH2D* fh2FFJetPtGen; //! Number of generated jets vs. jetPt and centrality
367 TProfile* fh1Xsec; //! pythia cross section and trials
368 TH1D* fh1Trials; //! sum of trials
370 AliCFContainer* fContainerEff; //! Container for efficiency determination
372 THnSparseD* fPtResolution[AliPID::kSPECIES]; //! Pt Resolution for the individual species
374 TObjArray* fOutputContainer; //! output data container
376 TObjArray* fPtResolutionContainer; //! output data container for pt resolution
378 AliAnalysisTaskPID(const AliAnalysisTaskPID&); // not implemented
379 AliAnalysisTaskPID& operator=(const AliAnalysisTaskPID&); // not implemented
381 ClassDef(AliAnalysisTaskPID, 16);
385 //_____________________________________________________________________________
386 inline Bool_t AliAnalysisTaskPID::FillEfficiencyContainer(const Double_t* values, AliAnalysisTaskPID::EffSteps step,
389 // Fill efficiency container at step "step" with the values
394 if (!fContainerEff) {
395 AliError("Efficiency container not initialised -> cannot be filled!");
399 fContainerEff->Fill(values, step, weight);
405 //_____________________________________________________________________________
406 inline Bool_t AliAnalysisTaskPID::FillGeneratedYield(const Double_t* values, Double_t weight)
408 // Fill histos with generated primary yields with provided values
413 if (!fhMCgeneratedYieldsPrimaries) {
414 AliError("Histo for generated primary yield not initialised -> cannot be filled!");
418 fhMCgeneratedYieldsPrimaries->Fill(values, weight);
424 //_____________________________________________________________________________
425 inline Bool_t AliAnalysisTaskPID::FillGenJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm)
427 if (!fDoPID && !fDoEfficiency)
434 fh2FFJetPtGen->Fill(centralityPercentile, jetPt, 1. / norm);
436 fh2FFJetPtGen->Fill(centralityPercentile, jetPt);
442 //_____________________________________________________________________________
443 inline Bool_t AliAnalysisTaskPID::FillRecJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm)
445 if (!fDoPID && !fDoEfficiency)
452 fh2FFJetPtRec->Fill(centralityPercentile, jetPt, 1. / norm);
454 fh2FFJetPtRec->Fill(centralityPercentile, jetPt);
460 //_____________________________________________________________________________
461 inline Bool_t AliAnalysisTaskPID::FillPtResolution(Int_t mcID, const Double_t* values)
463 // Fill histos with pT resolution with provided values
465 if (!fDoPtResolution || mcID < 0 || mcID >= AliPID::kSPECIES)
468 if (!fPtResolution[mcID]) {
469 AliError(Form("Histo for pT resolution (species: %s) not initialised -> cannot be filled!", AliPID::ParticleName(mcID)));
473 fPtResolution[mcID]->Fill(values);
479 //_____________________________________________________________________________
480 inline Bool_t AliAnalysisTaskPID::IncrementEventCounter(Double_t centralityPercentile, AliAnalysisTaskPID::EventCounterType type)
482 // Increment the number of events for the given centrality percentile and the corresponding counter type
484 if (type == kTriggerSel) {
485 if (!fhEventsTriggerSel) {
486 AliError("Histogram for number of events (kTriggerSel) not initialised -> cannot be incremented!");
490 fhEventsTriggerSel->Fill(centralityPercentile);
492 else if (type == kTriggerSelAndVtxCut) {
493 if (!fhEventsTriggerSelVtxCut) {
494 AliError("Histogram for number of events (kTriggerSelAndVtxCut) not initialised -> cannot be incremented!");
498 fhEventsTriggerSelVtxCut->Fill(centralityPercentile);
500 else if (type == kTriggerSelAndVtxCutAndZvtxCut) {
501 if (!fhEventsProcessed) {
502 AliError("Histogram for number of events (kTriggerSelAndVtxCutAndZvtxCut) not initialised -> cannot be incremented!");
506 fhEventsProcessed->Fill(centralityPercentile);
514 //_____________________________________________________________________________
515 inline Bool_t AliAnalysisTaskPID::SetEtaAbsCutRange(Double_t lowerLimit, Double_t upperLimit)
517 if (lowerLimit >= upperLimit) {
518 AliError(Form("Requested lower |eta| cut limit >= upper |eta| cut limit. Old eta cut range will be used (low %f, high %f).",
519 fEtaAbsCutLow, fEtaAbsCutUp));
523 fEtaAbsCutLow = lowerLimit;
524 fEtaAbsCutUp = upperLimit;
530 //_____________________________________________________________________________
531 inline Double_t AliAnalysisTaskPID::GetConvolutedGaussTransitionPar(Int_t index) const
533 if (index < 0 || index >= 3) {
534 printf("Invalid index %d!\n", index);
537 return fConvolutedGaussTransitionPars[index];
541 //_____________________________________________________________________________
542 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsTrackPt() const
544 if (!fFractionHists[AliPID::kPion])
547 return fFractionHists[AliPID::kPion]->GetNbinsX();
551 //_____________________________________________________________________________
552 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsJetPt() const
554 if (!fFractionHists[AliPID::kPion])
557 return fFractionHists[AliPID::kPion]->GetNbinsY();
561 //_____________________________________________________________________________
562 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsCentrality() const
564 if (!fFractionHists[AliPID::kPion])
567 return fFractionHists[AliPID::kPion]->GetNbinsZ();
571 //_____________________________________________________________________________
572 inline Double_t AliAnalysisTaskPID::GetCentralityPercentile(AliVEvent* evt) const
577 Double_t centralityPercentile = -1.;
578 if (fStoreCentralityPercentile)
579 centralityPercentile = evt->GetCentrality()->GetCentralityPercentile(fCentralityEstimator.Data());
581 return centralityPercentile;
585 //_____________________________________________________________________________
586 inline void AliAnalysisTaskPID::PostOutputData()
588 PostData(1, fOutputContainer);
591 PostData(2, fContainerEff);
594 PostData(3, fPtResolutionContainer);