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