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