]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/CheckESD.C
Commit message
[u/mrichter/AliRoot.git] / STEER / CheckESD.C
1 #if !defined( __CINT__) || defined(__MAKECINT__)
2 #include <TROOT.h>
3 #include <TFile.h>
4 #include <TError.h>
5 #include <TH1.h>
6 #include <TH2.h>
7 #include <TF1.h>
8 #include <TCanvas.h>
9 #include <TVector3.h>
10 #include <TPDGCode.h>
11 #include <TParticle.h>
12
13 #include "AliRunLoader.h"
14 #include "AliLoader.h"
15 #include "AliESDEvent.h"
16 #include "AliESDv0.h"
17 #include "AliESDcascade.h"
18 #include "AliESDMuonTrack.h"
19 #include "AliESDCaloCluster.h"
20 #include "AliRun.h"
21 #include "AliStack.h"
22 #include "AliHeader.h"
23 #include "AliGenEventHeader.h"
24 #include "AliPID.h"
25 #include "AliESDpid.h"
26 #endif
27
28 TH1F* CreateHisto(const char* name, const char* title, 
29                   Int_t nBins, Double_t xMin, Double_t xMax,
30                   const char* xLabel = NULL, const char* yLabel = NULL)
31 {
32 // create a histogram
33
34   TH1F* result = new TH1F(name, title, nBins, xMin, xMax);
35   result->SetOption("E");
36   if (xLabel) result->GetXaxis()->SetTitle(xLabel);
37   if (yLabel) result->GetYaxis()->SetTitle(yLabel);
38   result->SetMarkerStyle(kFullCircle);
39   return result;
40 }
41
42 TH1F* CreateEffHisto(TH1F* hGen, TH1F* hRec)
43 {
44 // create an efficiency histogram
45
46   Int_t nBins = hGen->GetNbinsX();
47   TH1F* hEff = (TH1F*) hGen->Clone("hEff");
48   hEff->SetTitle("");
49   hEff->SetStats(kFALSE);
50   hEff->SetMinimum(0.);
51   hEff->SetMaximum(110.);
52   hEff->GetYaxis()->SetTitle("#epsilon [%]");
53
54   for (Int_t iBin = 0; iBin <= nBins; iBin++) {
55     Double_t nGen = hGen->GetBinContent(iBin);
56     Double_t nRec = hRec->GetBinContent(iBin);
57     if (nGen > 0) {
58       Double_t eff = nRec/nGen;
59       hEff->SetBinContent(iBin, 100. * eff);
60       Double_t error = sqrt(eff*(1.-eff) / nGen);
61       if (error == 0) error = 0.0001;
62       hEff->SetBinError(iBin, 100. * error);
63     } else {
64       hEff->SetBinContent(iBin, -100.);
65       hEff->SetBinError(iBin, 0);
66     }
67   }
68
69   return hEff;
70 }
71
72 Bool_t FitHisto(TH1* histo, Double_t& res, Double_t& resError)
73 {
74 // fit a gaussian to a histogram
75
76   static TF1* fitFunc = new TF1("fitFunc", "gaus");
77   fitFunc->SetLineWidth(2);
78   fitFunc->SetFillStyle(0);
79   Double_t maxFitRange = 2;
80
81   if (histo->Integral() > 50) {
82     Float_t mean = histo->GetMean();
83     Float_t rms = histo->GetRMS();
84     fitFunc->SetRange(mean - maxFitRange*rms, mean + maxFitRange*rms);
85     fitFunc->SetParameters(mean, rms);
86     histo->Fit(fitFunc, "QRI0");
87     histo->GetFunction("fitFunc")->ResetBit(1<<9);
88     res = TMath::Abs(fitFunc->GetParameter(2));
89     resError = TMath::Abs(fitFunc->GetParError(2));
90     return kTRUE;
91   }
92
93   return kFALSE;
94 }
95
96
97 Bool_t CheckESD(const char* gAliceFileName = "galice.root", 
98                 const char* esdFileName = "AliESDs.root")
99 {
100 // check the content of the ESD
101  
102   // check values
103   Int_t    checkNGenLow = 1;
104
105   Double_t checkEffLow = 0.5;
106   Double_t checkEffSigma = 3;
107   Double_t checkFakeHigh = 0.5;
108   Double_t checkFakeSigma = 3;
109
110   Double_t checkResPtInvHigh = 5;
111   Double_t checkResPtInvSigma = 3;
112   Double_t checkResPhiHigh = 10;
113   Double_t checkResPhiSigma = 3;
114   Double_t checkResThetaHigh = 10;
115   Double_t checkResThetaSigma = 3;
116
117   Double_t checkPIDEffLow = 0.5;
118   Double_t checkPIDEffSigma = 3;
119   Double_t checkResTOFHigh = 500;
120   Double_t checkResTOFSigma = 3;
121
122   Double_t checkPHOSNLow = 5;
123   Double_t checkPHOSEnergyLow = 0.3;
124   Double_t checkPHOSEnergyHigh = 1.0;
125   Double_t checkEMCALNLow = 50;
126   Double_t checkEMCALEnergyLow = 0.05;
127   Double_t checkEMCALEnergyHigh = 1.0;
128
129   Double_t checkMUONNLow = 1;
130   Double_t checkMUONPtLow = 0.5;
131   Double_t checkMUONPtHigh = 10.;
132
133   Double_t cutPtV0 = 0.3;
134   Double_t checkV0EffLow = 0.02;
135   Double_t checkV0EffSigma = 3;
136   Double_t cutPtCascade = 0.5;
137   Double_t checkCascadeEffLow = 0.01;
138   Double_t checkCascadeEffSigma = 3;
139
140   // open run loader and load gAlice, kinematics and header
141   AliRunLoader* runLoader = AliRunLoader::Open(gAliceFileName);
142   if (!runLoader) {
143     Error("CheckESD", "getting run loader from file %s failed", 
144             gAliceFileName);
145     return kFALSE;
146   }
147   runLoader->LoadgAlice();
148   gAlice = runLoader->GetAliRun();
149   if (!gAlice) {
150     Error("CheckESD", "no galice object found");
151     return kFALSE;
152   }
153   runLoader->LoadKinematics();
154   runLoader->LoadHeader();
155
156   // open the ESD file
157   TFile* esdFile = TFile::Open(esdFileName);
158   if (!esdFile || !esdFile->IsOpen()) {
159     Error("CheckESD", "opening ESD file %s failed", esdFileName);
160     return kFALSE;
161   }
162   AliESDEvent * esd = new AliESDEvent;
163   TTree* tree = (TTree*) esdFile->Get("esdTree");
164   if (!tree) {
165     Error("CheckESD", "no ESD tree found");
166     return kFALSE;
167   }
168   esd->ReadFromTree(tree);
169
170   // PID
171
172   AliESDpid * pid = new AliESDpid(kTRUE);
173
174   // efficiency and resolution histograms
175   Int_t nBinsPt = 15;
176   Float_t minPt = 0.1;
177   Float_t maxPt = 3.1;
178   TH1F* hGen = CreateHisto("hGen", "generated tracks", 
179                            nBinsPt, minPt, maxPt, "p_{t} [GeV/c]", "N");
180   TH1F* hRec = CreateHisto("hRec", "reconstructed tracks", 
181                            nBinsPt, minPt, maxPt, "p_{t} [GeV/c]", "N");
182   Int_t nGen = 0;
183   Int_t nRec = 0;
184   Int_t nFake = 0;
185
186   TH1F* hResPtInv = CreateHisto("hResPtInv", "", 100, -10, 10, 
187            "(p_{t,rec}^{-1}-p_{t,sim}^{-1}) / p_{t,sim}^{-1} [%]", "N");
188   TH1F* hResPhi = CreateHisto("hResPhi", "", 100, -20, 20, 
189                               "#phi_{rec}-#phi_{sim} [mrad]", "N");
190   TH1F* hResTheta = CreateHisto("hResTheta", "", 100, -20, 20, 
191                                 "#theta_{rec}-#theta_{sim} [mrad]", "N");
192
193   // PID
194   Int_t partCode[AliPID::kSPECIES] = 
195     {kElectron, kMuonMinus, kPiPlus, kKPlus, kProton};
196   const char* partName[AliPID::kSPECIES+1] = 
197     {"electron", "muon", "pion", "kaon", "proton", "other"};
198   Double_t partFrac[AliPID::kSPECIES] = 
199     {0.01, 0.01, 0.85, 0.10, 0.05};
200   Int_t identified[AliPID::kSPECIES+1][AliPID::kSPECIES];
201   for (Int_t iGen = 0; iGen < AliPID::kSPECIES+1; iGen++) {
202     for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
203       identified[iGen][iRec] = 0;
204     }
205   }
206   Int_t nIdentified = 0;
207
208   // dE/dx and TOF
209   TH2F* hDEdxRight = new TH2F("hDEdxRight", "", 300, 0, 3, 100, 0, 400);
210   hDEdxRight->SetStats(kFALSE);
211   hDEdxRight->GetXaxis()->SetTitle("p [GeV/c]");
212   hDEdxRight->GetYaxis()->SetTitle("dE/dx_{TPC}");
213   hDEdxRight->SetMarkerStyle(kFullCircle);
214   hDEdxRight->SetMarkerSize(0.4);
215   TH2F* hDEdxWrong = new TH2F("hDEdxWrong", "", 300, 0, 3, 100, 0, 400);
216   hDEdxWrong->SetStats(kFALSE);
217   hDEdxWrong->GetXaxis()->SetTitle("p [GeV/c]");
218   hDEdxWrong->GetYaxis()->SetTitle("dE/dx_{TPC}");
219   hDEdxWrong->SetMarkerStyle(kFullCircle);
220   hDEdxWrong->SetMarkerSize(0.4);
221   hDEdxWrong->SetMarkerColor(kRed);
222   TH1F* hResTOFRight = CreateHisto("hResTOFRight", "", 100, -1000, 1000, 
223                                    "t_{TOF}-t_{track} [ps]", "N");
224   TH1F* hResTOFWrong = CreateHisto("hResTOFWrong", "", 100, -1000, 1000, 
225                                    "t_{TOF}-t_{track} [ps]", "N");
226   hResTOFWrong->SetLineColor(kRed);
227
228   // calorimeters
229   TH1F* hEPHOS = CreateHisto("hEPHOS", "PHOS", 100, 0, 50, "E [GeV]", "N");
230   TH1F* hEEMCAL = CreateHisto("hEEMCAL", "EMCAL", 100, 0, 50, "E [GeV]", "N");
231
232   // muons
233   TH1F* hPtMUON = CreateHisto("hPtMUON", "MUON", 100, 0, 20, 
234                               "p_{t} [GeV/c]", "N");
235
236   // V0s and cascades
237   TH1F* hMassK0 = CreateHisto("hMassK0", "K^{0}", 100, 0.4, 0.6, 
238                               "M(#pi^{+}#pi^{-}) [GeV/c^{2}]", "N");
239   TH1F* hMassLambda = CreateHisto("hMassLambda", "#Lambda", 100, 1.0, 1.2, 
240                                   "M(p#pi^{-}) [GeV/c^{2}]", "N");
241   TH1F* hMassLambdaBar = CreateHisto("hMassLambdaBar", "#bar{#Lambda}", 
242                                      100, 1.0, 1.2, 
243                                      "M(#bar{p}#pi^{+}) [GeV/c^{2}]", "N");
244   Int_t nGenV0s = 0;
245   Int_t nRecV0s = 0;
246   TH1F* hMassXi = CreateHisto("hMassXi", "#Xi", 100, 1.2, 1.5, 
247                               "M(#Lambda#pi) [GeV/c^{2}]", "N");
248   TH1F* hMassOmega = CreateHisto("hMassOmega", "#Omega", 100, 1.5, 1.8, 
249                                  "M(#LambdaK) [GeV/c^{2}]", "N");
250   Int_t nGenCascades = 0;
251   Int_t nRecCascades = 0;
252
253   // loop over events
254   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
255     runLoader->GetEvent(iEvent);
256
257     // select simulated primary particles, V0s and cascades
258     AliStack* stack = runLoader->Stack();
259     Int_t nParticles = stack->GetNtrack();
260     TArrayF vertex(3);
261     runLoader->GetHeader()->GenEventHeader()->PrimaryVertex(vertex);
262     TObjArray selParticles;
263     TObjArray selV0s;
264     TObjArray selCascades;
265     for (Int_t iParticle = 0; iParticle < nParticles; iParticle++) {
266       TParticle* particle = stack->Particle(iParticle);
267       if (!particle) continue;
268       if (particle->Pt() < 0.001) continue;
269       if (TMath::Abs(particle->Eta()) > 0.9) continue;
270       TVector3 dVertex(particle->Vx() - vertex[0], 
271                        particle->Vy() - vertex[1],
272                        particle->Vz() - vertex[2]);
273       if (dVertex.Mag() > 0.0001) continue;
274
275       switch (TMath::Abs(particle->GetPdgCode())) {
276       case kElectron:
277       case kMuonMinus:
278       case kPiPlus:
279       case kKPlus:
280       case kProton: {
281         if (particle->Pt() > minPt) {
282           selParticles.Add(particle);
283           nGen++;
284           hGen->Fill(particle->Pt());
285         }
286         break;
287       }
288       case kK0Short:
289       case kLambda0: {
290         if (particle->Pt() > cutPtV0) {
291           nGenV0s++;
292           selV0s.Add(particle);
293         }
294         break;
295       }
296       case kXiMinus:
297       case kOmegaMinus: {
298         if (particle->Pt() > cutPtCascade) {
299           nGenCascades++;
300           selCascades.Add(particle);
301         }
302         break;
303       }
304       default: break;
305       }
306     }
307
308     // get the event summary data
309     tree->GetEvent(iEvent);
310     if (!esd) {
311       Error("CheckESD", "no ESD object found for event %d", iEvent);
312       return kFALSE;
313     }
314
315     // PID for MC
316     pid->MakePID(esd,kTRUE);
317
318     // loop over tracks
319     for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
320       AliESDtrack* track = esd->GetTrack(iTrack);
321
322       // select tracks of selected particles
323       Int_t label = TMath::Abs(track->GetLabel());
324       if (label > stack->GetNtrack()) continue;     // background
325       TParticle* particle = stack->Particle(label);
326       if (!selParticles.Contains(particle)) continue;
327       if ((track->GetStatus() & AliESDtrack::kITSrefit) == 0) continue;
328       if (track->GetConstrainedChi2() > 1e9) continue;
329       selParticles.Remove(particle);   // don't count multiple tracks
330
331       nRec++;
332       hRec->Fill(particle->Pt());
333       if (track->GetLabel() < 0) nFake++;
334
335       // resolutions
336       hResPtInv->Fill(100. * (TMath::Abs(track->GetSigned1Pt()) - 1./particle->Pt()) * 
337                       particle->Pt());
338       hResPhi->Fill(1000. * (track->Phi() - particle->Phi()));
339       hResTheta->Fill(1000. * (track->Theta() - particle->Theta()));
340
341       // PID
342       if ((track->GetStatus() & AliESDtrack::kESDpid) == 0) continue;
343       Int_t iGen = 5;
344       for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
345         if (TMath::Abs(particle->GetPdgCode()) == partCode[i]) iGen = i;
346       }
347       Double_t probability[AliPID::kSPECIES];
348       track->GetESDpid(probability);
349       Double_t pMax = 0;
350       Int_t iRec = 0;
351       for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
352         probability[i] *= partFrac[i];
353         if (probability[i] > pMax) {
354           pMax = probability[i];
355           iRec = i;
356         }
357       }
358       identified[iGen][iRec]++;
359       if (iGen == iRec) nIdentified++;
360
361       // dE/dx and TOF
362       Double_t time[AliPID::kSPECIESC];
363       track->GetIntegratedTimes(time);
364       if (iGen == iRec) {
365         hDEdxRight->Fill(particle->P(), track->GetTPCsignal());
366         if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
367           hResTOFRight->Fill(track->GetTOFsignal() - time[iRec]);
368         }
369       } else {
370         hDEdxWrong->Fill(particle->P(), track->GetTPCsignal());
371         if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
372           hResTOFWrong->Fill(track->GetTOFsignal() - time[iRec]);
373         }
374       }
375     }
376
377     // loop over muon tracks
378     {
379     for (Int_t iTrack = 0; iTrack < esd->GetNumberOfMuonTracks(); iTrack++) {
380       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
381       Double_t ptInv = TMath::Abs(muonTrack->GetInverseBendingMomentum());
382       if (ptInv > 0.001) {
383         hPtMUON->Fill(1./ptInv);
384       }
385     }
386     }
387
388     // loop over V0s
389     for (Int_t iV0 = 0; iV0 < esd->GetNumberOfV0s(); iV0++) {
390       AliESDv0* v0 = esd->GetV0(iV0);
391       if (v0->GetOnFlyStatus()) continue;
392       v0->ChangeMassHypothesis(kK0Short);
393       hMassK0->Fill(v0->GetEffMass());
394       v0->ChangeMassHypothesis(kLambda0);
395       hMassLambda->Fill(v0->GetEffMass());
396       v0->ChangeMassHypothesis(kLambda0Bar);
397       hMassLambdaBar->Fill(v0->GetEffMass());
398
399       Int_t negLabel = TMath::Abs(esd->GetTrack(v0->GetNindex())->GetLabel());
400       if (negLabel > stack->GetNtrack()) continue;     // background
401       Int_t negMother = stack->Particle(negLabel)->GetMother(0);
402       if (negMother < 0) continue;
403       Int_t posLabel = TMath::Abs(esd->GetTrack(v0->GetPindex())->GetLabel());
404       if (posLabel > stack->GetNtrack()) continue;     // background
405       Int_t posMother = stack->Particle(posLabel)->GetMother(0);
406       if (negMother != posMother) continue;
407       TParticle* particle = stack->Particle(negMother);
408       if (!selV0s.Contains(particle)) continue;
409       selV0s.Remove(particle);
410       nRecV0s++;
411     }
412
413     // loop over Cascades
414     for (Int_t iCascade = 0; iCascade < esd->GetNumberOfCascades(); 
415          iCascade++) {
416       AliESDcascade* cascade = esd->GetCascade(iCascade);
417       Double_t v0q;
418       cascade->ChangeMassHypothesis(v0q,kXiMinus);
419       hMassXi->Fill(cascade->GetEffMassXi());
420       cascade->ChangeMassHypothesis(v0q,kOmegaMinus);
421       hMassOmega->Fill(cascade->GetEffMassXi());
422
423       Int_t negLabel = TMath::Abs(esd->GetTrack(cascade->GetNindex())
424                                   ->GetLabel());
425       if (negLabel > stack->GetNtrack()) continue;     // background
426       Int_t negMother = stack->Particle(negLabel)->GetMother(0);
427       if (negMother < 0) continue;
428       Int_t posLabel = TMath::Abs(esd->GetTrack(cascade->GetPindex())
429                                   ->GetLabel());
430       if (posLabel > stack->GetNtrack()) continue;     // background
431       Int_t posMother = stack->Particle(posLabel)->GetMother(0);
432       if (negMother != posMother) continue;
433       Int_t v0Mother = stack->Particle(negMother)->GetMother(0);
434       if (v0Mother < 0) continue;
435       Int_t bacLabel = TMath::Abs(esd->GetTrack(cascade->GetBindex())
436                                   ->GetLabel());
437       if (bacLabel > stack->GetNtrack()) continue;     // background
438       Int_t bacMother = stack->Particle(bacLabel)->GetMother(0);
439       if (v0Mother != bacMother) continue;
440       TParticle* particle = stack->Particle(v0Mother);
441       if (!selCascades.Contains(particle)) continue;
442       selCascades.Remove(particle);
443       nRecCascades++;
444     }
445
446     // loop over the clusters
447     {
448       for (Int_t iCluster=0; iCluster<esd->GetNumberOfCaloClusters(); iCluster++) {
449         AliESDCaloCluster * clust = esd->GetCaloCluster(iCluster);
450         if (clust->IsPHOS()) hEPHOS->Fill(clust->E());
451         if (clust->IsEMCAL()) hEEMCAL->Fill(clust->E());
452       }
453     }
454
455   }
456
457   // perform checks
458   if (nGen < checkNGenLow) {
459     Warning("CheckESD", "low number of generated particles: %d", Int_t(nGen));
460   }
461
462   TH1F* hEff = CreateEffHisto(hGen, hRec);
463
464   Info("CheckESD", "%d out of %d tracks reconstructed including %d "
465          "fake tracks", nRec, nGen, nFake);
466   if (nGen > 0) {
467     // efficiency
468     Double_t eff = nRec*1./nGen;
469     Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGen);
470     Double_t fake = nFake*1./nGen;
471     Double_t fakeError = TMath::Sqrt(fake*(1.-fake) / nGen);
472     Info("CheckESD", "eff = (%.1f +- %.1f) %%  fake = (%.1f +- %.1f) %%",
473          100.*eff, 100.*effError, 100.*fake, 100.*fakeError);
474
475     if (eff < checkEffLow - checkEffSigma*effError) {
476       Warning("CheckESD", "low efficiency: (%.1f +- %.1f) %%", 
477               100.*eff, 100.*effError);
478     }
479     if (fake > checkFakeHigh + checkFakeSigma*fakeError) {
480       Warning("CheckESD", "high fake: (%.1f +- %.1f) %%", 
481               100.*fake, 100.*fakeError);
482     }
483
484     // resolutions
485     Double_t res, resError;
486     if (FitHisto(hResPtInv, res, resError)) {
487       Info("CheckESD", "relative inverse pt resolution = (%.1f +- %.1f) %%",
488            res, resError);
489       if (res > checkResPtInvHigh + checkResPtInvSigma*resError) {
490         Warning("CheckESD", "bad pt resolution: (%.1f +- %.1f) %%", 
491                 res, resError);
492       }
493     }
494
495     if (FitHisto(hResPhi, res, resError)) {
496       Info("CheckESD", "phi resolution = (%.1f +- %.1f) mrad", res, resError);
497       if (res > checkResPhiHigh + checkResPhiSigma*resError) {
498         Warning("CheckESD", "bad phi resolution: (%.1f +- %.1f) mrad", 
499                 res, resError);
500       }
501     }
502
503     if (FitHisto(hResTheta, res, resError)) {
504       Info("CheckESD", "theta resolution = (%.1f +- %.1f) mrad", 
505            res, resError);
506       if (res > checkResThetaHigh + checkResThetaSigma*resError) {
507         Warning("CheckESD", "bad theta resolution: (%.1f +- %.1f) mrad", 
508                 res, resError);
509       }
510     }
511
512     // PID
513     if (nRec > 0) {
514       Double_t eff = nIdentified*1./nRec;
515       Double_t effError = TMath::Sqrt(eff*(1.-eff) / nRec);
516       Info("CheckESD", "PID eff = (%.1f +- %.1f) %%", 
517            100.*eff, 100.*effError);
518       if (eff < checkPIDEffLow - checkPIDEffSigma*effError) {
519         Warning("CheckESD", "low PID efficiency: (%.1f +- %.1f) %%", 
520                 100.*eff, 100.*effError);
521       }
522     }
523
524     printf("%9s:", "gen\\rec");
525     for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
526       printf("%9s", partName[iRec]);
527     }
528     printf("\n");
529     for (Int_t iGen = 0; iGen < AliPID::kSPECIES+1; iGen++) {
530       printf("%9s:", partName[iGen]);
531       for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
532         printf("%9d", identified[iGen][iRec]);
533       }
534       printf("\n");
535     }
536
537     if (FitHisto(hResTOFRight, res, resError)) {
538       Info("CheckESD", "TOF resolution = (%.1f +- %.1f) ps", res, resError);
539       if (res > checkResTOFHigh + checkResTOFSigma*resError) {
540         Warning("CheckESD", "bad TOF resolution: (%.1f +- %.1f) ps", 
541                 res, resError);
542       }
543     }
544
545     // calorimeters
546     if (hEPHOS->Integral() < checkPHOSNLow) {
547       Warning("CheckESD", "low number of PHOS particles: %d", 
548               Int_t(hEPHOS->Integral()));
549     } else {
550       Double_t mean = hEPHOS->GetMean();
551       if (mean < checkPHOSEnergyLow) {
552         Warning("CheckESD", "low mean PHOS energy: %.1f GeV", mean);
553       } else if (mean > checkPHOSEnergyHigh) {
554         Warning("CheckESD", "high mean PHOS energy: %.1f GeV", mean);
555       }
556     }
557
558     if (hEEMCAL->Integral() < checkEMCALNLow) {
559       Warning("CheckESD", "low number of EMCAL particles: %d", 
560               Int_t(hEEMCAL->Integral()));
561     } else {
562       Double_t mean = hEEMCAL->GetMean();
563       if (mean < checkEMCALEnergyLow) {
564         Warning("CheckESD", "low mean EMCAL energy: %.1f GeV", mean);
565       } else if (mean > checkEMCALEnergyHigh) {
566         Warning("CheckESD", "high mean EMCAL energy: %.1f GeV", mean);
567       }
568     }
569
570     // muons
571     if (hPtMUON->Integral() < checkMUONNLow) {
572       Warning("CheckESD", "low number of MUON particles: %d", 
573               Int_t(hPtMUON->Integral()));
574     } else {
575       Double_t mean = hPtMUON->GetMean();
576       if (mean < checkMUONPtLow) {
577         Warning("CheckESD", "low mean MUON pt: %.1f GeV/c", mean);
578       } else if (mean > checkMUONPtHigh) {
579         Warning("CheckESD", "high mean MUON pt: %.1f GeV/c", mean);
580       }
581     }
582
583     // V0s
584     if (nGenV0s > 0) {
585       Double_t eff = nRecV0s*1./nGenV0s;
586       Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGenV0s);
587       if (effError == 0) effError = checkV0EffLow / TMath::Sqrt(1.*nGenV0s);
588       Info("CheckESD", "V0 eff = (%.1f +- %.1f) %%", 
589            100.*eff, 100.*effError);
590       if (eff < checkV0EffLow - checkV0EffSigma*effError) {
591         Warning("CheckESD", "low V0 efficiency: (%.1f +- %.1f) %%", 
592                 100.*eff, 100.*effError);
593       }
594     }
595
596     // Cascades
597     if (nGenCascades > 0) {
598       Double_t eff = nRecCascades*1./nGenCascades;
599       Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGenCascades);
600       if (effError == 0) effError = checkV0EffLow / 
601                            TMath::Sqrt(1.*nGenCascades);
602       Info("CheckESD", "Cascade eff = (%.1f +- %.1f) %%", 
603            100.*eff, 100.*effError);
604       if (eff < checkCascadeEffLow - checkCascadeEffSigma*effError) {
605         Warning("CheckESD", "low Cascade efficiency: (%.1f +- %.1f) %%", 
606                 100.*eff, 100.*effError);
607       }
608     }
609   }
610
611   // draw the histograms if not in batch mode
612   if (!gROOT->IsBatch()) {
613     new TCanvas;
614     hEff->DrawCopy();
615     new TCanvas;
616     hResPtInv->DrawCopy("E");
617     new TCanvas;
618     hResPhi->DrawCopy("E");
619     new TCanvas;
620     hResTheta->DrawCopy("E");
621     new TCanvas;
622     hDEdxRight->DrawCopy();
623     hDEdxWrong->DrawCopy("SAME");
624     new TCanvas;
625     hResTOFRight->DrawCopy("E");
626     hResTOFWrong->DrawCopy("SAME");
627     new TCanvas;
628     hEPHOS->DrawCopy("E");
629     new TCanvas;
630     hEEMCAL->DrawCopy("E");
631     new TCanvas;
632     hPtMUON->DrawCopy("E");
633     new TCanvas;
634     hMassK0->DrawCopy("E");
635     new TCanvas;
636     hMassLambda->DrawCopy("E");
637     new TCanvas;
638     hMassLambdaBar->DrawCopy("E");
639     new TCanvas;
640     hMassXi->DrawCopy("E");
641     new TCanvas;
642     hMassOmega->DrawCopy("E");
643   }
644
645   // write the output histograms to a file
646   TFile* outputFile = TFile::Open("check.root", "recreate");
647   if (!outputFile || !outputFile->IsOpen()) {
648     Error("CheckESD", "opening output file check.root failed");
649     return kFALSE;
650   }
651   hEff->Write();
652   hResPtInv->Write();
653   hResPhi->Write();
654   hResTheta->Write();
655   hDEdxRight->Write();
656   hDEdxWrong->Write();
657   hResTOFRight->Write();
658   hResTOFWrong->Write();
659   hEPHOS->Write();
660   hEEMCAL->Write();
661   hPtMUON->Write();
662   hMassK0->Write();
663   hMassLambda->Write();
664   hMassLambdaBar->Write();
665   hMassXi->Write();
666   hMassOmega->Write();
667   outputFile->Close();
668   delete outputFile;
669
670   // clean up
671   delete hGen;
672   delete hRec;
673   delete hEff;
674   delete hResPtInv;
675   delete hResPhi;
676   delete hResTheta;
677   delete hDEdxRight;
678   delete hDEdxWrong;
679   delete hResTOFRight;
680   delete hResTOFWrong;
681   delete hEPHOS;
682   delete hEEMCAL;
683   delete hPtMUON;
684   delete hMassK0;
685   delete hMassLambda;
686   delete hMassLambdaBar;
687   delete hMassXi;
688   delete hMassOmega;
689
690   delete esd;
691   esdFile->Close();
692   delete esdFile;
693
694   runLoader->UnloadHeader();
695   runLoader->UnloadKinematics();
696   delete runLoader;
697   delete pid;
698
699   // result of check
700   Info("CheckESD", "check of ESD was successfull");
701   return kTRUE;
702 }