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;
32 #include "THnSparse.h"
36 #include "AliCentrality.h"
37 #include "AliCFContainer.h"
40 #include "AliAnalysisTaskPIDV0base.h"
42 class AliAnalysisTaskPID : public AliAnalysisTaskPIDV0base {
45 AliAnalysisTaskPID(const char *name);
46 virtual ~AliAnalysisTaskPID();
48 virtual void UserCreateOutputObjects();
49 virtual void UserExec(Option_t *option);
50 virtual void Terminate(const Option_t*);
52 enum ErrorCode { kNoErrors = 1, kWarning = 0, kError = -1};
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 };
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 };
60 enum genYieldAxes { kGenYieldMCID = 0, kGenYieldPt = 1, kGenYieldCentrality = 2, kGenYieldJetPt = 3, kGenYieldZ = 4, kGenYieldXi = 5,
61 kGenYieldCharge = 6, kGenYieldNumAxes = 7 };
63 enum ptResolutionAxes { kPtResJetPt = 0, kPtResGenPt = 1, kPtResRecPt = 2, kPtResCharge = 3, kPtResCentrality = 4, kPtResNumAxes = 5 };
65 enum qaSharedClsAxes { kQASharedClsJetPt = 0, kQASharedClsPt = 1, kQASharedClsNumSharedCls = 2, kQASharedClsPadRow = 3,
66 kQASharedClsNumAxes = 4 };
68 enum dEdxCheckAxes { kDeDxCheckPID = 0, kDeDxCheckP = 1, kDeDxCheckJetPt = 2, kDeDxCheckEtaAbs = 3 , kDeDxCheckDeDx = 4,
69 kDeDxCheckNumAxes = 5 };
71 enum efficiencyAxes { kEffMCID = 0, kEffTrackPt = 1, kEffTrackEta = 2, kEffTrackCharge = 3, kEffCentrality = 4, kEffJetPt = 5,
72 kEffZ = 6, kEffXi = 7, kEffNumAxes = 8 };
74 enum EffSteps { kStepGenWithGenCuts = 0, kStepRecWithGenCuts = 1, kStepRecWithGenCutsMeasuredObs = 2,
75 kStepRecWithRecCutsMeasuredObs = 3, kStepRecWithRecCutsMeasuredObsPrimaries = 4,
76 kStepRecWithRecCutsMeasuredObsStrangenessScaled = 5, kStepRecWithRecCutsPrimaries = 6, kNumSteps = 7};
78 enum TOFpidInfo { kNoTOFinfo = -2, kNoTOFpid = -1, kTOFpion = 0, kTOFkaon = 1, kTOFproton = 2, kNumTOFspecies = 3,
79 kNumTOFpidInfoBins = 5 };
81 enum EventCounterType { kTriggerSel = 0, kTriggerSelAndVtxCut = 1, kTriggerSelAndVtxCutAndZvtxCutNoPileUpRejection = 2,
82 kTriggerSelAndVtxCutAndZvtxCut = 3 };
84 static Int_t PDGtoMCID(Int_t pdg);
86 static void GetJetTrackObservables(Double_t trackPt, Double_t jetPt, Double_t& z, Double_t& xi);
88 static Double_t GetMCStrangenessFactorCMS(Int_t motherPDG, Double_t motherGenPt);
89 static Double_t GetMCStrangenessFactorCMS(AliMCEvent* mcEvent, AliMCParticle* daughter);
91 static Bool_t IsSecondaryWithStrangeMotherMC(AliMCEvent* mcEvent, Int_t partLabel);
93 virtual void ConfigureTaskForCurrentEvent(AliVEvent* event);
95 Int_t GetIndexOfChargeAxisData() const
96 { return fStoreAdditionalJetInformation ? kDataCharge : kDataCharge - fgkNumJetAxes; };
97 Int_t GetIndexOfChargeAxisGen() const
98 { return fStoreAdditionalJetInformation ? kGenCharge : kGenCharge - fgkNumJetAxes; };
99 Int_t GetIndexOfChargeAxisGenYield() const
100 { return fStoreAdditionalJetInformation ? kGenYieldCharge : kGenYieldCharge - fgkNumJetAxes; };
102 Int_t GetIndexOfTOFpidInfoAxisData() const
103 { return fStoreAdditionalJetInformation ? kDataTOFpidInfo : kDataTOFpidInfo - fgkNumJetAxes; };
104 Int_t GetIndexOfTOFpidInfoAxisGen() const
105 { return fStoreAdditionalJetInformation ? kGenTOFpidInfo : kGenTOFpidInfo - fgkNumJetAxes; };
107 Bool_t FillXsec(Double_t xsection)
108 { if (!fh1Xsec) return kFALSE; fh1Xsec->Fill("<#sigma>", xsection); return kTRUE; };
109 Bool_t FillPythiaTrials(Double_t avgTrials)
110 { if (!fh1Trials) return kFALSE; fh1Trials->Fill("#sum{ntrials}", avgTrials); return kTRUE; };
112 Bool_t FillEfficiencyContainer(const Double_t* values, EffSteps step, Double_t weight = 1.0);
114 Bool_t FillGeneratedYield(const Double_t* values, Double_t weight = 1.0);
115 Bool_t FillPtResolution(Int_t mcID, const Double_t* values);
116 Bool_t FillGenJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm = -1.);
117 Bool_t FillRecJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm = -1.);
119 Bool_t IncrementEventCounter(Double_t centralityPercentile, EventCounterType type);
121 void PostOutputData();
123 void PrintSettings(Bool_t printSystematicsSettings = kFALSE) const;
125 void PrintSystematicsSettings() const;
127 Bool_t ProcessTrack(const AliVTrack* track, Int_t particlePDGcode, Double_t centralityPercentile, Double_t jetPt) ;
129 ErrorCode GenerateDetectorResponse(ErrorCode errCode, Double_t mean, Double_t sigma, Double_t* responses,
131 Bool_t usePureGaus = kFALSE);
132 ErrorCode SetParamsForConvolutedGaus(Double_t gausMean, Double_t gausSigma);
134 const TString GetCentralityEstimator() const { return fCentralityEstimator; };
135 void SetCentralityEstimator(TString estimator) { fCentralityEstimator = estimator; };
137 Double_t GetCentralityPercentile(AliVEvent* evt) const;
139 inline Double_t GetConvolutedGaussTransitionPar(Int_t index) const;
141 Bool_t SetConvolutedGaussLambdaParameter(Double_t lambda);
143 Bool_t GetInputFromOtherTask() const { return fInputFromOtherTask; };
144 void SetInputFromOtherTask(Bool_t flag) { fInputFromOtherTask = flag; };
146 Bool_t GetDoPID() const { return fDoPID; };
147 void SetDoPID(Bool_t flag) { fDoPID = flag; };
149 Bool_t GetDoEfficiency() const { return fDoEfficiency; };
150 void SetDoEfficiency(Bool_t flag) { fDoEfficiency = flag; };
152 Bool_t GetDoPtResolution() const { return fDoPtResolution; };
153 void SetDoPtResolution(Bool_t flag) { fDoPtResolution = flag; };
155 Bool_t GetDoDeDxCheck() const { return fDoDeDxCheck; };
156 void SetDoDeDxCheck(Bool_t flag) { fDoDeDxCheck = flag; };
158 Bool_t GetStoreCentralityPercentile() const { return fStoreCentralityPercentile; };
159 void SetStoreCentralityPercentile(Bool_t flag) { fStoreCentralityPercentile = flag; };
161 Bool_t GetStoreAdditionalJetInformation() const { return fStoreAdditionalJetInformation; };
162 void SetStoreAdditionalJetInformation(Bool_t flag) { fStoreAdditionalJetInformation = flag; };
164 Bool_t GetUseMCidForGeneration() const { return fUseMCidForGeneration; };
165 void SetUseMCidForGeneration(Bool_t flag) { fUseMCidForGeneration = flag; };
167 Bool_t GetUseConvolutedGaus() const { return fUseConvolutedGaus; };
168 void SetUseConvolutedGaus(Bool_t flag) { fUseConvolutedGaus = flag; };
170 Double_t GetAccuracyNonGaussianTail() const { return fAccuracyNonGaussianTail; };
171 void SetAccuracyNonGaussianTail(Double_t value) { fAccuracyNonGaussianTail = value; };
173 Bool_t GetTakeIntoAccountMuons() const { return fTakeIntoAccountMuons; };
174 void SetTakeIntoAccountMuons(Bool_t flag) { fTakeIntoAccountMuons = flag; };
176 Int_t GetTOFmode() const { return fTOFmode; };
177 void SetTOFmode(Int_t tofMode) { fTOFmode = tofMode; };
179 Bool_t GetUseTPCDefaultPriors() const { return fTPCDefaultPriors; };
180 void SetUseTPCDefaultPriors(Bool_t flag) { fTPCDefaultPriors = flag; };
182 Bool_t GetUsePriors() const { return fUsePriors; };
183 void SetUsePriors(Bool_t flag) { fUsePriors = flag; };
185 Bool_t GetUseITS() const { return fUseITS; };
186 void SetUseITS(Bool_t flag) { fUseITS = flag; };
188 Bool_t GetUseTOF() const { return fUseTOF; };
189 void SetUseTOF(Bool_t flag) { fUseTOF = flag; };
191 Double_t GetEtaAbsCutLow() const { return fEtaAbsCutLow; };
192 Double_t GetEtaAbsCutUp() const { return fEtaAbsCutUp; };
193 Bool_t SetEtaAbsCutRange(Double_t lowerLimit, Double_t upperLimit);
195 Bool_t IsInAcceptedEtaRange(Double_t etaAbs) const { return (etaAbs >= fEtaAbsCutLow && etaAbs <= fEtaAbsCutUp); };
197 AliAnalysisTaskPIDV0base::PileUpRejectionType GetPileUpRejectionType() const { return fPileUpRejectionType; };
198 void SetPileUpRejectionType(AliAnalysisTaskPIDV0base::PileUpRejectionType newType) { fPileUpRejectionType = newType; };
200 Double_t GetSystematicScalingSplinesThreshold() const { return fSystematicScalingSplinesThreshold; };
201 void SetSystematicScalingSplinesThreshold(Double_t threshold) { fSystematicScalingSplinesThreshold = threshold; };
203 Double_t GetSystematicScalingSplinesBelowThreshold() const { return fSystematicScalingSplinesBelowThreshold; };
204 void SetSystematicScalingSplinesBelowThreshold(Double_t scaleFactor)
205 { fSystematicScalingSplinesBelowThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
207 Double_t GetSystematicScalingSplinesAboveThreshold() const { return fSystematicScalingSplinesAboveThreshold; };
208 void SetSystematicScalingSplinesAboveThreshold(Double_t scaleFactor)
209 { fSystematicScalingSplinesAboveThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
211 Double_t GetSystematicScalingEtaCorrectionMomentumThr() const { return fSystematicScalingEtaCorrectionMomentumThr; };
212 void SetSystematicScalingEtaCorrectionMomentumThr(Double_t threshold) { fSystematicScalingEtaCorrectionMomentumThr = threshold; };
214 Double_t GetSystematicScalingEtaCorrectionLowMomenta() const { return fSystematicScalingEtaCorrectionLowMomenta; };
215 void SetSystematicScalingEtaCorrectionLowMomenta(Double_t scaleFactor)
216 { fSystematicScalingEtaCorrectionLowMomenta = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
218 Double_t GetSystematicScalingEtaCorrectionHighMomenta() const { return fSystematicScalingEtaCorrectionHighMomenta; };
219 void SetSystematicScalingEtaCorrectionHighMomenta(Double_t scaleFactor)
220 { fSystematicScalingEtaCorrectionHighMomenta = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
222 Double_t GetSystematicScalingEtaSigmaParaThreshold() const { return fSystematicScalingEtaSigmaParaThreshold; };
223 void SetSystematicScalingEtaSigmaParaThreshold(Double_t threshold) { fSystematicScalingEtaSigmaParaThreshold = threshold; };
225 Double_t GetSystematicScalingEtaSigmaParaBelowThreshold() const { return fSystematicScalingEtaSigmaParaBelowThreshold; };
226 void SetSystematicScalingEtaSigmaParaBelowThreshold(Double_t scaleFactor)
227 { fSystematicScalingEtaSigmaParaBelowThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
229 Double_t GetSystematicScalingEtaSigmaParaAboveThreshold() const { return fSystematicScalingEtaSigmaParaAboveThreshold; };
230 void SetSystematicScalingEtaSigmaParaAboveThreshold(Double_t scaleFactor)
231 { fSystematicScalingEtaSigmaParaAboveThreshold = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
233 Double_t GetSystematicScalingMultCorrection() const { return fSystematicScalingMultCorrection; };
234 void SetSystematicScalingMultCorrection(Double_t scaleFactor)
235 { fSystematicScalingMultCorrection = scaleFactor; CheckDoAnyStematicStudiesOnTheExpectedSignal(); };
237 Double_t GetMaxEtaVariation(Double_t dEdxSplines);
238 Bool_t CalculateMaxEtaVariationMapFromPIDResponse();
240 void CleanupParticleFractionHistos();
241 Bool_t GetParticleFraction(Double_t trackPt, Double_t jetPt, Double_t multiplicity,
242 AliPID::EParticleType species, Double_t& fraction, Double_t& fractionErrorStat,
243 Double_t& fractionErrorSys) const;
244 Bool_t GetParticleFractions(Double_t trackPt, Double_t jetPt, Double_t centralityPercentile,
245 Double_t* prob, Int_t smearSpeciesByError, Int_t takeIntoAccountSpeciesSysError,
246 Bool_t uniformSystematicError = kFALSE) const;
247 const TH3D* GetParticleFractionHisto(Int_t species, Bool_t sysError = kFALSE) const;
248 Bool_t SetParticleFractionHisto(const TH3D* hist, Int_t species, Bool_t sysError = kFALSE);
249 Int_t GetParticleFractionHistoNbinsTrackPt() const;
250 Int_t GetParticleFractionHistoNbinsJetPt() const;
251 Int_t GetParticleFractionHistoNbinsCentrality() const;
252 Bool_t SetParticleFractionHistosFromFile(const TString filePathName, Bool_t sysError = kFALSE);
253 Int_t GetRandomParticleTypeAccordingToParticleFractions(Double_t trackPt, Double_t jetPt,
254 Double_t centralityPercentile,
256 Bool_t takeIntoAccountSysError = kFALSE) const;
258 TOFpidInfo GetTOFType(const AliVTrack* track, Int_t tofMode) const;
261 void CheckDoAnyStematicStudiesOnTheExpectedSignal();
262 Double_t ConvolutedGaus(const Double_t* xx, const Double_t* par) const;
263 inline Double_t FastGaus(Double_t x, Double_t mean, Double_t sigma) const;
264 inline Double_t FastNormalisedGaus(Double_t x, Double_t mean, Double_t sigma) const;
265 Int_t FindBinWithinRange(TAxis* axis, Double_t value) const;
266 Int_t FindFirstBinAboveIn3dSubset(const TH3* hist, Double_t threshold, Int_t yValue, Int_t zValue) const;
267 Int_t FindLastBinAboveIn3dSubset(const TH3* hist, Double_t threshold, Int_t yValue, Int_t zValue) const;
268 virtual void SetUpGenHist(THnSparse* hist, Double_t* binsPt, Double_t* binsDeltaPrime, Double_t* binsCent, Double_t* binsJetPt) const;
269 virtual void SetUpGenYieldHist(THnSparse* hist, Double_t* binsPt, Double_t* binsCent, Double_t* binsJetPt) const;
270 virtual void SetUpHist(THnSparse* hist, Double_t* binsPt, Double_t* binsDeltaPrime, Double_t* binsCent, Double_t* binsJetPt) const;
271 virtual void SetUpPtResHist(THnSparse* hist, Double_t* binsPt, Double_t* binsJetPt, Double_t* binsCent) const;
272 virtual void SetUpSharedClsHist(THnSparse* hist, Double_t* binsPt, Double_t* binsJetPt) const;
273 virtual void SetUpDeDxCheckHist(THnSparse* hist, const Double_t* binsPt, const Double_t* binsJetPt, const Double_t* binsEtaAbs) const;
274 virtual void SetUpPIDcombined();
276 static const Int_t fgkNumJetAxes; // Number of additional axes for jets
277 static const Double_t fgkEpsilon; // Double_t threshold above zero
278 static const Int_t fgkMaxNumGenEntries; // Maximum number of generated detector responses per track and delta(Prime) and associated species
280 static const Double_t fgkOneOverSqrt2; // = 1. / TMath::Sqrt2();
283 Int_t fRun; // Current run number
284 AliPIDCombined* fPIDcombined; //! PID combined object
286 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
288 Bool_t fDoPID; // Only do PID processing (and post the output), if flag is set to kTRUE
289 Bool_t fDoEfficiency; // Only do efficiency processing (and post the output), if flag is set to kTRUE
290 Bool_t fDoPtResolution; // Only do pT resolution processing (and post the output), if flag is set to kTRUE
291 Bool_t fDoDeDxCheck; // Only check dEdx, if flag set to kTRUE
293 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
294 Bool_t fStoreAdditionalJetInformation; // If set to kTRUE, additional jet information like jetPt, z, xi will be stored in the THnSparses
296 Bool_t fTakeIntoAccountMuons; // Also take into account muons for the generation of the expected response and the most probable PID
297 Bool_t fUseITS; // Use ITS for PID combined probabilities
298 Bool_t fUseTOF; // Use TOF for PID combined probabilities
299 Bool_t fUsePriors; // Use priors for PID combined probabilities
300 Bool_t fTPCDefaultPriors; // Use TPC default priors for PID combined probabilities, if priors are enabled
302 Bool_t fUseMCidForGeneration; // If MC, use MCid instead of PIDcombined to generate the signals
304 Bool_t fUseConvolutedGaus; // Use convoluted gaus to generate detector response instead of pure gaus
305 const Int_t fkConvolutedGausNPar; // Number of parameters for convoluted gaus
306 Double_t fAccuracyNonGaussianTail; // Accuracy of the non-gaussian tail (fraction of the maximum)
307 const Double_t fkDeltaPrimeLowLimit; // Lower deltaPrime limit
308 const Double_t fkDeltaPrimeUpLimit; // Upper deltaPrime limit
309 TF1* fConvolutedGausDeltaPrime; // Gaus convoluted with exponential tail to generate detector response (deltaPrime)
311 Int_t fTOFmode; // TOF mode used for TOF PID info (affects num sigma inclusion/exclusion)
312 Double_t fConvolutedGaussTransitionPars[3]; // Parameter for transition from gaussian parameters to asymmetric shape
313 static const Double_t fgkSigmaReferenceForTransitionPars; // Reference sigma chosen to calculate transition parameters
315 Double_t fEtaAbsCutLow; // Lower cut value on |eta|
316 Double_t fEtaAbsCutUp; // Upper cut value on |eta|
318 AliAnalysisTaskPIDV0base::PileUpRejectionType fPileUpRejectionType; // Which pile-up rejection is used (if any)
320 // For systematic studies
321 Bool_t fDoAnySystematicStudiesOnTheExpectedSignal; // Internal flag indicating whether any systematic studies are going to be performed
322 Double_t fSystematicScalingSplinesThreshold; // beta-gamma threshold for the systematic spline scale factor
323 Double_t fSystematicScalingSplinesBelowThreshold; // Systematic scale factor for the splines (1. = no systematics) below threshold
324 Double_t fSystematicScalingSplinesAboveThreshold; // Systematic scale factor for the splines (1. = no systematics) above threshold
325 Double_t fSystematicScalingEtaCorrectionMomentumThr; // Momentum threshold for the systematic scale factor for the eta correction (separates low-p from high-p
326 Double_t fSystematicScalingEtaCorrectionLowMomenta; // Systematic scale factor for the eta correction (1. = no systematics) at low momenta
327 Double_t fSystematicScalingEtaCorrectionHighMomenta; // Systematic scale factor for the eta correction (1. = no systematics) at high momenta
329 Double_t fSystematicScalingEtaSigmaParaThreshold; // dEdx threshold for the systematic scale factor for the parametrisation of the relative signal width
330 Double_t fSystematicScalingEtaSigmaParaBelowThreshold; // Systematic scale factor for the parametrisation of the relative signal width (1. = no systematics) below threshold
331 Double_t fSystematicScalingEtaSigmaParaAboveThreshold; // Systematic scale factor for the parametrisation of the relative signal width (1. = no systematics) above threshold
332 Double_t fSystematicScalingMultCorrection; // Systematic scale factor for the multiplicity correction (1. = no systematics)
334 TH3D* fFractionHists[AliPID::kSPECIES]; // 3D histos of particle fraction as function with trackPt, jetPt (-1 for inclusive spectra), centralityPercentile (-1 for pp)
335 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)
337 TString fCentralityEstimator; // Estimator for the centrality (e.g. V0A, V0M)
339 THnSparseD* fhPIDdataAll; //! Data histo
341 // Generated response information
342 THnSparseD* fhGenEl; //! Generated response for el
343 THnSparseD* fhGenKa; //! Generated response for ka
344 THnSparseD* fhGenPi; //! Generated response for pi
345 THnSparseD* fhGenMu; //! Generated response for mu
346 THnSparseD* fhGenPr; //! Generated response for pr
348 // Generated responses for a single track
349 Double_t* fGenRespElDeltaPrimeEl; //! Generated responses for a single track
350 Double_t* fGenRespElDeltaPrimeKa; //! Generated responses for a single track
351 Double_t* fGenRespElDeltaPrimePi; //! Generated responses for a single track
352 Double_t* fGenRespElDeltaPrimePr; //! Generated responses for a single track
353 Double_t* fGenRespKaDeltaPrimeEl; //! Generated responses for a single track
354 Double_t* fGenRespKaDeltaPrimeKa; //! Generated responses for a single track
355 Double_t* fGenRespKaDeltaPrimePi; //! Generated responses for a single track
356 Double_t* fGenRespKaDeltaPrimePr; //! Generated responses for a single track
357 Double_t* fGenRespPiDeltaPrimeEl; //! Generated responses for a single track
358 Double_t* fGenRespPiDeltaPrimeKa; //! Generated responses for a single track
359 Double_t* fGenRespPiDeltaPrimePi; //! Generated responses for a single track
360 Double_t* fGenRespPiDeltaPrimePr; //! Generated responses for a single track
361 Double_t* fGenRespMuDeltaPrimeEl; //! Generated responses for a single track
362 Double_t* fGenRespMuDeltaPrimeKa; //! Generated responses for a single track
363 Double_t* fGenRespMuDeltaPrimePi; //! Generated responses for a single track
364 Double_t* fGenRespMuDeltaPrimePr; //! Generated responses for a single track
365 Double_t* fGenRespPrDeltaPrimeEl; //! Generated responses for a single track
366 Double_t* fGenRespPrDeltaPrimeKa; //! Generated responses for a single track
367 Double_t* fGenRespPrDeltaPrimePi; //! Generated responses for a single track
368 Double_t* fGenRespPrDeltaPrimePr; //! Generated responses for a single track
370 Double_t* fGenRespElDeltaEl; //! Generated responses for a single track
371 Double_t* fGenRespElDeltaKa; //! Generated responses for a single track
372 Double_t* fGenRespElDeltaPi; //! Generated responses for a single track
373 Double_t* fGenRespElDeltaPr; //! Generated responses for a single track
374 Double_t* fGenRespKaDeltaEl; //! Generated responses for a single track
375 Double_t* fGenRespKaDeltaKa; //! Generated responses for a single track
376 Double_t* fGenRespKaDeltaPi; //! Generated responses for a single track
377 Double_t* fGenRespKaDeltaPr; //! Generated responses for a single track
378 Double_t* fGenRespPiDeltaEl; //! Generated responses for a single track
379 Double_t* fGenRespPiDeltaKa; //! Generated responses for a single track
380 Double_t* fGenRespPiDeltaPi; //! Generated responses for a single track
381 Double_t* fGenRespPiDeltaPr; //! Generated responses for a single track
382 Double_t* fGenRespMuDeltaEl; //! Generated responses for a single track
383 Double_t* fGenRespMuDeltaKa; //! Generated responses for a single track
384 Double_t* fGenRespMuDeltaPi; //! Generated responses for a single track
385 Double_t* fGenRespMuDeltaPr; //! Generated responses for a single track
386 Double_t* fGenRespPrDeltaEl; //! Generated responses for a single track
387 Double_t* fGenRespPrDeltaKa; //! Generated responses for a single track
388 Double_t* fGenRespPrDeltaPi; //! Generated responses for a single track
389 Double_t* fGenRespPrDeltaPr; //! Generated responses for a single track
392 TAxis* fDeltaPrimeAxis; //! Axis holding the deltaPrime binning
393 TH1D* fhMaxEtaVariation; //! Histo holding the maximum deviation of the eta correction factor from unity vs. 1/dEdx(splines)
395 TH1D* fhEventsProcessed; //! Histo holding the number of processed events (i.e. passing trigger selection, vtx and zvtx cuts and (if enabled) pile-up rejection)
396 TH1D* fhEventsTriggerSel; //! Histo holding the number of events passing trigger selection
397 TH1D* fhEventsTriggerSelVtxCut; //! Histo holding the number of events passing trigger selection and vtx cut
398 TH1D* fhEventsProcessedNoPileUpRejection; //! Histo holding the number of processed events before pile-up rejection
400 THnSparseD* fhMCgeneratedYieldsPrimaries; //! Histo holding the generated (no reco, no cuts) primary particle yields in considered eta range
402 TH2D* fh2FFJetPtRec; //! Number of reconstructed jets vs. jetPt and centrality
403 TH2D* fh2FFJetPtGen; //! Number of generated jets vs. jetPt and centrality
405 TProfile* fh1Xsec; //! pythia cross section and trials
406 TH1D* fh1Trials; //! sum of trials
408 AliCFContainer* fContainerEff; //! Container for efficiency determination
410 THnSparseD* fPtResolution[AliPID::kSPECIES]; //! Pt Resolution for the individual species
411 THnSparseD* fQASharedCls; //! QA for shared clusters
413 THnSparseD* fDeDxCheck; //! dEdx check
415 TObjArray* fOutputContainer; //! output data container
417 TObjArray* fQAContainer; //! output data container for QA
419 AliAnalysisTaskPID(const AliAnalysisTaskPID&); // not implemented
420 AliAnalysisTaskPID& operator=(const AliAnalysisTaskPID&); // not implemented
422 ClassDef(AliAnalysisTaskPID, 21);
426 //_____________________________________________________________________________
427 inline Bool_t AliAnalysisTaskPID::FillEfficiencyContainer(const Double_t* values, AliAnalysisTaskPID::EffSteps step,
430 // Fill efficiency container at step "step" with the values
435 if (!fContainerEff) {
436 AliError("Efficiency container not initialised -> cannot be filled!");
440 fContainerEff->Fill(values, step, weight);
446 //_____________________________________________________________________________
447 inline Bool_t AliAnalysisTaskPID::FillGeneratedYield(const Double_t* values, Double_t weight)
449 // Fill histos with generated primary yields with provided values
454 if (!fhMCgeneratedYieldsPrimaries) {
455 AliError("Histo for generated primary yield not initialised -> cannot be filled!");
459 fhMCgeneratedYieldsPrimaries->Fill(values, weight);
465 //_____________________________________________________________________________
466 inline Bool_t AliAnalysisTaskPID::FillGenJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm)
468 if (!fDoPID && !fDoEfficiency)
475 fh2FFJetPtGen->Fill(centralityPercentile, jetPt, 1. / norm);
477 fh2FFJetPtGen->Fill(centralityPercentile, jetPt);
483 //_____________________________________________________________________________
484 inline Bool_t AliAnalysisTaskPID::FillRecJets(Double_t centralityPercentile, Double_t jetPt, Double_t norm)
486 if (!fDoPID && !fDoEfficiency)
493 fh2FFJetPtRec->Fill(centralityPercentile, jetPt, 1. / norm);
495 fh2FFJetPtRec->Fill(centralityPercentile, jetPt);
501 //_____________________________________________________________________________
502 inline Bool_t AliAnalysisTaskPID::FillPtResolution(Int_t mcID, const Double_t* values)
504 // Fill histos with pT resolution with provided values
506 if (!fDoPtResolution || mcID < 0 || mcID >= AliPID::kSPECIES)
509 if (!fPtResolution[mcID]) {
510 AliError(Form("Histo for pT resolution (species: %s) not initialised -> cannot be filled!", AliPID::ParticleName(mcID)));
514 fPtResolution[mcID]->Fill(values);
520 //_____________________________________________________________________________
521 inline Bool_t AliAnalysisTaskPID::IncrementEventCounter(Double_t centralityPercentile, AliAnalysisTaskPID::EventCounterType type)
523 // Increment the number of events for the given centrality percentile and the corresponding counter type
525 if (type == kTriggerSel) {
526 if (!fhEventsTriggerSel) {
527 AliError("Histogram for number of events (kTriggerSel) not initialised -> cannot be incremented!");
531 fhEventsTriggerSel->Fill(centralityPercentile);
533 else if (type == kTriggerSelAndVtxCut) {
534 if (!fhEventsTriggerSelVtxCut) {
535 AliError("Histogram for number of events (kTriggerSelAndVtxCut) not initialised -> cannot be incremented!");
539 fhEventsTriggerSelVtxCut->Fill(centralityPercentile);
541 else if (type == kTriggerSelAndVtxCutAndZvtxCut) {
542 if (!fhEventsProcessed) {
543 AliError("Histogram for number of events (kTriggerSelAndVtxCutAndZvtxCut) not initialised -> cannot be incremented!");
547 fhEventsProcessed->Fill(centralityPercentile);
549 else if (type == kTriggerSelAndVtxCutAndZvtxCutNoPileUpRejection) {
550 if (!fhEventsProcessedNoPileUpRejection) {
551 AliError("Histogram for number of events (kTriggerSelAndVtxCutAndZvtxCutNoPileUpRejection) not initialised -> cannot be incremented!");
555 fhEventsProcessedNoPileUpRejection->Fill(centralityPercentile);
563 //_____________________________________________________________________________
564 inline Bool_t AliAnalysisTaskPID::SetEtaAbsCutRange(Double_t lowerLimit, Double_t upperLimit)
566 if (lowerLimit >= upperLimit) {
567 AliError(Form("Requested lower |eta| cut limit >= upper |eta| cut limit. Old eta cut range will be used (low %f, high %f).",
568 fEtaAbsCutLow, fEtaAbsCutUp));
572 fEtaAbsCutLow = lowerLimit;
573 fEtaAbsCutUp = upperLimit;
579 //_____________________________________________________________________________
580 inline Double_t AliAnalysisTaskPID::GetConvolutedGaussTransitionPar(Int_t index) const
582 if (index < 0 || index >= 3) {
583 printf("Invalid index %d!\n", index);
586 return fConvolutedGaussTransitionPars[index];
590 //_____________________________________________________________________________
591 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsTrackPt() const
593 if (!fFractionHists[AliPID::kPion])
596 return fFractionHists[AliPID::kPion]->GetNbinsX();
600 //_____________________________________________________________________________
601 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsJetPt() const
603 if (!fFractionHists[AliPID::kPion])
606 return fFractionHists[AliPID::kPion]->GetNbinsY();
610 //_____________________________________________________________________________
611 inline Int_t AliAnalysisTaskPID::GetParticleFractionHistoNbinsCentrality() const
613 if (!fFractionHists[AliPID::kPion])
616 return fFractionHists[AliPID::kPion]->GetNbinsZ();
620 //_____________________________________________________________________________
621 inline Double_t AliAnalysisTaskPID::GetCentralityPercentile(AliVEvent* evt) const
626 Double_t centralityPercentile = -1.;
627 if (fStoreCentralityPercentile)
628 centralityPercentile = evt->GetCentrality()->GetCentralityPercentile(fCentralityEstimator.Data());
630 return centralityPercentile;
634 //_____________________________________________________________________________
635 inline void AliAnalysisTaskPID::PostOutputData()
637 PostData(1, fOutputContainer);
640 PostData(2, fContainerEff);
642 if (fDoPtResolution || fDoDeDxCheck)
643 PostData(3, fQAContainer);