1 #if !defined( __CINT__) || defined(__MAKECINT__)
11 #include <TParticle.h>
13 #include "AliRunLoader.h"
14 #include "AliLoader.h"
15 #include "AliESDEvent.h"
17 #include "AliESDcascade.h"
18 #include "AliESDMuonTrack.h"
19 #include "AliESDCaloCluster.h"
22 #include "AliHeader.h"
23 #include "AliGenEventHeader.h"
27 TH1F* CreateHisto(const char* name, const char* title,
28 Int_t nBins, Double_t xMin, Double_t xMax,
29 const char* xLabel = NULL, const char* yLabel = NULL)
33 TH1F* result = new TH1F(name, title, nBins, xMin, xMax);
34 result->SetOption("E");
35 if (xLabel) result->GetXaxis()->SetTitle(xLabel);
36 if (yLabel) result->GetYaxis()->SetTitle(yLabel);
37 result->SetMarkerStyle(kFullCircle);
41 TH1F* CreateEffHisto(TH1F* hGen, TH1F* hRec)
43 // create an efficiency histogram
45 Int_t nBins = hGen->GetNbinsX();
46 TH1F* hEff = (TH1F*) hGen->Clone("hEff");
48 hEff->SetStats(kFALSE);
50 hEff->SetMaximum(110.);
51 hEff->GetYaxis()->SetTitle("#epsilon [%]");
53 for (Int_t iBin = 0; iBin <= nBins; iBin++) {
54 Double_t nGen = hGen->GetBinContent(iBin);
55 Double_t nRec = hRec->GetBinContent(iBin);
57 Double_t eff = nRec/nGen;
58 hEff->SetBinContent(iBin, 100. * eff);
59 Double_t error = sqrt(eff*(1.-eff) / nGen);
60 if (error == 0) error = 0.0001;
61 hEff->SetBinError(iBin, 100. * error);
63 hEff->SetBinContent(iBin, -100.);
64 hEff->SetBinError(iBin, 0);
71 Bool_t FitHisto(TH1* histo, Double_t& res, Double_t& resError)
73 // fit a gaussian to a histogram
75 static TF1* fitFunc = new TF1("fitFunc", "gaus");
76 fitFunc->SetLineWidth(2);
77 fitFunc->SetFillStyle(0);
78 Double_t maxFitRange = 2;
80 if (histo->Integral() > 50) {
81 Float_t mean = histo->GetMean();
82 Float_t rms = histo->GetRMS();
83 fitFunc->SetRange(mean - maxFitRange*rms, mean + maxFitRange*rms);
84 fitFunc->SetParameters(mean, rms);
85 histo->Fit(fitFunc, "QRI0");
86 histo->GetFunction("fitFunc")->ResetBit(1<<9);
87 res = TMath::Abs(fitFunc->GetParameter(2));
88 resError = TMath::Abs(fitFunc->GetParError(2));
96 Bool_t CheckESD(const char* gAliceFileName = "galice.root",
97 const char* esdFileName = "AliESDs.root")
99 // check the content of the ESD
102 Int_t checkNGenLow = 1;
104 Double_t checkEffLow = 0.5;
105 Double_t checkEffSigma = 3;
106 Double_t checkFakeHigh = 0.5;
107 Double_t checkFakeSigma = 3;
109 Double_t checkResPtInvHigh = 5;
110 Double_t checkResPtInvSigma = 3;
111 Double_t checkResPhiHigh = 10;
112 Double_t checkResPhiSigma = 3;
113 Double_t checkResThetaHigh = 10;
114 Double_t checkResThetaSigma = 3;
116 Double_t checkPIDEffLow = 0.5;
117 Double_t checkPIDEffSigma = 3;
118 Double_t checkResTOFHigh = 500;
119 Double_t checkResTOFSigma = 3;
121 Double_t checkPHOSNLow = 5;
122 Double_t checkPHOSEnergyLow = 0.3;
123 Double_t checkPHOSEnergyHigh = 1.0;
124 Double_t checkEMCALNLow = 50;
125 Double_t checkEMCALEnergyLow = 0.05;
126 Double_t checkEMCALEnergyHigh = 1.0;
128 Double_t checkMUONNLow = 1;
129 Double_t checkMUONPtLow = 0.5;
130 Double_t checkMUONPtHigh = 10.;
132 Double_t cutPtV0 = 0.3;
133 Double_t checkV0EffLow = 0.02;
134 Double_t checkV0EffSigma = 3;
135 Double_t cutPtCascade = 0.5;
136 Double_t checkCascadeEffLow = 0.01;
137 Double_t checkCascadeEffSigma = 3;
139 // open run loader and load gAlice, kinematics and header
140 AliRunLoader* runLoader = AliRunLoader::Open(gAliceFileName);
142 Error("CheckESD", "getting run loader from file %s failed",
146 runLoader->LoadgAlice();
147 gAlice = runLoader->GetAliRun();
149 Error("CheckESD", "no galice object found");
152 runLoader->LoadKinematics();
153 runLoader->LoadHeader();
156 TFile* esdFile = TFile::Open(esdFileName);
157 if (!esdFile || !esdFile->IsOpen()) {
158 Error("CheckESD", "opening ESD file %s failed", esdFileName);
161 AliESDEvent * esd = new AliESDEvent;
162 TTree* tree = (TTree*) esdFile->Get("esdTree");
164 Error("CheckESD", "no ESD tree found");
167 esd->ReadFromTree(tree);
169 // efficiency and resolution histograms
173 TH1F* hGen = CreateHisto("hGen", "generated tracks",
174 nBinsPt, minPt, maxPt, "p_{t} [GeV/c]", "N");
175 TH1F* hRec = CreateHisto("hRec", "reconstructed tracks",
176 nBinsPt, minPt, maxPt, "p_{t} [GeV/c]", "N");
181 TH1F* hResPtInv = CreateHisto("hResPtInv", "", 100, -10, 10,
182 "(p_{t,rec}^{-1}-p_{t,sim}^{-1}) / p_{t,sim}^{-1} [%]", "N");
183 TH1F* hResPhi = CreateHisto("hResPhi", "", 100, -20, 20,
184 "#phi_{rec}-#phi_{sim} [mrad]", "N");
185 TH1F* hResTheta = CreateHisto("hResTheta", "", 100, -20, 20,
186 "#theta_{rec}-#theta_{sim} [mrad]", "N");
189 Int_t partCode[AliPID::kSPECIES] =
190 {kElectron, kMuonMinus, kPiPlus, kKPlus, kProton};
191 const char* partName[AliPID::kSPECIES+1] =
192 {"electron", "muon", "pion", "kaon", "proton", "other"};
193 Double_t partFrac[AliPID::kSPECIES] =
194 {0.01, 0.01, 0.85, 0.10, 0.05};
195 Int_t identified[AliPID::kSPECIES+1][AliPID::kSPECIES];
196 for (Int_t iGen = 0; iGen < AliPID::kSPECIES+1; iGen++) {
197 for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
198 identified[iGen][iRec] = 0;
201 Int_t nIdentified = 0;
204 TH2F* hDEdxRight = new TH2F("hDEdxRight", "", 300, 0, 3, 100, 0, 400);
205 hDEdxRight->SetStats(kFALSE);
206 hDEdxRight->GetXaxis()->SetTitle("p [GeV/c]");
207 hDEdxRight->GetYaxis()->SetTitle("dE/dx_{TPC}");
208 hDEdxRight->SetMarkerStyle(kFullCircle);
209 hDEdxRight->SetMarkerSize(0.4);
210 TH2F* hDEdxWrong = new TH2F("hDEdxWrong", "", 300, 0, 3, 100, 0, 400);
211 hDEdxWrong->SetStats(kFALSE);
212 hDEdxWrong->GetXaxis()->SetTitle("p [GeV/c]");
213 hDEdxWrong->GetYaxis()->SetTitle("dE/dx_{TPC}");
214 hDEdxWrong->SetMarkerStyle(kFullCircle);
215 hDEdxWrong->SetMarkerSize(0.4);
216 hDEdxWrong->SetMarkerColor(kRed);
217 TH1F* hResTOFRight = CreateHisto("hResTOFRight", "", 100, -1000, 1000,
218 "t_{TOF}-t_{track} [ps]", "N");
219 TH1F* hResTOFWrong = CreateHisto("hResTOFWrong", "", 100, -1000, 1000,
220 "t_{TOF}-t_{track} [ps]", "N");
221 hResTOFWrong->SetLineColor(kRed);
224 TH1F* hEPHOS = CreateHisto("hEPHOS", "PHOS", 100, 0, 50, "E [GeV]", "N");
225 TH1F* hEEMCAL = CreateHisto("hEEMCAL", "EMCAL", 100, 0, 50, "E [GeV]", "N");
228 TH1F* hPtMUON = CreateHisto("hPtMUON", "MUON", 100, 0, 20,
229 "p_{t} [GeV/c]", "N");
232 TH1F* hMassK0 = CreateHisto("hMassK0", "K^{0}", 100, 0.4, 0.6,
233 "M(#pi^{+}#pi^{-}) [GeV/c^{2}]", "N");
234 TH1F* hMassLambda = CreateHisto("hMassLambda", "#Lambda", 100, 1.0, 1.2,
235 "M(p#pi^{-}) [GeV/c^{2}]", "N");
236 TH1F* hMassLambdaBar = CreateHisto("hMassLambdaBar", "#bar{#Lambda}",
238 "M(#bar{p}#pi^{+}) [GeV/c^{2}]", "N");
241 TH1F* hMassXi = CreateHisto("hMassXi", "#Xi", 100, 1.2, 1.5,
242 "M(#Lambda#pi) [GeV/c^{2}]", "N");
243 TH1F* hMassOmega = CreateHisto("hMassOmega", "#Omega", 100, 1.5, 1.8,
244 "M(#LambdaK) [GeV/c^{2}]", "N");
245 Int_t nGenCascades = 0;
246 Int_t nRecCascades = 0;
249 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
250 runLoader->GetEvent(iEvent);
252 // select simulated primary particles, V0s and cascades
253 AliStack* stack = runLoader->Stack();
254 Int_t nParticles = stack->GetNtrack();
256 runLoader->GetHeader()->GenEventHeader()->PrimaryVertex(vertex);
257 TObjArray selParticles;
259 TObjArray selCascades;
260 for (Int_t iParticle = 0; iParticle < nParticles; iParticle++) {
261 TParticle* particle = stack->Particle(iParticle);
262 if (!particle) continue;
263 if (particle->Pt() < 0.001) continue;
264 if (TMath::Abs(particle->Eta()) > 0.9) continue;
265 TVector3 dVertex(particle->Vx() - vertex[0],
266 particle->Vy() - vertex[1],
267 particle->Vz() - vertex[2]);
268 if (dVertex.Mag() > 0.0001) continue;
270 switch (TMath::Abs(particle->GetPdgCode())) {
276 if (particle->Pt() > minPt) {
277 selParticles.Add(particle);
279 hGen->Fill(particle->Pt());
285 if (particle->Pt() > cutPtV0) {
287 selV0s.Add(particle);
293 if (particle->Pt() > cutPtCascade) {
295 selCascades.Add(particle);
303 // get the event summary data
304 tree->GetEvent(iEvent);
306 Error("CheckESD", "no ESD object found for event %d", iEvent);
311 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
312 AliESDtrack* track = esd->GetTrack(iTrack);
314 // select tracks of selected particles
315 Int_t label = TMath::Abs(track->GetLabel());
316 if (label > stack->GetNtrack()) continue; // background
317 TParticle* particle = stack->Particle(label);
318 if (!selParticles.Contains(particle)) continue;
319 if ((track->GetStatus() & AliESDtrack::kITSrefit) == 0) continue;
320 if (track->GetConstrainedChi2() > 1e9) continue;
321 selParticles.Remove(particle); // don't count multiple tracks
324 hRec->Fill(particle->Pt());
325 if (track->GetLabel() < 0) nFake++;
328 hResPtInv->Fill(100. * (TMath::Abs(track->GetSigned1Pt()) - 1./particle->Pt()) *
330 hResPhi->Fill(1000. * (track->Phi() - particle->Phi()));
331 hResTheta->Fill(1000. * (track->Theta() - particle->Theta()));
334 if ((track->GetStatus() & AliESDtrack::kESDpid) == 0) continue;
336 for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
337 if (TMath::Abs(particle->GetPdgCode()) == partCode[i]) iGen = i;
339 Double_t probability[AliPID::kSPECIES];
340 track->GetESDpid(probability);
343 for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
344 probability[i] *= partFrac[i];
345 if (probability[i] > pMax) {
346 pMax = probability[i];
350 identified[iGen][iRec]++;
351 if (iGen == iRec) nIdentified++;
354 Double_t time[AliPID::kSPECIES];
355 track->GetIntegratedTimes(time);
357 hDEdxRight->Fill(particle->P(), track->GetTPCsignal());
358 if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
359 hResTOFRight->Fill(track->GetTOFsignal() - time[iRec]);
362 hDEdxWrong->Fill(particle->P(), track->GetTPCsignal());
363 if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
364 hResTOFWrong->Fill(track->GetTOFsignal() - time[iRec]);
369 // loop over muon tracks
371 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfMuonTracks(); iTrack++) {
372 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
373 Double_t ptInv = TMath::Abs(muonTrack->GetInverseBendingMomentum());
375 hPtMUON->Fill(1./ptInv);
381 for (Int_t iV0 = 0; iV0 < esd->GetNumberOfV0s(); iV0++) {
382 AliESDv0* v0 = esd->GetV0(iV0);
383 if (v0->GetOnFlyStatus()) continue;
384 v0->ChangeMassHypothesis(kK0Short);
385 hMassK0->Fill(v0->GetEffMass());
386 v0->ChangeMassHypothesis(kLambda0);
387 hMassLambda->Fill(v0->GetEffMass());
388 v0->ChangeMassHypothesis(kLambda0Bar);
389 hMassLambdaBar->Fill(v0->GetEffMass());
391 Int_t negLabel = TMath::Abs(esd->GetTrack(v0->GetNindex())->GetLabel());
392 if (negLabel > stack->GetNtrack()) continue; // background
393 Int_t negMother = stack->Particle(negLabel)->GetMother(0);
394 if (negMother < 0) continue;
395 Int_t posLabel = TMath::Abs(esd->GetTrack(v0->GetPindex())->GetLabel());
396 if (posLabel > stack->GetNtrack()) continue; // background
397 Int_t posMother = stack->Particle(posLabel)->GetMother(0);
398 if (negMother != posMother) continue;
399 TParticle* particle = stack->Particle(negMother);
400 if (!selV0s.Contains(particle)) continue;
401 selV0s.Remove(particle);
405 // loop over Cascades
406 for (Int_t iCascade = 0; iCascade < esd->GetNumberOfCascades();
408 AliESDcascade* cascade = esd->GetCascade(iCascade);
410 cascade->ChangeMassHypothesis(v0q,kXiMinus);
411 hMassXi->Fill(cascade->GetEffMassXi());
412 cascade->ChangeMassHypothesis(v0q,kOmegaMinus);
413 hMassOmega->Fill(cascade->GetEffMassXi());
415 Int_t negLabel = TMath::Abs(esd->GetTrack(cascade->GetNindex())
417 if (negLabel > stack->GetNtrack()) continue; // background
418 Int_t negMother = stack->Particle(negLabel)->GetMother(0);
419 if (negMother < 0) continue;
420 Int_t posLabel = TMath::Abs(esd->GetTrack(cascade->GetPindex())
422 if (posLabel > stack->GetNtrack()) continue; // background
423 Int_t posMother = stack->Particle(posLabel)->GetMother(0);
424 if (negMother != posMother) continue;
425 Int_t v0Mother = stack->Particle(negMother)->GetMother(0);
426 if (v0Mother < 0) continue;
427 Int_t bacLabel = TMath::Abs(esd->GetTrack(cascade->GetBindex())
429 if (bacLabel > stack->GetNtrack()) continue; // background
430 Int_t bacMother = stack->Particle(bacLabel)->GetMother(0);
431 if (v0Mother != bacMother) continue;
432 TParticle* particle = stack->Particle(v0Mother);
433 if (!selCascades.Contains(particle)) continue;
434 selCascades.Remove(particle);
438 // loop over the clusters
440 for (Int_t iCluster=0; iCluster<esd->GetNumberOfCaloClusters(); iCluster++) {
441 AliESDCaloCluster * clust = esd->GetCaloCluster(iCluster);
442 if (clust->IsPHOS()) hEPHOS->Fill(clust->E());
443 if (clust->IsEMCAL()) hEEMCAL->Fill(clust->E());
450 if (nGen < checkNGenLow) {
451 Warning("CheckESD", "low number of generated particles: %d", Int_t(nGen));
454 TH1F* hEff = CreateEffHisto(hGen, hRec);
456 Info("CheckESD", "%d out of %d tracks reconstructed including %d "
457 "fake tracks", nRec, nGen, nFake);
460 Double_t eff = nRec*1./nGen;
461 Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGen);
462 Double_t fake = nFake*1./nGen;
463 Double_t fakeError = TMath::Sqrt(fake*(1.-fake) / nGen);
464 Info("CheckESD", "eff = (%.1f +- %.1f) %% fake = (%.1f +- %.1f) %%",
465 100.*eff, 100.*effError, 100.*fake, 100.*fakeError);
467 if (eff < checkEffLow - checkEffSigma*effError) {
468 Warning("CheckESD", "low efficiency: (%.1f +- %.1f) %%",
469 100.*eff, 100.*effError);
471 if (fake > checkFakeHigh + checkFakeSigma*fakeError) {
472 Warning("CheckESD", "high fake: (%.1f +- %.1f) %%",
473 100.*fake, 100.*fakeError);
477 Double_t res, resError;
478 if (FitHisto(hResPtInv, res, resError)) {
479 Info("CheckESD", "relative inverse pt resolution = (%.1f +- %.1f) %%",
481 if (res > checkResPtInvHigh + checkResPtInvSigma*resError) {
482 Warning("CheckESD", "bad pt resolution: (%.1f +- %.1f) %%",
487 if (FitHisto(hResPhi, res, resError)) {
488 Info("CheckESD", "phi resolution = (%.1f +- %.1f) mrad", res, resError);
489 if (res > checkResPhiHigh + checkResPhiSigma*resError) {
490 Warning("CheckESD", "bad phi resolution: (%.1f +- %.1f) mrad",
495 if (FitHisto(hResTheta, res, resError)) {
496 Info("CheckESD", "theta resolution = (%.1f +- %.1f) mrad",
498 if (res > checkResThetaHigh + checkResThetaSigma*resError) {
499 Warning("CheckESD", "bad theta resolution: (%.1f +- %.1f) mrad",
506 Double_t eff = nIdentified*1./nRec;
507 Double_t effError = TMath::Sqrt(eff*(1.-eff) / nRec);
508 Info("CheckESD", "PID eff = (%.1f +- %.1f) %%",
509 100.*eff, 100.*effError);
510 if (eff < checkPIDEffLow - checkPIDEffSigma*effError) {
511 Warning("CheckESD", "low PID efficiency: (%.1f +- %.1f) %%",
512 100.*eff, 100.*effError);
516 printf("%9s:", "gen\\rec");
517 for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
518 printf("%9s", partName[iRec]);
521 for (Int_t iGen = 0; iGen < AliPID::kSPECIES+1; iGen++) {
522 printf("%9s:", partName[iGen]);
523 for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
524 printf("%9d", identified[iGen][iRec]);
529 if (FitHisto(hResTOFRight, res, resError)) {
530 Info("CheckESD", "TOF resolution = (%.1f +- %.1f) ps", res, resError);
531 if (res > checkResTOFHigh + checkResTOFSigma*resError) {
532 Warning("CheckESD", "bad TOF resolution: (%.1f +- %.1f) ps",
538 if (hEPHOS->Integral() < checkPHOSNLow) {
539 Warning("CheckESD", "low number of PHOS particles: %d",
540 Int_t(hEPHOS->Integral()));
542 Double_t mean = hEPHOS->GetMean();
543 if (mean < checkPHOSEnergyLow) {
544 Warning("CheckESD", "low mean PHOS energy: %.1f GeV", mean);
545 } else if (mean > checkPHOSEnergyHigh) {
546 Warning("CheckESD", "high mean PHOS energy: %.1f GeV", mean);
550 if (hEEMCAL->Integral() < checkEMCALNLow) {
551 Warning("CheckESD", "low number of EMCAL particles: %d",
552 Int_t(hEEMCAL->Integral()));
554 Double_t mean = hEEMCAL->GetMean();
555 if (mean < checkEMCALEnergyLow) {
556 Warning("CheckESD", "low mean EMCAL energy: %.1f GeV", mean);
557 } else if (mean > checkEMCALEnergyHigh) {
558 Warning("CheckESD", "high mean EMCAL energy: %.1f GeV", mean);
563 if (hPtMUON->Integral() < checkMUONNLow) {
564 Warning("CheckESD", "low number of MUON particles: %d",
565 Int_t(hPtMUON->Integral()));
567 Double_t mean = hPtMUON->GetMean();
568 if (mean < checkMUONPtLow) {
569 Warning("CheckESD", "low mean MUON pt: %.1f GeV/c", mean);
570 } else if (mean > checkMUONPtHigh) {
571 Warning("CheckESD", "high mean MUON pt: %.1f GeV/c", mean);
577 Double_t eff = nRecV0s*1./nGenV0s;
578 Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGenV0s);
579 if (effError == 0) effError = checkV0EffLow / TMath::Sqrt(1.*nGenV0s);
580 Info("CheckESD", "V0 eff = (%.1f +- %.1f) %%",
581 100.*eff, 100.*effError);
582 if (eff < checkV0EffLow - checkV0EffSigma*effError) {
583 Warning("CheckESD", "low V0 efficiency: (%.1f +- %.1f) %%",
584 100.*eff, 100.*effError);
589 if (nGenCascades > 0) {
590 Double_t eff = nRecCascades*1./nGenCascades;
591 Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGenCascades);
592 if (effError == 0) effError = checkV0EffLow /
593 TMath::Sqrt(1.*nGenCascades);
594 Info("CheckESD", "Cascade eff = (%.1f +- %.1f) %%",
595 100.*eff, 100.*effError);
596 if (eff < checkCascadeEffLow - checkCascadeEffSigma*effError) {
597 Warning("CheckESD", "low Cascade efficiency: (%.1f +- %.1f) %%",
598 100.*eff, 100.*effError);
603 // draw the histograms if not in batch mode
604 if (!gROOT->IsBatch()) {
608 hResPtInv->DrawCopy("E");
610 hResPhi->DrawCopy("E");
612 hResTheta->DrawCopy("E");
614 hDEdxRight->DrawCopy();
615 hDEdxWrong->DrawCopy("SAME");
617 hResTOFRight->DrawCopy("E");
618 hResTOFWrong->DrawCopy("SAME");
620 hEPHOS->DrawCopy("E");
622 hEEMCAL->DrawCopy("E");
624 hPtMUON->DrawCopy("E");
626 hMassK0->DrawCopy("E");
628 hMassLambda->DrawCopy("E");
630 hMassLambdaBar->DrawCopy("E");
632 hMassXi->DrawCopy("E");
634 hMassOmega->DrawCopy("E");
637 // write the output histograms to a file
638 TFile* outputFile = TFile::Open("check.root", "recreate");
639 if (!outputFile || !outputFile->IsOpen()) {
640 Error("CheckESD", "opening output file check.root failed");
649 hResTOFRight->Write();
650 hResTOFWrong->Write();
655 hMassLambda->Write();
656 hMassLambdaBar->Write();
678 delete hMassLambdaBar;
686 runLoader->UnloadHeader();
687 runLoader->UnloadKinematics();
691 Info("CheckESD", "check of ESD was successfull");