]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/CheckESD.C
Possibility to compile the macro
[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, 5, "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       Double_t p[3];
332       track->GetConstrainedPxPyPz(p);
333       TVector3 pTrack(p);
334       hResPtInv->Fill(100. * (1./pTrack.Pt() - 1./particle->Pt()) * 
335                       particle->Pt());
336       hResPhi->Fill(1000. * (pTrack.Phi() - particle->Phi()));
337       hResTheta->Fill(1000. * (pTrack.Theta() - particle->Theta()));
338
339       // PID
340       if ((track->GetStatus() & AliESDtrack::kESDpid) == 0) continue;
341       Int_t iGen = 5;
342       for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
343         if (TMath::Abs(particle->GetPdgCode()) == partCode[i]) iGen = i;
344       }
345       Double_t probability[AliPID::kSPECIES];
346       track->GetESDpid(probability);
347       Double_t pMax = 0;
348       Int_t iRec = 0;
349       for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
350         probability[i] *= partFrac[i];
351         if (probability[i] > pMax) {
352           pMax = probability[i];
353           iRec = i;
354         }
355       }
356       identified[iGen][iRec]++;
357       if (iGen == iRec) nIdentified++;
358
359       // dE/dx and TOF
360       Double_t time[AliPID::kSPECIES];
361       track->GetIntegratedTimes(time);
362       if (iGen == iRec) {
363         hDEdxRight->Fill(pTrack.Mag(), track->GetTPCsignal());
364         if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
365           hResTOFRight->Fill(track->GetTOFsignal() - time[iRec]);
366         }
367       } else {
368         hDEdxWrong->Fill(pTrack.Mag(), track->GetTPCsignal());
369         if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
370           hResTOFWrong->Fill(track->GetTOFsignal() - time[iRec]);
371         }
372       }
373     }
374
375     // loop over muon tracks
376     {
377     for (Int_t iTrack = 0; iTrack < esd->GetNumberOfMuonTracks(); iTrack++) {
378       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
379       Double_t ptInv = TMath::Abs(muonTrack->GetInverseBendingMomentum());
380       if (ptInv > 0.001) {
381         hPtMUON->Fill(1./ptInv);
382       }
383     }
384     }
385
386     // loop over V0s
387     for (Int_t iV0 = 0; iV0 < esd->GetNumberOfV0s(); iV0++) {
388       AliESDv0* v0 = esd->GetV0(iV0);
389       if (v0->GetOnFlyStatus()) continue;
390       v0->ChangeMassHypothesis(kK0Short);
391       hMassK0->Fill(v0->GetEffMass());
392       v0->ChangeMassHypothesis(kLambda0);
393       hMassLambda->Fill(v0->GetEffMass());
394       v0->ChangeMassHypothesis(kLambda0Bar);
395       hMassLambdaBar->Fill(v0->GetEffMass());
396
397       Int_t negLabel = TMath::Abs(esd->GetTrack(v0->GetNindex())->GetLabel());
398       if (negLabel > stack->GetNtrack()) continue;     // background
399       Int_t negMother = stack->Particle(negLabel)->GetMother(0);
400       if (negMother < 0) continue;
401       Int_t posLabel = TMath::Abs(esd->GetTrack(v0->GetPindex())->GetLabel());
402       if (posLabel > stack->GetNtrack()) continue;     // background
403       Int_t posMother = stack->Particle(posLabel)->GetMother(0);
404       if (negMother != posMother) continue;
405       TParticle* particle = stack->Particle(negMother);
406       if (!selV0s.Contains(particle)) continue;
407       selV0s.Remove(particle);
408       nRecV0s++;
409     }
410
411     // loop over Cascades
412     for (Int_t iCascade = 0; iCascade < esd->GetNumberOfCascades(); 
413          iCascade++) {
414       AliESDcascade* cascade = esd->GetCascade(iCascade);
415       Double_t v0q;
416       cascade->ChangeMassHypothesis(v0q,kXiMinus);
417       hMassXi->Fill(cascade->GetEffMass());
418       cascade->ChangeMassHypothesis(v0q,kOmegaMinus);
419       hMassOmega->Fill(cascade->GetEffMass());
420
421       Int_t negLabel = TMath::Abs(esd->GetTrack(cascade->GetNindex())
422                                   ->GetLabel());
423       if (negLabel > stack->GetNtrack()) continue;     // background
424       Int_t negMother = stack->Particle(negLabel)->GetMother(0);
425       if (negMother < 0) continue;
426       Int_t posLabel = TMath::Abs(esd->GetTrack(cascade->GetPindex())
427                                   ->GetLabel());
428       if (posLabel > stack->GetNtrack()) continue;     // background
429       Int_t posMother = stack->Particle(posLabel)->GetMother(0);
430       if (negMother != posMother) continue;
431       Int_t v0Mother = stack->Particle(negMother)->GetMother(0);
432       if (v0Mother < 0) continue;
433       Int_t bacLabel = TMath::Abs(esd->GetTrack(cascade->GetBindex())
434                                   ->GetLabel());
435       if (bacLabel > stack->GetNtrack()) continue;     // background
436       Int_t bacMother = stack->Particle(bacLabel)->GetMother(0);
437       if (v0Mother != bacMother) continue;
438       TParticle* particle = stack->Particle(v0Mother);
439       if (!selCascades.Contains(particle)) continue;
440       selCascades.Remove(particle);
441       nRecCascades++;
442     }
443
444     // loop over the PHOS clusters
445     {
446     Int_t firstPHOSCluster = esd->GetFirstPHOSCluster();
447     Int_t lastPHOSCluster  = firstPHOSCluster + esd->GetNumberOfPHOSClusters();
448     for (Int_t iCluster=firstPHOSCluster; iCluster<lastPHOSCluster; iCluster++)
449       hEPHOS->Fill(esd->GetCaloCluster(iCluster)->E());
450     }
451
452     // loop over the EMCAL clusters
453     {
454     Int_t firstEMCALCluster = esd->GetFirstEMCALCluster();
455     Int_t lastEMCALCluster  = firstEMCALCluster + esd->GetNumberOfEMCALClusters();
456     for (Int_t iCluster=firstEMCALCluster; iCluster<lastEMCALCluster; iCluster++)
457       hEEMCAL->Fill(esd->GetCaloCluster(iCluster)->E());
458     }
459   }
460
461   // perform checks
462   if (nGen < checkNGenLow) {
463     Warning("CheckESD", "low number of generated particles: %d", Int_t(nGen));
464   }
465
466   TH1F* hEff = CreateEffHisto(hGen, hRec);
467
468   Info("CheckESD", "%d out of %d tracks reconstructed including %d "
469          "fake tracks", nRec, nGen, nFake);
470   if (nGen > 0) {
471     // efficiency
472     Double_t eff = nRec*1./nGen;
473     Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGen);
474     Double_t fake = nFake*1./nGen;
475     Double_t fakeError = TMath::Sqrt(fake*(1.-fake) / nGen);
476     Info("CheckESD", "eff = (%.1f +- %.1f) %%  fake = (%.1f +- %.1f) %%",
477          100.*eff, 100.*effError, 100.*fake, 100.*fakeError);
478
479     if (eff < checkEffLow - checkEffSigma*effError) {
480       Warning("CheckESD", "low efficiency: (%.1f +- %.1f) %%", 
481               100.*eff, 100.*effError);
482     }
483     if (fake > checkFakeHigh + checkFakeSigma*fakeError) {
484       Warning("CheckESD", "high fake: (%.1f +- %.1f) %%", 
485               100.*fake, 100.*fakeError);
486     }
487
488     // resolutions
489     Double_t res, resError;
490     if (FitHisto(hResPtInv, res, resError)) {
491       Info("CheckESD", "relative inverse pt resolution = (%.1f +- %.1f) %%",
492            res, resError);
493       if (res > checkResPtInvHigh + checkResPtInvSigma*resError) {
494         Warning("CheckESD", "bad pt resolution: (%.1f +- %.1f) %%", 
495                 res, resError);
496       }
497     }
498
499     if (FitHisto(hResPhi, res, resError)) {
500       Info("CheckESD", "phi resolution = (%.1f +- %.1f) mrad", res, resError);
501       if (res > checkResPhiHigh + checkResPhiSigma*resError) {
502         Warning("CheckESD", "bad phi resolution: (%.1f +- %.1f) mrad", 
503                 res, resError);
504       }
505     }
506
507     if (FitHisto(hResTheta, res, resError)) {
508       Info("CheckESD", "theta resolution = (%.1f +- %.1f) mrad", 
509            res, resError);
510       if (res > checkResThetaHigh + checkResThetaSigma*resError) {
511         Warning("CheckESD", "bad theta resolution: (%.1f +- %.1f) mrad", 
512                 res, resError);
513       }
514     }
515
516     // PID
517     if (nRec > 0) {
518       Double_t eff = nIdentified*1./nRec;
519       Double_t effError = TMath::Sqrt(eff*(1.-eff) / nRec);
520       Info("CheckESD", "PID eff = (%.1f +- %.1f) %%", 
521            100.*eff, 100.*effError);
522       if (eff < checkPIDEffLow - checkPIDEffSigma*effError) {
523         Warning("CheckESD", "low PID efficiency: (%.1f +- %.1f) %%", 
524                 100.*eff, 100.*effError);
525       }
526     }
527
528     printf("%9s:", "gen\\rec");
529     for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
530       printf("%9s", partName[iRec]);
531     }
532     printf("\n");
533     for (Int_t iGen = 0; iGen < AliPID::kSPECIES+1; iGen++) {
534       printf("%9s:", partName[iGen]);
535       for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
536         printf("%9d", identified[iGen][iRec]);
537       }
538       printf("\n");
539     }
540
541     if (FitHisto(hResTOFRight, res, resError)) {
542       Info("CheckESD", "TOF resolution = (%.1f +- %.1f) ps", res, resError);
543       if (res > checkResTOFHigh + checkResTOFSigma*resError) {
544         Warning("CheckESD", "bad TOF resolution: (%.1f +- %.1f) ps", 
545                 res, resError);
546       }
547     }
548
549     // calorimeters
550     if (hEPHOS->Integral() < checkPHOSNLow) {
551       Warning("CheckESD", "low number of PHOS particles: %d", 
552               Int_t(hEPHOS->Integral()));
553     } else {
554       Double_t mean = hEPHOS->GetMean();
555       if (mean < checkPHOSEnergyLow) {
556         Warning("CheckESD", "low mean PHOS energy: %.1f GeV", mean);
557       } else if (mean > checkPHOSEnergyHigh) {
558         Warning("CheckESD", "high mean PHOS energy: %.1f GeV", mean);
559       }
560     }
561
562     if (hEEMCAL->Integral() < checkEMCALNLow) {
563       Warning("CheckESD", "low number of EMCAL particles: %d", 
564               Int_t(hEEMCAL->Integral()));
565     } else {
566       Double_t mean = hEEMCAL->GetMean();
567       if (mean < checkEMCALEnergyLow) {
568         Warning("CheckESD", "low mean EMCAL energy: %.1f GeV", mean);
569       } else if (mean > checkEMCALEnergyHigh) {
570         Warning("CheckESD", "high mean EMCAL energy: %.1f GeV", mean);
571       }
572     }
573
574     // muons
575     if (hPtMUON->Integral() < checkMUONNLow) {
576       Warning("CheckESD", "low number of MUON particles: %d", 
577               Int_t(hPtMUON->Integral()));
578     } else {
579       Double_t mean = hPtMUON->GetMean();
580       if (mean < checkMUONPtLow) {
581         Warning("CheckESD", "low mean MUON pt: %.1f GeV/c", mean);
582       } else if (mean > checkMUONPtHigh) {
583         Warning("CheckESD", "high mean MUON pt: %.1f GeV/c", mean);
584       }
585     }
586
587     // V0s
588     if (nGenV0s > 0) {
589       Double_t eff = nRecV0s*1./nGenV0s;
590       Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGenV0s);
591       if (effError == 0) effError = checkV0EffLow / TMath::Sqrt(1.*nGenV0s);
592       Info("CheckESD", "V0 eff = (%.1f +- %.1f) %%", 
593            100.*eff, 100.*effError);
594       if (eff < checkV0EffLow - checkV0EffSigma*effError) {
595         Warning("CheckESD", "low V0 efficiency: (%.1f +- %.1f) %%", 
596                 100.*eff, 100.*effError);
597       }
598     }
599
600     // Cascades
601     if (nGenCascades > 0) {
602       Double_t eff = nRecCascades*1./nGenCascades;
603       Double_t effError = TMath::Sqrt(eff*(1.-eff) / nGenCascades);
604       if (effError == 0) effError = checkV0EffLow / 
605                            TMath::Sqrt(1.*nGenCascades);
606       Info("CheckESD", "Cascade eff = (%.1f +- %.1f) %%", 
607            100.*eff, 100.*effError);
608       if (eff < checkCascadeEffLow - checkCascadeEffSigma*effError) {
609         Warning("CheckESD", "low Cascade efficiency: (%.1f +- %.1f) %%", 
610                 100.*eff, 100.*effError);
611       }
612     }
613   }
614
615   // draw the histograms if not in batch mode
616   if (!gROOT->IsBatch()) {
617     new TCanvas;
618     hEff->DrawCopy();
619     new TCanvas;
620     hResPtInv->DrawCopy("E");
621     new TCanvas;
622     hResPhi->DrawCopy("E");
623     new TCanvas;
624     hResTheta->DrawCopy("E");
625     new TCanvas;
626     hDEdxRight->DrawCopy();
627     hDEdxWrong->DrawCopy("SAME");
628     new TCanvas;
629     hResTOFRight->DrawCopy("E");
630     hResTOFWrong->DrawCopy("SAME");
631     new TCanvas;
632     hEPHOS->DrawCopy("E");
633     new TCanvas;
634     hEEMCAL->DrawCopy("E");
635     new TCanvas;
636     hPtMUON->DrawCopy("E");
637     new TCanvas;
638     hMassK0->DrawCopy("E");
639     new TCanvas;
640     hMassLambda->DrawCopy("E");
641     new TCanvas;
642     hMassLambdaBar->DrawCopy("E");
643     new TCanvas;
644     hMassXi->DrawCopy("E");
645     new TCanvas;
646     hMassOmega->DrawCopy("E");
647   }
648
649   // write the output histograms to a file
650   TFile* outputFile = TFile::Open("check.root", "recreate");
651   if (!outputFile || !outputFile->IsOpen()) {
652     Error("CheckESD", "opening output file check.root failed");
653     return kFALSE;
654   }
655   hEff->Write();
656   hResPtInv->Write();
657   hResPhi->Write();
658   hResTheta->Write();
659   hDEdxRight->Write();
660   hDEdxWrong->Write();
661   hResTOFRight->Write();
662   hResTOFWrong->Write();
663   hEPHOS->Write();
664   hEEMCAL->Write();
665   hPtMUON->Write();
666   hMassK0->Write();
667   hMassLambda->Write();
668   hMassLambdaBar->Write();
669   hMassXi->Write();
670   hMassOmega->Write();
671   outputFile->Close();
672   delete outputFile;
673
674   // clean up
675   delete hGen;
676   delete hRec;
677   delete hEff;
678   delete hResPtInv;
679   delete hResPhi;
680   delete hResTheta;
681   delete hDEdxRight;
682   delete hDEdxWrong;
683   delete hResTOFRight;
684   delete hResTOFWrong;
685   delete hEPHOS;
686   delete hEEMCAL;
687   delete hPtMUON;
688   delete hMassK0;
689   delete hMassLambda;
690   delete hMassLambdaBar;
691   delete hMassXi;
692   delete hMassOmega;
693
694   delete esd;
695   esdFile->Close();
696   delete esdFile;
697
698   runLoader->UnloadHeader();
699   runLoader->UnloadKinematics();
700   delete runLoader;
701
702   // result of check
703   Info("CheckESD", "check of ESD was successfull");
704   return kTRUE;
705 }