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