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