]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/CheckESD.C
Deleting array created with TString::Tokenize
[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
27 TH1F* CreateHisto(const char* name, const char* title, 
28                   Int_t nBins, Double_t xMin, Double_t xMax,
29                   const char* xLabel = NULL, const char* yLabel = NULL)
30 {
31 // create a histogram
32
33   TH1F* result = new TH1F(name, title, nBins, xMin, xMax);
34   result->SetOption("E");
35   if (xLabel) result->GetXaxis()->SetTitle(xLabel);
36   if (yLabel) result->GetYaxis()->SetTitle(yLabel);
37   result->SetMarkerStyle(kFullCircle);
38   return result;
39 }
40
41 TH1F* CreateEffHisto(TH1F* hGen, TH1F* hRec)
42 {
43 // create an efficiency histogram
44
45   Int_t nBins = hGen->GetNbinsX();
46   TH1F* hEff = (TH1F*) hGen->Clone("hEff");
47   hEff->SetTitle("");
48   hEff->SetStats(kFALSE);
49   hEff->SetMinimum(0.);
50   hEff->SetMaximum(110.);
51   hEff->GetYaxis()->SetTitle("#epsilon [%]");
52
53   for (Int_t iBin = 0; iBin <= nBins; iBin++) {
54     Double_t nGen = hGen->GetBinContent(iBin);
55     Double_t nRec = hRec->GetBinContent(iBin);
56     if (nGen > 0) {
57       Double_t eff = nRec/nGen;
58       hEff->SetBinContent(iBin, 100. * eff);
59       Double_t error = sqrt(eff*(1.-eff) / nGen);
60       if (error == 0) error = 0.0001;
61       hEff->SetBinError(iBin, 100. * error);
62     } else {
63       hEff->SetBinContent(iBin, -100.);
64       hEff->SetBinError(iBin, 0);
65     }
66   }
67
68   return hEff;
69 }
70
71 Bool_t FitHisto(TH1* histo, Double_t& res, Double_t& resError)
72 {
73 // fit a gaussian to a histogram
74
75   static TF1* fitFunc = new TF1("fitFunc", "gaus");
76   fitFunc->SetLineWidth(2);
77   fitFunc->SetFillStyle(0);
78   Double_t maxFitRange = 2;
79
80   if (histo->Integral() > 50) {
81     Float_t mean = histo->GetMean();
82     Float_t rms = histo->GetRMS();
83     fitFunc->SetRange(mean - maxFitRange*rms, mean + maxFitRange*rms);
84     fitFunc->SetParameters(mean, rms);
85     histo->Fit(fitFunc, "QRI0");
86     histo->GetFunction("fitFunc")->ResetBit(1<<9);
87     res = TMath::Abs(fitFunc->GetParameter(2));
88     resError = TMath::Abs(fitFunc->GetParError(2));
89     return kTRUE;
90   }
91
92   return kFALSE;
93 }
94
95
96 Bool_t CheckESD(const char* gAliceFileName = "galice.root", 
97                 const char* esdFileName = "AliESDs.root")
98 {
99 // check the content of the ESD
100  
101   // check values
102   Int_t    checkNGenLow = 1;
103
104   Double_t checkEffLow = 0.5;
105   Double_t checkEffSigma = 3;
106   Double_t checkFakeHigh = 0.5;
107   Double_t checkFakeSigma = 3;
108
109   Double_t checkResPtInvHigh = 5;
110   Double_t checkResPtInvSigma = 3;
111   Double_t checkResPhiHigh = 10;
112   Double_t checkResPhiSigma = 3;
113   Double_t checkResThetaHigh = 10;
114   Double_t checkResThetaSigma = 3;
115
116   Double_t checkPIDEffLow = 0.5;
117   Double_t checkPIDEffSigma = 3;
118   Double_t checkResTOFHigh = 500;
119   Double_t checkResTOFSigma = 3;
120
121   Double_t checkPHOSNLow = 5;
122   Double_t checkPHOSEnergyLow = 0.3;
123   Double_t checkPHOSEnergyHigh = 1.0;
124   Double_t checkEMCALNLow = 50;
125   Double_t checkEMCALEnergyLow = 0.05;
126   Double_t checkEMCALEnergyHigh = 1.0;
127
128   Double_t checkMUONNLow = 1;
129   Double_t checkMUONPtLow = 0.5;
130   Double_t checkMUONPtHigh = 10.;
131
132   Double_t cutPtV0 = 0.3;
133   Double_t checkV0EffLow = 0.02;
134   Double_t checkV0EffSigma = 3;
135   Double_t cutPtCascade = 0.5;
136   Double_t checkCascadeEffLow = 0.01;
137   Double_t checkCascadeEffSigma = 3;
138
139   // open run loader and load gAlice, kinematics and header
140   AliRunLoader* runLoader = AliRunLoader::Open(gAliceFileName);
141   if (!runLoader) {
142     Error("CheckESD", "getting run loader from file %s failed", 
143             gAliceFileName);
144     return kFALSE;
145   }
146   runLoader->LoadgAlice();
147   gAlice = runLoader->GetAliRun();
148   if (!gAlice) {
149     Error("CheckESD", "no galice object found");
150     return kFALSE;
151   }
152   runLoader->LoadKinematics();
153   runLoader->LoadHeader();
154
155   // open the ESD file
156   TFile* esdFile = TFile::Open(esdFileName);
157   if (!esdFile || !esdFile->IsOpen()) {
158     Error("CheckESD", "opening ESD file %s failed", esdFileName);
159     return kFALSE;
160   }
161   AliESDEvent * esd = new AliESDEvent;
162   TTree* tree = (TTree*) esdFile->Get("esdTree");
163   if (!tree) {
164     Error("CheckESD", "no ESD tree found");
165     return kFALSE;
166   }
167   esd->ReadFromTree(tree);
168
169   // efficiency and resolution histograms
170   Int_t nBinsPt = 15;
171   Float_t minPt = 0.1;
172   Float_t maxPt = 3.1;
173   TH1F* hGen = CreateHisto("hGen", "generated tracks", 
174                            nBinsPt, minPt, maxPt, "p_{t} [GeV/c]", "N");
175   TH1F* hRec = CreateHisto("hRec", "reconstructed tracks", 
176                            nBinsPt, minPt, maxPt, "p_{t} [GeV/c]", "N");
177   Int_t nGen = 0;
178   Int_t nRec = 0;
179   Int_t nFake = 0;
180
181   TH1F* hResPtInv = CreateHisto("hResPtInv", "", 100, -10, 10, 
182            "(p_{t,rec}^{-1}-p_{t,sim}^{-1}) / p_{t,sim}^{-1} [%]", "N");
183   TH1F* hResPhi = CreateHisto("hResPhi", "", 100, -20, 20, 
184                               "#phi_{rec}-#phi_{sim} [mrad]", "N");
185   TH1F* hResTheta = CreateHisto("hResTheta", "", 100, -20, 20, 
186                                 "#theta_{rec}-#theta_{sim} [mrad]", "N");
187
188   // PID
189   Int_t partCode[AliPID::kSPECIES] = 
190     {kElectron, kMuonMinus, kPiPlus, kKPlus, kProton};
191   const char* partName[AliPID::kSPECIES+1] = 
192     {"electron", "muon", "pion", "kaon", "proton", "other"};
193   Double_t partFrac[AliPID::kSPECIES] = 
194     {0.01, 0.01, 0.85, 0.10, 0.05};
195   Int_t identified[AliPID::kSPECIES+1][AliPID::kSPECIES];
196   for (Int_t iGen = 0; iGen < AliPID::kSPECIES+1; iGen++) {
197     for (Int_t iRec = 0; iRec < AliPID::kSPECIES; iRec++) {
198       identified[iGen][iRec] = 0;
199     }
200   }
201   Int_t nIdentified = 0;
202
203   // dE/dx and TOF
204   TH2F* hDEdxRight = new TH2F("hDEdxRight", "", 300, 0, 3, 100, 0, 400);
205   hDEdxRight->SetStats(kFALSE);
206   hDEdxRight->GetXaxis()->SetTitle("p [GeV/c]");
207   hDEdxRight->GetYaxis()->SetTitle("dE/dx_{TPC}");
208   hDEdxRight->SetMarkerStyle(kFullCircle);
209   hDEdxRight->SetMarkerSize(0.4);
210   TH2F* hDEdxWrong = new TH2F("hDEdxWrong", "", 300, 0, 3, 100, 0, 400);
211   hDEdxWrong->SetStats(kFALSE);
212   hDEdxWrong->GetXaxis()->SetTitle("p [GeV/c]");
213   hDEdxWrong->GetYaxis()->SetTitle("dE/dx_{TPC}");
214   hDEdxWrong->SetMarkerStyle(kFullCircle);
215   hDEdxWrong->SetMarkerSize(0.4);
216   hDEdxWrong->SetMarkerColor(kRed);
217   TH1F* hResTOFRight = CreateHisto("hResTOFRight", "", 100, -1000, 1000, 
218                                    "t_{TOF}-t_{track} [ps]", "N");
219   TH1F* hResTOFWrong = CreateHisto("hResTOFWrong", "", 100, -1000, 1000, 
220                                    "t_{TOF}-t_{track} [ps]", "N");
221   hResTOFWrong->SetLineColor(kRed);
222
223   // calorimeters
224   TH1F* hEPHOS = CreateHisto("hEPHOS", "PHOS", 100, 0, 50, "E [GeV]", "N");
225   TH1F* hEEMCAL = CreateHisto("hEEMCAL", "EMCAL", 100, 0, 50, "E [GeV]", "N");
226
227   // muons
228   TH1F* hPtMUON = CreateHisto("hPtMUON", "MUON", 100, 0, 20, 
229                               "p_{t} [GeV/c]", "N");
230
231   // V0s and cascades
232   TH1F* hMassK0 = CreateHisto("hMassK0", "K^{0}", 100, 0.4, 0.6, 
233                               "M(#pi^{+}#pi^{-}) [GeV/c^{2}]", "N");
234   TH1F* hMassLambda = CreateHisto("hMassLambda", "#Lambda", 100, 1.0, 1.2, 
235                                   "M(p#pi^{-}) [GeV/c^{2}]", "N");
236   TH1F* hMassLambdaBar = CreateHisto("hMassLambdaBar", "#bar{#Lambda}", 
237                                      100, 1.0, 1.2, 
238                                      "M(#bar{p}#pi^{+}) [GeV/c^{2}]", "N");
239   Int_t nGenV0s = 0;
240   Int_t nRecV0s = 0;
241   TH1F* hMassXi = CreateHisto("hMassXi", "#Xi", 100, 1.2, 1.5, 
242                               "M(#Lambda#pi) [GeV/c^{2}]", "N");
243   TH1F* hMassOmega = CreateHisto("hMassOmega", "#Omega", 100, 1.5, 1.8, 
244                                  "M(#LambdaK) [GeV/c^{2}]", "N");
245   Int_t nGenCascades = 0;
246   Int_t nRecCascades = 0;
247
248   // loop over events
249   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
250     runLoader->GetEvent(iEvent);
251
252     // select simulated primary particles, V0s and cascades
253     AliStack* stack = runLoader->Stack();
254     Int_t nParticles = stack->GetNtrack();
255     TArrayF vertex(3);
256     runLoader->GetHeader()->GenEventHeader()->PrimaryVertex(vertex);
257     TObjArray selParticles;
258     TObjArray selV0s;
259     TObjArray selCascades;
260     for (Int_t iParticle = 0; iParticle < nParticles; iParticle++) {
261       TParticle* particle = stack->Particle(iParticle);
262       if (!particle) continue;
263       if (particle->Pt() < 0.001) continue;
264       if (TMath::Abs(particle->Eta()) > 0.9) continue;
265       TVector3 dVertex(particle->Vx() - vertex[0], 
266                        particle->Vy() - vertex[1],
267                        particle->Vz() - vertex[2]);
268       if (dVertex.Mag() > 0.0001) continue;
269
270       switch (TMath::Abs(particle->GetPdgCode())) {
271       case kElectron:
272       case kMuonMinus:
273       case kPiPlus:
274       case kKPlus:
275       case kProton: {
276         if (particle->Pt() > minPt) {
277           selParticles.Add(particle);
278           nGen++;
279           hGen->Fill(particle->Pt());
280         }
281         break;
282       }
283       case kK0Short:
284       case kLambda0: {
285         if (particle->Pt() > cutPtV0) {
286           nGenV0s++;
287           selV0s.Add(particle);
288         }
289         break;
290       }
291       case kXiMinus:
292       case kOmegaMinus: {
293         if (particle->Pt() > cutPtCascade) {
294           nGenCascades++;
295           selCascades.Add(particle);
296         }
297         break;
298       }
299       default: break;
300       }
301     }
302
303     // get the event summary data
304     tree->GetEvent(iEvent);
305     if (!esd) {
306       Error("CheckESD", "no ESD object found for event %d", iEvent);
307       return kFALSE;
308     }
309
310     // loop over tracks
311     for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
312       AliESDtrack* track = esd->GetTrack(iTrack);
313
314       // select tracks of selected particles
315       Int_t label = TMath::Abs(track->GetLabel());
316       if (label > stack->GetNtrack()) continue;     // background
317       TParticle* particle = stack->Particle(label);
318       if (!selParticles.Contains(particle)) continue;
319       if ((track->GetStatus() & AliESDtrack::kITSrefit) == 0) continue;
320       if (track->GetConstrainedChi2() > 1e9) continue;
321       selParticles.Remove(particle);   // don't count multiple tracks
322
323       nRec++;
324       hRec->Fill(particle->Pt());
325       if (track->GetLabel() < 0) nFake++;
326
327       // resolutions
328       hResPtInv->Fill(100. * (TMath::Abs(track->GetSigned1Pt()) - 1./particle->Pt()) * 
329                       particle->Pt());
330       hResPhi->Fill(1000. * (track->Phi() - particle->Phi()));
331       hResTheta->Fill(1000. * (track->Theta() - particle->Theta()));
332
333       // PID
334       if ((track->GetStatus() & AliESDtrack::kESDpid) == 0) continue;
335       Int_t iGen = 5;
336       for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
337         if (TMath::Abs(particle->GetPdgCode()) == partCode[i]) iGen = i;
338       }
339       Double_t probability[AliPID::kSPECIES];
340       track->GetESDpid(probability);
341       Double_t pMax = 0;
342       Int_t iRec = 0;
343       for (Int_t i = 0; i < AliPID::kSPECIES; i++) {
344         probability[i] *= partFrac[i];
345         if (probability[i] > pMax) {
346           pMax = probability[i];
347           iRec = i;
348         }
349       }
350       identified[iGen][iRec]++;
351       if (iGen == iRec) nIdentified++;
352
353       // dE/dx and TOF
354       Double_t time[AliPID::kSPECIES];
355       track->GetIntegratedTimes(time);
356       if (iGen == iRec) {
357         hDEdxRight->Fill(particle->P(), track->GetTPCsignal());
358         if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
359           hResTOFRight->Fill(track->GetTOFsignal() - time[iRec]);
360         }
361       } else {
362         hDEdxWrong->Fill(particle->P(), track->GetTPCsignal());
363         if ((track->GetStatus() & AliESDtrack::kTOFpid) != 0) {
364           hResTOFWrong->Fill(track->GetTOFsignal() - time[iRec]);
365         }
366       }
367     }
368
369     // loop over muon tracks
370     {
371     for (Int_t iTrack = 0; iTrack < esd->GetNumberOfMuonTracks(); iTrack++) {
372       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
373       Double_t ptInv = TMath::Abs(muonTrack->GetInverseBendingMomentum());
374       if (ptInv > 0.001) {
375         hPtMUON->Fill(1./ptInv);
376       }
377     }
378     }
379
380     // loop over V0s
381     for (Int_t iV0 = 0; iV0 < esd->GetNumberOfV0s(); iV0++) {
382       AliESDv0* v0 = esd->GetV0(iV0);
383       if (v0->GetOnFlyStatus()) continue;
384       v0->ChangeMassHypothesis(kK0Short);
385       hMassK0->Fill(v0->GetEffMass());
386       v0->ChangeMassHypothesis(kLambda0);
387       hMassLambda->Fill(v0->GetEffMass());
388       v0->ChangeMassHypothesis(kLambda0Bar);
389       hMassLambdaBar->Fill(v0->GetEffMass());
390
391       Int_t negLabel = TMath::Abs(esd->GetTrack(v0->GetNindex())->GetLabel());
392       if (negLabel > stack->GetNtrack()) continue;     // background
393       Int_t negMother = stack->Particle(negLabel)->GetMother(0);
394       if (negMother < 0) continue;
395       Int_t posLabel = TMath::Abs(esd->GetTrack(v0->GetPindex())->GetLabel());
396       if (posLabel > stack->GetNtrack()) continue;     // background
397       Int_t posMother = stack->Particle(posLabel)->GetMother(0);
398       if (negMother != posMother) continue;
399       TParticle* particle = stack->Particle(negMother);
400       if (!selV0s.Contains(particle)) continue;
401       selV0s.Remove(particle);
402       nRecV0s++;
403     }
404
405     // loop over Cascades
406     for (Int_t iCascade = 0; iCascade < esd->GetNumberOfCascades(); 
407          iCascade++) {
408       AliESDcascade* cascade = esd->GetCascade(iCascade);
409       Double_t v0q;
410       cascade->ChangeMassHypothesis(v0q,kXiMinus);
411       hMassXi->Fill(cascade->GetEffMassXi());
412       cascade->ChangeMassHypothesis(v0q,kOmegaMinus);
413       hMassOmega->Fill(cascade->GetEffMassXi());
414
415       Int_t negLabel = TMath::Abs(esd->GetTrack(cascade->GetNindex())
416                                   ->GetLabel());
417       if (negLabel > stack->GetNtrack()) continue;     // background
418       Int_t negMother = stack->Particle(negLabel)->GetMother(0);
419       if (negMother < 0) continue;
420       Int_t posLabel = TMath::Abs(esd->GetTrack(cascade->GetPindex())
421                                   ->GetLabel());
422       if (posLabel > stack->GetNtrack()) continue;     // background
423       Int_t posMother = stack->Particle(posLabel)->GetMother(0);
424       if (negMother != posMother) continue;
425       Int_t v0Mother = stack->Particle(negMother)->GetMother(0);
426       if (v0Mother < 0) continue;
427       Int_t bacLabel = TMath::Abs(esd->GetTrack(cascade->GetBindex())
428                                   ->GetLabel());
429       if (bacLabel > stack->GetNtrack()) continue;     // background
430       Int_t bacMother = stack->Particle(bacLabel)->GetMother(0);
431       if (v0Mother != bacMother) continue;
432       TParticle* particle = stack->Particle(v0Mother);
433       if (!selCascades.Contains(particle)) continue;
434       selCascades.Remove(particle);
435       nRecCascades++;
436     }
437
438     // loop over the clusters
439     {
440       for (Int_t iCluster=0; iCluster<esd->GetNumberOfCaloClusters(); iCluster++) {
441         AliESDCaloCluster * clust = esd->GetCaloCluster(iCluster);
442         if (clust->IsPHOS()) hEPHOS->Fill(clust->E());
443         if (clust->IsEMCAL()) hEEMCAL->Fill(clust->E());
444       }
445     }
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 < AliPID::kSPECIES; iRec++) {
518       printf("%9s", partName[iRec]);
519     }
520     printf("\n");
521     for (Int_t iGen = 0; iGen < AliPID::kSPECIES+1; iGen++) {
522       printf("%9s:", partName[iGen]);
523       for (Int_t iRec = 0; iRec < AliPID::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 }