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