]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/FMDrec/AliFMDQAChecker.cxx
Dummy implementation
[u/mrichter/AliRoot.git] / FMD / FMDrec / AliFMDQAChecker.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 //__________________________________________________________________
16 //
17 // Yves?
18 // What 
19 // is 
20 // this 
21 // class 
22 // supposed 
23 // to
24 // do?
25 //__________________________________________________________________
26 //
27 // --- ROOT system ---
28 #include <TClass.h>
29 #include <TH1F.h> 
30 #include <TF1.h> 
31 #include <TH2.h> 
32 #include <TH1I.h> 
33 #include <TIterator.h> 
34 #include <TKey.h> 
35 #include <TFile.h> 
36 #include <iostream>
37 #include <TCanvas.h>
38 #include <TStyle.h>
39 #include <TLatex.h>
40 #include <TFitResult.h>
41 #include <TParameter.h>
42 #include <TMacro.h>
43 #include <TPaveText.h>
44 #include <TVirtualFitter.h>
45
46 // --- AliRoot header files ---
47 #include "AliLog.h"
48 #include "AliQAv1.h"
49 #include "AliQAChecker.h"
50 #include "AliFMDQAChecker.h"
51 #include "AliFMDQADataMakerRec.h"
52 #include "AliRecoParam.h"
53 #include <AliCDBManager.h>
54 #include <AliCDBEntry.h>
55 #include <AliCDBId.h>
56 #include <AliQAThresholds.h>
57
58 ClassImp(AliFMDQAChecker)
59 #if 0
60 ; // This is for Emacs! - do not delete
61 #endif
62
63 namespace {
64   void addFitsMacro(TList* l) {
65     TMacro* m = new TMacro("fits");
66     m->AddLine("void fits() {");
67     m->AddLine("  if (!gPad) { Printf(\"No gPad\"); return; }");
68     m->AddLine("  TList* lp = gPad->GetListOfPrimitives();");
69     m->AddLine("  if (!lp) return;");
70     m->AddLine("  TObject* po = 0;");
71     m->AddLine("  TIter next(lp);");
72     m->AddLine("  while ((po = next())) {");
73     m->AddLine("    if (!po->IsA()->InheritsFrom(TH1::Class())) continue;");
74     m->AddLine("    TH1*   htmp = dynamic_cast<TH1*>(po);");
75     m->AddLine("    TList* lf   = htmp->GetListOfFunctions();");
76     m->AddLine("    TObject* pso = (lf ? lf->FindObject(\"stats\") : 0);");
77     m->AddLine("    if (!pso) continue;");
78     m->AddLine("    TPaveStats* ps = static_cast<TPaveStats*>(pso);");
79     m->AddLine("    ps->SetOptFit(111);");
80     m->AddLine("    UShort_t qual = htmp->GetUniqueID();");
81     m->AddLine("    ps->SetFillColor(qual >= 3 ? kRed-4 : qual >= 2 ? kOrange-4 : qual >= 1 ? kYellow-4 : kGreen-4);");
82     // m->AddLine("    lf->Remove(lf->FindObject(\"fits\"));");
83     // m->AddLine("    ps->Paint();");
84     m->AddLine("    break;");
85     m->AddLine("  }");
86     // m->AddLine("  gPad->Modified(); gPad->Update(); gPad->cd();");
87     m->AddLine("}");
88
89     TObject* old = l->FindObject(m->GetName());
90     if (old) l->Remove(old);
91     l->Add(m);
92   }
93   
94   const Double_t kROErrorsLabelY    = .30;
95   
96   const Int_t    kConvolutionSteps  = 100;
97   const Double_t kConvolutionNSigma = 5;
98
99   // 
100   // The shift of the most probable value for the ROOT function TMath::Landau 
101   //
102   const Double_t  kMpShift  = -0.22278298;
103   // 
104   // Integration normalisation 
105   //
106   const Double_t  kInvSq2Pi = 1. / TMath::Sqrt(2*TMath::Pi());
107
108   Double_t landau(Double_t x, Double_t delta, Double_t xi)
109   {
110     // 
111     // Calculate the shifted Landau
112     // @f[
113     //    f'_{L}(x;\Delta,\xi) = f_L(x;\Delta+0.22278298\xi)
114     // @f]
115     //
116     // where @f$ f_{L}@f$ is the ROOT implementation of the Landau
117     // distribution (known to have @f$ \Delta_{p}=-0.22278298@f$ for
118     // @f$\Delta=0,\xi=1@f$.
119     // 
120     // Parameters:
121     //    x      Where to evaluate @f$ f'_{L}@f$ 
122     //    delta  Most probable value 
123     //    xi     The 'width' of the distribution 
124     //
125     // Return:
126     //    @f$ f'_{L}(x;\Delta,\xi) @f$
127     //
128     return TMath::Landau(x, delta - xi * kMpShift, xi);
129   }
130   Double_t landauGaus(Double_t x, Double_t delta, Double_t xi,
131                       Double_t sigma, Double_t sigmaN)
132   {
133     // 
134     // Calculate the value of a Landau convolved with a Gaussian 
135     // 
136     // @f[ 
137     // f(x;\Delta,\xi,\sigma') = \frac{1}{\sigma' \sqrt{2 \pi}}
138     //    \int_{-\infty}^{+\infty} d\Delta' f'_{L}(x;\Delta',\xi)
139     //    \exp{-\frac{(\Delta-\Delta')^2}{2\sigma'^2}}
140     // @f]
141     // 
142     // where @f$ f'_{L}@f$ is the Landau distribution, @f$ \Delta@f$
143     // the energy loss, @f$ \xi@f$ the width of the Landau, and @f$
144     // \sigma'^2=\sigma^2-\sigma_n^2 @f$.  Here, @f$\sigma@f$ is the
145     // variance of the Gaussian, and @f$\sigma_n@f$ is a parameter
146     // modelling noise in the detector.
147     //
148     // Note that this function uses the constants kConvolutionSteps and
149     // kConvolutionNSigma
150     // 
151     // References: 
152     //  - <a href="http://dx.doi.org/10.1016/0168-583X(84)90472-5">Nucl.Instrum.Meth.B1:16</a>
153     //  - <a href="http://dx.doi.org/10.1103/PhysRevA.28.615">Phys.Rev.A28:615</a>
154     //  - <a href="http://root.cern.ch/root/htmldoc/tutorials/fit/langaus.C.html">ROOT implementation</a>
155     // 
156     // Parameters:
157     //    x         where to evaluate @f$ f@f$
158     //    delta     @f$ \Delta@f$ of @f$ f(x;\Delta,\xi,\sigma')@f$
159     //    xi        @f$ \xi@f$ of @f$ f(x;\Delta,\xi,\sigma')@f$
160     //    sigma     @f$ \sigma@f$ of @f$\sigma'^2=\sigma^2-\sigma_n^2 @f$
161     //    sigma_n   @f$ \sigma_n@f$ of @f$\sigma'^2=\sigma^2-\sigma_n^2 @f$
162     // 
163     // Return:
164     //    @f$ f@f$ evaluated at @f$ x@f$.  
165     //
166     Double_t deltap = delta - xi * kMpShift;
167     Double_t sigma2 = sigmaN*sigmaN + sigma*sigma;
168     Double_t sigma1 = sigmaN == 0 ? sigma : TMath::Sqrt(sigma2);
169     Double_t xlow   = x - kConvolutionNSigma * sigma1;
170     Double_t xhigh  = x + kConvolutionNSigma * sigma1;
171     Double_t step   = (xhigh - xlow) / kConvolutionSteps;
172     Double_t sum    = 0;
173   
174     for (Int_t i = 0; i <= kConvolutionSteps/2; i++) { 
175       Double_t x1 = xlow  + (i - .5) * step;
176       Double_t x2 = xhigh - (i - .5) * step;
177     
178       sum += TMath::Landau(x1, deltap, xi, kTRUE) * TMath::Gaus(x, x1, sigma1);
179       sum += TMath::Landau(x2, deltap, xi, kTRUE) * TMath::Gaus(x, x2, sigma1);
180     }
181     return step * sum * kInvSq2Pi / sigma1;
182   }
183
184   // 
185   // Utility function to use in TF1 defintition 
186   //
187   Double_t landauGaus1(Double_t* xp, Double_t* pp) 
188   {
189     Double_t x        = xp[0];
190     Double_t constant = pp[0];
191     Double_t delta    = pp[1];
192     Double_t xi       = pp[2];
193     Double_t sigma    = pp[3];
194     Double_t sigmaN   = pp[4];
195
196     return constant * landauGaus(x, delta, xi, sigma, sigmaN);
197   }
198
199   //____________________________________________________________________
200   TF1* makeLandauGaus(const char* , 
201                       Double_t  c=1, 
202                       Double_t  delta=.5,  Double_t xi=0.07, 
203                       Double_t  sigma=.1,  Double_t sigmaN=-1, 
204                       Double_t  xmin=0,    Double_t xmax=15)
205   {
206     // 
207     // Generate a TF1 object of @f$ f_I@f$ 
208     // 
209     // Parameters:
210     //    c        Constant
211     //    delta    @f$ \Delta@f$ 
212     //    xi       @f$ \xi_1@f$        
213     //    sigma    @f$ \sigma_1@f$             
214     //    sigma_n  @f$ \sigma_n@f$             
215     //    xmin     Least value of range
216     //    xmax     Largest value of range
217     // 
218     // Return:
219     //    Newly allocated TF1 object
220     //
221     Int_t npar     = 5;
222     TF1*  func     = new TF1("landauGaus", 
223                              &landauGaus1,xmin,xmax,npar);
224     // func->SetLineStyle(((i-2) % 10)+2); // start at dashed
225     func->SetLineColor(kBlack); 
226     func->SetLineWidth(2);
227     func->SetNpx(500);
228     func->SetParNames("C","#Delta_{p}","#xi", "#sigma", "#sigma_{n}");
229
230     // Set the initial parameters from the seed fit 
231     func->SetParameter(0,      c);       
232     func->SetParameter(1,  delta);   
233     func->SetParameter(2,     xi);      
234     func->SetParameter(3,  sigma);   
235     func->SetParameter(4, sigmaN); 
236
237     func->SetParLimits(1, 0,    xmax);
238     func->SetParLimits(2, 0,    xmax);
239     func->SetParLimits(3, 0.01, 1);
240     
241     if (sigmaN < 0) func->FixParameter(4, 0);
242     else            func->SetParLimits(4, 0, xmax);
243
244     return func;
245   }
246 }
247   
248 //__________________________________________________________________
249 AliFMDQAChecker::AliFMDQAChecker() 
250   : AliQACheckerBase("FMD","FMD Quality Assurance Checker") ,
251     fDoScale(false),
252     fDidExternal(false), 
253     fShowFitResults(true),
254     fELossLowCut(0.2), 
255     fELossNRMS(3), 
256     fELossBadChi2Nu(10), 
257     fELossFkupChi2Nu(100), 
258     fELossMinEntries(1000),
259     fELossMaxEntries(-1),
260     fELossGoodParError(0.1),
261     fELossMinSharing(0.1),
262     fROErrorsBad(0.3), 
263     fROErrorsFkup(0.5),
264     fMaxNProblem(10),
265     fMaxNBad(10),
266     fNoFits(false)
267 {
268 }          
269
270 //__________________________________________________________________
271 void
272 AliFMDQAChecker::ProcessExternalParams()
273 {
274   ProcessExternalParam("ELossLowCut",           fELossLowCut);
275   ProcessExternalParam("ELossNRMS",             fELossNRMS);
276   ProcessExternalParam("ELossBadChi2Nu",        fELossBadChi2Nu);
277   ProcessExternalParam("ELossFkupChi2Nu",       fELossFkupChi2Nu);
278   ProcessExternalParam("ELossGoodParError",     fELossGoodParError);
279   ProcessExternalParam("ROErrorsBad",           fROErrorsBad);
280   ProcessExternalParam("ROErrorsFkup",          fROErrorsFkup);
281   ProcessExternalParam("ELossMinSharing",       fELossMinSharing);
282   Double_t tmp = 0;
283   ProcessExternalParam("CommonScale", tmp);
284   fDoScale = tmp > 0; tmp = fELossMinEntries;
285   ProcessExternalParam("ELossMinEntries", tmp);
286   fELossMinEntries = tmp; tmp = fELossMaxEntries;
287   ProcessExternalParam("ELossMaxEntries", tmp);
288   fELossMaxEntries = tmp; tmp = fMaxNProblem;
289   ProcessExternalParam("MaxNProblem", tmp);
290   fMaxNProblem = tmp; tmp = 0;
291   fELossMaxEntries = tmp; tmp = fMaxNBad;
292   ProcessExternalParam("MaxNBad", tmp);
293   fMaxNBad = tmp; tmp = 0;
294   ProcessExternalParam("NoFits", tmp);
295   fNoFits = tmp > 0; tmp = 0;
296
297   GetThresholds();
298
299   fDidExternal = true;
300 }
301 //__________________________________________________________________
302 void
303 AliFMDQAChecker::ProcessExternalParam(const char* name, Double_t& v)
304 {
305   TObject* o = fExternParamList->FindObject(name);
306   if (!o) return; 
307   TParameter<double>* p = static_cast<TParameter<double>*>(o);
308   v = p->GetVal();
309   AliDebugF(3, "External parameter: %-20s=%lf", name, v);
310 }
311
312 //__________________________________________________________________
313 void
314 AliFMDQAChecker::GetThresholds()
315 {
316   const char*    path   = "GRP/Calib/QAThresholds";
317   AliCDBManager* cdbMan = AliCDBManager::Instance();
318   AliCDBEntry*   cdbEnt = cdbMan->Get(path);
319   if (!cdbEnt) { 
320     AliWarningF("Failed to get CDB entry at %s", path);
321     return;
322   }
323   
324   TObjArray*     cdbObj = static_cast<TObjArray*>(cdbEnt->GetObject());
325   if (!cdbObj) { 
326     AliWarningF("Failed to get CDB object at %s", path);
327     return;
328   }
329   
330   TObject*       fmdObj = cdbObj->FindObject("FMD");
331   if (!fmdObj) { 
332     AliWarningF("Failed to get FMD object at from CDB %s", path);
333     return;
334   }
335   
336   AliQAThresholds* qaThr = static_cast<AliQAThresholds*>(fmdObj);
337   Int_t nThr = qaThr->GetSize();
338   for (Int_t i = 0; i < nThr; i++) { 
339     TObject* thr = qaThr->GetThreshold(i);
340     if (!thr) continue; 
341
342     TParameter<double>* d = dynamic_cast<TParameter<double>*>(thr);
343     if (!d) { 
344       AliWarningF("Parameter %s not of type double", thr->GetName());
345       continue;
346     }
347     Double_t val = d->GetVal();
348     TString  name(thr->GetName());
349     if      (name.EqualTo("ELossBadChi2Nu"))     fELossBadChi2Nu    = val;
350     else if (name.EqualTo("ELossFkupChi2Nu"))    fELossFkupChi2Nu   = val;
351     else if (name.EqualTo("ELossGoodParError"))  fELossGoodParError = val;
352     else if (name.EqualTo("ROErrorsBad"))        fROErrorsBad       = val;
353     else if (name.EqualTo("ROErrorsFkup"))       fROErrorsFkup      = val;    
354     else if (name.EqualTo("MaxNProblem"))        fMaxNProblem       = val;
355     else if (name.EqualTo("MaxNBad"))            fMaxNBad           = val;
356     AliDebugF(3, "Threshold %s=%f", name.Data(), val);
357   }
358 }
359
360 //__________________________________________________________________
361 AliQAv1::QABIT_t
362 AliFMDQAChecker::Quality2Bit(UShort_t value) const
363 {
364   AliQAv1::QABIT_t  ret   = AliQAv1::kINFO; // Assume success 
365   if      (value >= kWhatTheFk) ret = AliQAv1::kFATAL;
366   else if (value >= kBad)       ret = AliQAv1::kERROR;
367   else if (value >= kProblem)   ret = AliQAv1::kWARNING;
368   
369   return ret;
370 }
371
372 //__________________________________________________________________
373 void
374 AliFMDQAChecker::SetQA(AliQAv1::ALITASK_t index, Double_t* values) const
375 {
376   AliQAv1 * qa = AliQAv1::Instance(index) ;
377
378   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
379     // Check if specie is defined 
380     if (!qa->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)))
381       continue ;
382     
383     // No checker is implemented, set all QA to Fatal
384     if (!values) { 
385       qa->Set(AliQAv1::kFATAL, specie) ; 
386       continue;
387     }
388     
389     UShort_t          value = values[specie];
390     AliQAv1::QABIT_t  ret   = Quality2Bit(value);
391     
392     qa->Set(ret, AliRecoParam::ConvertIndex(specie));
393     AliDebugF(3, "Quality of %s: %d -> %d", 
394               AliRecoParam::GetEventSpecieName(specie), value, ret);
395   }
396 }
397
398 //__________________________________________________________________
399 UShort_t
400 AliFMDQAChecker::BasicCheck(TH1* hist) const
401 {
402   if (hist->GetEntries() <= 0) return kOK;
403   return (hist->GetMean() > 0 ? kOK : kProblem);
404 }
405
406 //__________________________________________________________________
407 UShort_t
408 AliFMDQAChecker::CheckOne(AliQAv1::ALITASK_t          what,
409                           AliRecoParam::EventSpecie_t specie, 
410                           TH1*                        hist) const
411 {
412   if(what == AliQAv1::kESD) return CheckESD(specie, hist);
413   if(what == AliQAv1::kRAW) return CheckRaw(specie, hist);
414   if(what == AliQAv1::kSIM) return CheckSim(specie, hist);
415   if(what == AliQAv1::kREC) return CheckRec(specie, hist);
416   return 0;
417 }
418 //__________________________________________________________________
419 UShort_t 
420 AliFMDQAChecker::CheckESD(AliRecoParam::EventSpecie_t /* specie*/, 
421                           TH1*                        hist) const
422 {
423   return BasicCheck(hist);
424 }
425 namespace {
426   Double_t Chi2Scale(TH1* h, Double_t base=10000) 
427   {
428     return 1. / TMath::Max(1., h->GetEntries() / base);
429   }
430   void AddLine(TObjArray* lines, 
431                Double_t x1, Double_t x2, Double_t x3, 
432                Double_t y, Double_t dy,
433                const char* name, Double_t val, Double_t lim, 
434                Bool_t ok, Int_t color)
435   {
436     TString n; n.Form("%s:", name);
437     TLatex* ltx = new TLatex(x1, y, n);
438     ltx->SetNDC(true);
439     ltx->SetTextSize(dy-0.01);
440     ltx->SetTextColor(color);
441     lines->Add(ltx);
442     
443     n.Form("%7.3f", val);
444     ltx = new TLatex(x2, y, n);
445     ltx->SetNDC(true);
446     ltx->SetTextSize(dy-0.01);
447     ltx->SetTextColor(color);
448     lines->Add(ltx);
449     
450     if (lim < 0) n = "(ignored)";
451     else  n.Form("%c %4.2f", ok ? '<' : '>', lim);
452     ltx = new TLatex(x3, y, n);
453     ltx->SetNDC(true);
454     ltx->SetTextSize(dy-0.01);
455     ltx->SetTextColor(color);
456     lines->Add(ltx);
457   }
458 }
459
460 //__________________________________________________________________
461 UShort_t
462 AliFMDQAChecker::CheckFit(TH1* hist, const TFitResultPtr& res, 
463                           Double_t low, Double_t high, Int_t& color) const
464 {
465   color = kGreen+4;
466
467   // Check if there's indeed a result - if not, flag as OK
468   if (!res.Get()) return 0;
469
470   UShort_t   ret   = 0;
471   Int_t      nPar  = res->NPar();
472   Double_t   dy    = .06;
473   Double_t   x     = .2;
474   Double_t   x2    = .3;
475   Double_t   x3    = .4;
476   Double_t   y     = .9-dy;
477   Double_t   chi2  = res->Chi2();
478   Int_t      nu    = res->Ndf();
479   Double_t   s     = Chi2Scale(hist,fELossMinEntries);
480   Double_t   red   = (nu == 0 ? fELossFkupChi2Nu : chi2 / nu);
481   TObjArray* lines = 0;
482   // TLatex*    lRed  = 0;
483   TLatex*    ltx   = 0;
484   Int_t      chi2Check = 0;
485   Double_t   chi2Lim   = fELossBadChi2Nu;
486   if (AliDebugLevel() > 0) 
487     printf("FIT: %s, 1, %d, %f, %f\n", hist->GetName(),
488             Int_t(hist->GetEntries()), red, s * red);
489   red *= s;
490   if (red > fELossBadChi2Nu) { // || res->Prob() < .01) { 
491     // AliWarningF("Fit gave chi^2/nu=%f/%d=%f>%f (%f)", 
492     //             res->Chi2(), res->Ndf(), red, fELossBadChi2Nu, 
493     //             fELossFkupChi2Nu);
494     // res->Print();
495     chi2Check++;
496     if (red > fELossFkupChi2Nu) { 
497       chi2Check++;
498       chi2Lim = fELossFkupChi2Nu;
499     }
500   }
501   ret += chi2Check;
502
503   if (fShowFitResults) { 
504     lines = new TObjArray(nPar+3);
505     lines->SetName("lines");
506     lines->SetOwner(true);
507
508     AddLine(lines, x, x2, x3, y, dy, "#chi^{2}/#nu", red, chi2Lim,
509             chi2Check < 1, chi2Check < 1 ? color : 
510             chi2Check < 2 ? kOrange+2 : kRed+2);
511     
512     Double_t x1 = .85;
513     Double_t y1 = .5;
514
515     // y1 -= dy;
516     ltx = new TLatex(x1, y1, Form("Fit range: [%6.2f,%6.2f]", low, high));
517     ltx->SetTextColor(kGray+3);
518     ltx->SetTextSize(dy-.01);
519     ltx->SetTextAlign(31);
520     ltx->SetNDC(true);
521     lines->Add(ltx);
522
523     y1 -= dy;
524     ltx = new TLatex(x1, y1, Form("Entries: %d (%d)", 
525                                   Int_t(hist->GetEffectiveEntries()),
526                                   fELossMaxEntries));
527     ltx->SetTextColor(kGray+3);
528     ltx->SetTextSize(dy-.01);
529     ltx->SetTextAlign(31);
530     ltx->SetNDC(true);
531     lines->Add(ltx);
532
533     y1 -= dy;
534     ltx = new TLatex(x1, y1, Form("%s: %f #pm %f", 
535                                   res->ParName(1).c_str(),
536                                   res->Parameter(1),
537                                   res->ParError(1)));
538     ltx->SetTextColor(kGray+3);
539     ltx->SetTextSize(dy-.01);
540     ltx->SetTextAlign(31);
541     ltx->SetNDC(true);
542     lines->Add(ltx);
543   }
544   
545   // Now check the relative error on the fit parameters 
546   Int_t parsOk = 0;
547   for (Int_t i = 0; i < nPar; i++) { 
548     if (res->IsParameterFixed(i)) continue; 
549     Double_t thr = fELossGoodParError;
550     Double_t pv  = res->Parameter(i);
551     Double_t pe  = res->ParError(i);
552     Double_t rel = (pv == 0 ? 100 : pe / pv);
553     Bool_t   ok  = (i == 3) || (rel < thr);
554     if (lines) {
555       y -= dy;
556       AddLine(lines, x, x2, x3, y, dy,Form("#delta%s/%s", 
557                                            res->ParName(i).c_str(),
558                                            res->ParName(i).c_str()),
559               rel, (i == 3 ? -1 : thr), ok, ok ? color : kOrange+2);
560     }
561     if (i == 3) continue; // Skip sigma 
562     if (ok) parsOk++;
563   }
564   if (parsOk > 0) 
565     ret = TMath::Max(ret-(parsOk-1),0);
566   if (ret > 1) color = kRed+2;
567   if (ret > 0) color = kOrange+2;
568
569   if (lines) {
570     TList*   lf  = hist->GetListOfFunctions();
571     TObject* old = lf->FindObject(lines->GetName());
572     if (old) {
573       lf->Remove(old);
574       delete old;
575     }
576     lf->Add(lines);
577   }
578   hist->SetStats(false);
579     
580   return ret;
581 }
582
583 //__________________________________________________________________
584 UShort_t
585 AliFMDQAChecker::CheckRaw(AliRecoParam::EventSpecie_t specie, 
586                           TH1*                        hist) const
587 {
588   Int_t ret = BasicCheck(hist);
589   TString name(hist->GetName());
590   if (name.Contains("readouterrors", TString::kIgnoreCase)) { 
591     // Check the mean number of errors per event 
592     TH2*  roErrors = static_cast<TH2*>(hist);
593     Int_t nY       = roErrors->GetNbinsY();
594
595     TLatex* ltx = new TLatex(.15, .9, Form("Thresholds: %5.2f,%5.2f",
596                                            fROErrorsBad, fROErrorsFkup));
597     ltx->SetName("thresholds");
598     ltx->SetTextColor(kGray+3);
599     ltx->SetNDC();
600
601     TList* ll = hist->GetListOfFunctions(); 
602     TObject* old = ll->FindObject(ltx->GetName());
603     if (old) { 
604       ll->Remove(old);
605       delete old;
606     }
607     ll->Add(ltx);
608     
609     for (Int_t i = 1; i <= 3; i++) {
610       Double_t sum = 0;
611       Int_t    cnt = 0;
612       for (Int_t j = 1; j <= nY; j++) {
613         Int_t n =  roErrors->GetBinContent(i, j);
614         sum     += n * roErrors->GetYaxis()->GetBinCenter(j);
615         cnt     += n;
616       }
617       Double_t mean = (cnt <= 0 ? 0 : sum / cnt);
618       Double_t x    = ((i-.5) * (1-0.1-0.1) / 3 + 0.1);
619       
620       ltx = new TLatex(x, kROErrorsLabelY, Form("Mean: %6.3f", mean));
621       ltx->SetName(Form("FMD%d", i));
622       ltx->SetNDC();
623       ltx->SetTextAngle(90);
624       ltx->SetTextColor(kGreen+4);
625       old = ll->FindObject(ltx->GetName());
626       if (old) { 
627         ll->Remove(old);
628         delete old;
629       }
630       ll->Add(ltx);
631
632       if (mean > fROErrorsBad) { 
633         AliWarningF("Mean of readout errors for FMD%d = %f > %f (%f)", 
634                     i, mean, fROErrorsBad, fROErrorsFkup);
635         ret++;
636         ltx->SetTextColor(kOrange+2);
637         if (mean > fROErrorsFkup) {
638           ret++;
639           ltx->SetTextColor(kRed+2);
640         }
641       }
642     }
643   }
644   else if (name.Contains("eloss",TString::kIgnoreCase)) { 
645     // If we' asked to not fit the data, return immediately
646     if (fNoFits) return ret;
647     // Do not fit cosmic or calibration data 
648     if (specie == AliRecoParam::kCosmic || 
649         specie == AliRecoParam::kCalib) return ret;
650     // Do not fit `expert' histograms 
651     if (hist->TestBit(AliQAv1::GetExpertBit())) return ret;
652     // Do not fit histograms with too little data 
653     if (hist->GetEntries() < fELossMinEntries) return ret;
654
655     // Try to fit a function to the histogram 
656     Double_t xMin  = hist->GetXaxis()->GetXmin();
657     Double_t xMax  = hist->GetXaxis()->GetXmax();
658
659     hist->GetXaxis()->SetRangeUser(fELossLowCut, xMax);
660     Int_t    bMaxY = hist->GetMaximumBin();
661     Double_t xMaxY = hist->GetXaxis()->GetBinCenter(bMaxY);
662     Double_t rms   = hist->GetRMS();
663     Double_t low   = hist->GetXaxis()->GetBinCenter(bMaxY-4);
664     hist->GetXaxis()->SetRangeUser(0.2, xMaxY+(fELossNRMS+1)*rms);
665     rms  = hist->GetRMS();
666     hist->GetXaxis()->SetRange(0,-1);
667     TF1*          func = makeLandauGaus(name);
668     func->SetParameter(1, xMaxY);
669     func->SetLineColor(kGreen+4);
670     // func->SetLineStyle(2);
671     Double_t high = xMax; // xMaxY+fELossNRMS*rms;
672     if (fELossNRMS > 0) high = xMaxY+fELossNRMS*rms;
673     
674     // Check we don't have an empty fit range 
675     if (low >= high) return ret;
676
677     // Check that we have enough counts in the fit range 
678     Int_t bLow  = hist->FindBin(low);
679     Int_t bHigh = hist->FindBin(high);
680     if (bLow >= bHigh || hist->Integral(bLow, bHigh) < fELossMinEntries)
681       return ret;
682
683     // Set our fit function 
684     TString fitOpt("QS");
685     TFitResultPtr res   = hist->Fit(func, fitOpt, "", low, high);
686     Int_t         color = func->GetLineColor();
687     UShort_t      qual  = CheckFit(hist, res, low, high, color);
688
689     // Make sure we save the function in the full range of the histogram
690     func = hist->GetFunction("landauGaus");
691     if (fELossNRMS <= 0) func->SetRange(xMin, xMax);
692     // func->SetParent(hist);
693     func->Save(xMin, xMax, 0, 0, 0, 0);
694     func->SetLineColor(color);
695
696     fitOpt.Append("+");
697     res = hist->Fit("pol2", fitOpt, "", fELossMinSharing, low-0.05);
698     func = hist->GetFunction("pol2");
699     Double_t   s     = Chi2Scale(hist,fELossMinEntries*100);
700     Double_t   chi2  = (!res.Get() ? 0 : res->Chi2());
701     Int_t      nu    = (!res.Get() ? 1 : res->Ndf());
702     Double_t   red   = s * (nu == 0 ? fELossFkupChi2Nu : chi2 / nu);
703     if (AliDebugLevel()) 
704       printf("FIT: %s, 2, %d, %f, %f\n", hist->GetName(),
705              Int_t(hist->GetEntries()), red, s * red);
706     red *= s;
707     if (red > fELossFkupChi2Nu) func->SetLineColor(kRed);
708     else                        func->SetLineColor(kGreen+4);
709
710     // Now check if this histogram should be cleared or not 
711     if (fELossMaxEntries > 0 && hist->GetEntries() > fELossMaxEntries)
712       hist->SetBit(AliFMDQADataMakerRec::kResetBit);
713     if (qual > 0) { 
714       func->SetLineWidth(3);
715       func->SetLineStyle(1);
716       if (qual > 1) 
717         func->SetLineWidth(4);
718     }
719     ret += qual;
720   }
721
722   return ret;
723 }
724 //__________________________________________________________________
725 UShort_t
726 AliFMDQAChecker::CheckSim(AliRecoParam::EventSpecie_t /* specie*/, 
727                           TH1*                        hist) const
728 {
729   // 
730   // Check simulated hits 
731   // 
732   return BasicCheck(hist);
733 }
734 //__________________________________________________________________
735 UShort_t
736 AliFMDQAChecker::CheckRec(AliRecoParam::EventSpecie_t /* specie*/, 
737                           TH1*                        hist) const
738 {
739   // 
740   // Check reconstructed data 
741   // 
742   return BasicCheck(hist);
743 }
744
745 //__________________________________________________________________
746 void AliFMDQAChecker::AddStatusPave(TH1* hist, Int_t qual, 
747                                     Double_t xl, Double_t yl, 
748                                     Double_t xh, Double_t yh) const
749 {
750   //
751   // Add a status pave to a plot
752   // 
753   if (xh < 0) xh = gStyle->GetStatX();
754   if (xl < 0) xl = xh - gStyle->GetStatW(); 
755   if (yh < 0) yh = gStyle->GetStatY();
756   if (yl < 0) yl = xl - gStyle->GetStatH(); 
757   
758   TPaveText* text = new TPaveText(xl, yl, xh, yh, "brNDC");
759   Int_t   bg  = kGreen-10;
760   Int_t   fg  = kBlack;
761   TString msg = "OK";
762   if      (qual >= kWhatTheFk) { bg = kRed+1; fg = kWhite; msg = "Argh!"; }
763   else if (qual >= kBad)       { bg = kRed-3; fg = kWhite; msg = "Bad"; }
764   else if (qual >= kProblem)   { bg = kOrange-4; msg = "Warning"; }
765   text->AddText(msg);
766   text->SetTextFont(62);
767   text->SetTextColor(fg);
768   text->SetFillColor(bg);
769
770   TList*   ll  = hist->GetListOfFunctions();
771   TObject* old = ll->FindObject(text->GetName());
772   if (old) { 
773     ll->Remove(old);
774     delete old;
775   }
776   ll->Add(text);
777 }
778
779 //__________________________________________________________________
780 void AliFMDQAChecker::Check(Double_t*                   rv, 
781                             AliQAv1::ALITASK_t          what, 
782                             TObjArray**                 list, 
783                             const AliDetectorRecoParam* /*t*/) 
784 {
785   // 
786   // Member function called to do the actual checking
787   //
788   // Parameters: 
789   //    rv   Array of return values. 
790   //    what What to check 
791   //    list Array of arrays of histograms.  There's one arrat for
792   //         each 'specie'
793   //    t    Reconstruction parameters - not used. 
794   //
795   // The bounds defined for RV are 
796   // 
797   //   FATAL:      [-1,   0.0]
798   //   ERROR:      (0.0,  0.002]
799   //   WARNING:    (0.002,0.5]
800   //   INFO:       (0.5,  1.0]
801   //
802   // Double_t* rv = new Double_t[AliRecoParam::kNSpecies] ; 
803
804   if (!fDidExternal) ProcessExternalParams();
805
806   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
807     // Int_t count   = 0;
808     rv[specie]    = 0.; 
809
810     if (!AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
811       continue ;
812     
813     if(!list[specie]) continue;
814     
815     TH1*     hist  = 0;
816     Int_t    nHist = list[specie]->GetEntriesFast();
817
818     // Find the status histogram if any 
819     TH2*  status  = 0;
820     Int_t istatus = AliFMDQADataMakerRec::GetHalfringIndex(4, 'i', 0, 0);
821     if (istatus < nHist) 
822       status = dynamic_cast<TH2*>(list[specie]->At(istatus));
823       
824     UShort_t ret   = 0;
825     for(Int_t i= 0; i< nHist; i++) {
826       if (!(hist = static_cast<TH1*>(list[specie]->At(i)))) continue;
827       if (hist == status) continue;
828       
829       Int_t qual = CheckOne(what, AliRecoParam::ConvertIndex(specie), hist);
830       hist->SetUniqueID(Quality2Bit(qual));
831       hist->SetStats(0);
832       AddStatusPave(hist, qual);
833       ret += qual;
834
835       if (!status) continue;
836
837       // Parse out the detector and ring, calculate the bin, and fill
838       // status histogram.
839       TString nme(hist->GetName());
840       Char_t cD   = nme[nme.Length()-2];
841       Char_t cR   = nme[nme.Length()-1];
842       Int_t  xbin = 0;
843       switch (cD) { 
844       case '1': xbin = 1; break;
845       case '2': xbin = 2 + ((cR == 'i' || cR == 'I') ? 0 : 1); break;
846       case '3': xbin = 4 + ((cR == 'i' || cR == 'I') ? 0 : 1); break;
847       }
848       if (xbin == 0) continue;
849       status->Fill(xbin, qual);
850                    
851     } // for (int i ...)
852     rv[specie] = ret;
853     // if      (ret > kWhatTheFk) rv[specie] = fLowTestValue[AliQAv1::kFATAL];
854     // else if (ret > kBad)       rv[specie] = fUpTestValue[AliQAv1::kERROR]; 
855     // else if (ret > kProblem)   rv[specie] = fUpTestValue[AliQAv1::kWARNING]; 
856     // else                       rv[specie] = fUpTestValue[AliQAv1::kINFO]; 
857     AliDebugF(3, "Combined sum is %d -> %f", ret, rv[specie]);
858
859     if (status) { 
860       Int_t nProblem = 0;
861       Int_t nBad     = 0;
862       for (Int_t i = 1; i < status->GetXaxis()->GetNbins(); i++) { 
863         nProblem += status->GetBinContent(i, 3);
864         nBad     += status->GetBinContent(i, 4);
865       }
866       Int_t qual = 0;
867       if (nProblem > fMaxNProblem) qual++;
868       if (nBad     > fMaxNBad)     qual += 2;
869       status->SetUniqueID(Quality2Bit(qual));
870       AddStatusPave(status, qual);
871     }
872     // if (count != 0) rv[specie] /= count;
873   }
874   // return rv;
875 }
876
877 namespace {
878   Int_t CheckForLog(TAxis*       axis,
879                     TVirtualPad* pad, 
880                     Int_t        xyz)
881   {
882     Int_t ret = 0;
883     TString t(axis->GetTitle());
884     if (!t.Contains("[log]", TString::kIgnoreCase)) return 0;
885     t.ReplaceAll("[log]", "");
886     switch (xyz) { 
887     case 1: pad->SetLogx(); ret |= 0x1; break;
888     case 2: pad->SetLogy(); ret |= 0x2; break;
889     case 3: pad->SetLogz(); ret |= 0x4; break;
890     }
891     axis->SetTitle(t);
892     return ret;
893   }
894   void RestoreLog(TAxis* axis, Bool_t log) 
895   {
896     if (!log) return;
897     TString t(axis->GetTitle());
898     t.Append("[log]");
899     axis->SetTitle(t);
900   }
901 }
902
903 namespace {
904   void FindMinMax(TH1* h, Double_t& min, Double_t& max)
905   {
906     Double_t tmin = 1e9;
907     Double_t tmax = 0;
908     for (Int_t i = 1; i <= h->GetNbinsX(); i++) { 
909       Double_t c = h->GetBinContent(i);
910       if (c < 1e-8) continue;
911       tmin = TMath::Min(tmin, c);
912       tmax = TMath::Max(tmax, c);
913     }
914     min = tmin;
915     max = tmax;
916   }
917 }
918
919 namespace { 
920   Int_t GetHalfringPad(TH1* h) {
921     TString nme(h->GetName());
922     Char_t cD   = nme[nme.Length()-2];
923     Char_t cR   = nme[nme.Length()-1];
924     Int_t  xbin = 0;
925     switch (cD) { 
926     case '1': xbin = 1; break;
927     case '2': xbin = ((cR == 'i' || cR == 'I') ? 2 : 5); break;
928     case '3': xbin = ((cR == 'i' || cR == 'I') ? 3 : 6); break;
929     }
930     return xbin;
931   }
932 }
933
934 //____________________________________________________________________________ 
935 void 
936 AliFMDQAChecker::MakeImage(TObjArray** list, 
937                            AliQAv1::TASKINDEX_t task, 
938                            AliQAv1::MODE_t mode) 
939 {
940   // makes the QA image for sim and rec
941   // 
942   // Parameters: 
943   //    task What to check 
944   //    list Array of arrays of histograms.  There's one array for
945   //         each 'specie'
946   //    t    Reconstruction parameters - not used. 
947   // 
948   Int_t    nImages = 0 ;
949   Double_t max     = 0;
950   Double_t min     = 10000;
951
952   // Loop over all species 
953   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
954     AliRecoParam::EventSpecie_t spe = AliRecoParam::ConvertIndex(specie);
955     if (!AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))
956         ->IsEventSpecieSet(spe)) 
957       continue;
958                                                                         
959     // If nothing is defined for this specie, go on. 
960     if(!list[specie] || list[specie]->GetEntriesFast() == 0) continue;
961
962     // Loop over the histograms and figure out how many histograms we
963     // have and the min/max 
964     TH1* hist  = 0;
965     Int_t nHist = list[specie]->GetEntriesFast();
966     for(Int_t i= 0; i< nHist; i++) {
967       hist = static_cast<TH1F*>(list[specie]->At(i));
968       if (hist && hist->TestBit(AliQAv1::GetImageBit())) {
969         nImages++; 
970         TString name(hist->GetName());
971         if (name.Contains("readouterrors", TString::kIgnoreCase) || 
972             name.Contains("status", TString::kIgnoreCase)) continue;
973
974         // Double_t hMax = hist->GetMaximum(); 
975         // hist->GetBinContent(hist->GetMaximumBin());
976         // Double_t hMin = hist->GetMinimum();
977         // hist->GetBinContent(hist->GetMinimumBin());
978         Double_t hMax, hMin;
979         FindMinMax(hist, hMin, hMax);
980         max = TMath::Max(max, hMax);
981         min = TMath::Min(min, hMin);
982         // AliInfoF("Min/max of %40s: %f/%f, global -> %f/%f", 
983         //          hist->GetName(), hMin, hMax, min, max);
984       }
985     }
986     break ; 
987   }
988   // AliInfoF("Global min/max=%f/%f", min, max);
989   min = TMath::Max(1e-1, min);
990   max = TMath::Max(1e5,  max);
991
992   // IF no images, go on. 
993   if (nImages == 0) {
994     AliDebug(AliQAv1::GetQADebugLevel(), 
995              Form("No histogram will be plotted for %s %s\n", GetName(), 
996                   AliQAv1::GetTaskName(task).Data()));
997     return;
998   }
999
1000   AliDebug(AliQAv1::GetQADebugLevel(), 
1001            Form("%d histograms will be plotted for %s %s\n", 
1002                 nImages, GetName(), AliQAv1::GetTaskName(task).Data()));  
1003   gStyle->SetOptStat(0);
1004   
1005   // Again loop over species and draw a canvas 
1006   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
1007     if (!AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))
1008         ->IsEventSpecieSet(specie)) continue;
1009
1010     // if Nothing here, go on
1011     if(!list[specie] || list[specie]->GetEntries() <= 0 || 
1012        nImages <= 0) continue;
1013
1014    // Form the title 
1015     const Char_t * title = Form("QA_%s_%s_%s", GetName(), 
1016                                 AliQAv1::GetTaskName(task).Data(), 
1017                                 AliRecoParam::GetEventSpecieName(specie)); 
1018     if (!fImage[specie]) fImage[specie] = new TCanvas(title, title) ;
1019     fImage[specie]->Clear() ; 
1020     fImage[specie]->SetTitle(title) ; 
1021     fImage[specie]->cd() ; 
1022
1023     // Put something in the canvas - even if empty 
1024     TPaveText someText(0.015, 0.015, 0.98, 0.98) ;
1025     someText.AddText(title) ;
1026     someText.SetFillColor(0);
1027     someText.SetFillStyle(0);
1028     someText.SetBorderSize(0);
1029     someText.SetTextColor(kRed+1);
1030     someText.Draw() ; 
1031     TString outName(Form("%s%s%d.%s", AliQAv1::GetImageFileName(), 
1032                          AliQAv1::GetModeName(mode), 
1033                          AliQAChecker::Instance()->GetRunNumber(), 
1034                          AliQAv1::GetImageFileFormat()));
1035     fImage[specie]->Print(outName, "ps") ; 
1036
1037     // Now set some parameters on the canvas 
1038     fImage[specie]->Clear(); 
1039     fImage[specie]->SetTopMargin(0.10);
1040     fImage[specie]->SetBottomMargin(0.15);
1041     fImage[specie]->SetLeftMargin(0.15);
1042     fImage[specie]->SetRightMargin(0.05);
1043     
1044     // Put title on top 
1045     const char* topT = Form("Mode: %s, Task: %s, Specie: %s, Run: %d",
1046                             AliQAv1::GetModeName(mode), 
1047                             AliQAv1::GetTaskName(task).Data(), 
1048                             AliRecoParam::GetEventSpecieName(specie),
1049                             AliQAChecker::Instance()->GetRunNumber());
1050     TLatex* topText = new TLatex(.5, .99, topT);
1051     topText->SetTextAlign(23);
1052     topText->SetTextSize(.038);
1053     topText->SetTextFont(42);
1054     topText->SetTextColor(kBlue+3);
1055     topText->SetNDC();
1056     topText->Draw();
1057
1058     // Find the status histogram if any 
1059     TH2*  status  = 0;
1060     Int_t istatus = AliFMDQADataMakerRec::GetHalfringIndex(4, 'i', 0, 0);
1061     if (istatus < list[specie]->GetEntriesFast()) 
1062       status = dynamic_cast<TH2*>(list[specie]->At(istatus));
1063
1064     // Divide canvas 
1065     // if (fDoScale) 
1066     TVirtualPad* plots = fImage[specie];
1067     TVirtualPad* stat  = 0;
1068     if (status) {
1069       // AliWarning("Drawing plots sub-pad");
1070       TPad* pM = new TPad("plots", "Plots Pad", 0, .2, 1., .9, 0, 0);
1071       fImage[specie]->cd();
1072       pM->Draw();
1073       plots = pM;
1074       // AliWarning("Drawing status sub-pad");
1075       TPad* pS = new TPad("status", "Status Pad", 0, 0, 1., .2, 0, 0);
1076       fImage[specie]->cd();
1077       pS->Draw();
1078       pS->SetLogz();
1079       stat = pS;
1080       // status->DrawCopy("colz");
1081     }
1082     // AliWarningF("fImage[specie]=%p, plots=%p", fImage[specie], plots);
1083     // plots->cd();
1084     Int_t nx = 3;
1085     Int_t ny = (nImages + .5) / nx;
1086     plots->Divide(nx, ny, 0, 0);
1087     // else fImage[specie]->Divide(nx, ny);
1088     
1089     
1090     // Loop over histograms 
1091     TH1*  hist  = 0;
1092     Int_t nHist = list[specie]->GetEntriesFast();
1093     for (Int_t i = 0; i < nHist; i++) { 
1094       hist = static_cast<TH1*>(list[specie]->At(i));
1095       if (!hist || !hist->TestBit(AliQAv1::GetImageBit())) continue;
1096       if (hist == status) continue;
1097       TString name(hist->GetName());
1098       Bool_t isROE = name.Contains("readouterrors", TString::kIgnoreCase);
1099
1100       // Go to sub-pad 
1101       TVirtualPad* pad = 0;
1102       if      (isROE) pad = plots->cd(4);
1103       else            pad = plots->cd(GetHalfringPad(hist));
1104       
1105       pad->SetRightMargin(0.01);
1106       if (!fDoScale) { 
1107         pad->SetLeftMargin(0.10);
1108         pad->SetBottomMargin(0.10);
1109       }
1110
1111       // Check for log scale 
1112       Int_t logOpts = 0;
1113       logOpts |= CheckForLog(hist->GetXaxis(), pad, 1);
1114       logOpts |= CheckForLog(hist->GetYaxis(), pad, 2);
1115       logOpts |= CheckForLog(hist->GetZaxis(), pad, 3);
1116
1117       // Figure out special cases 
1118       TString opt("");
1119       if (isROE) {
1120         pad->SetRightMargin(0.15);
1121         pad->SetBottomMargin(0.10);
1122         // pad->SetTopMargin(0.02);
1123         opt="COLZ";
1124       }
1125       else {
1126         pad->SetGridx();
1127         pad->SetGridy();
1128         if (fDoScale) { 
1129           hist->SetMinimum(min);
1130           hist->SetMaximum(max);
1131         }
1132         else { 
1133           hist->SetMinimum();
1134           hist->SetMaximum();
1135         }
1136       }
1137       // Draw (As a copy)
1138       hist->DrawCopy(opt);
1139       
1140       // Special cases 
1141       if (!name.Contains("readouterrors", TString::kIgnoreCase)) {
1142         gStyle->SetOptTitle(0);
1143         TPad* insert = new TPad("insert", "Zoom", 
1144                                 .5,.5, .99, .95, 0, 0, 0);
1145         insert->SetTopMargin(0.01);
1146         insert->SetRightMargin(0.01);
1147         insert->SetFillColor(0);
1148         insert->SetBorderSize(1);
1149         insert->SetBorderMode(0);
1150         insert->Draw();
1151         insert->cd();
1152         if (logOpts & 0x1) insert->SetLogx();
1153         if (logOpts & 0x2) insert->SetLogy();
1154         if (logOpts & 0x4) insert->SetLogz();
1155         hist->GetXaxis()->SetRange(1, hist->GetNbinsX()/8);
1156         TH1* copy = hist->DrawCopy(opt);
1157         copy->GetXaxis()->SetNdivisions(408, false);
1158         // Restore full range 
1159         hist->GetXaxis()->SetRange(0, 0);
1160         gStyle->SetOptTitle(1);
1161       }
1162       pad->cd();
1163       // Possibly restore the log options 
1164       RestoreLog(hist->GetXaxis(), logOpts & 0x1);
1165       RestoreLog(hist->GetYaxis(), logOpts & 0x2);
1166       RestoreLog(hist->GetZaxis(), logOpts & 0x4);
1167     }
1168     if (status && stat) {
1169       stat->cd();
1170       status->DrawCopy("BOX TEXT");
1171     }
1172     // Print to a post-script file 
1173     fImage[specie]->Print(outName, "ps");
1174     if (AliDebugLevel() > 0) 
1175       fImage[specie]->Print(Form("%s_%d.png", 
1176                                  AliRecoParam::GetEventSpecieName(specie), 
1177                                  AliQAChecker::Instance()->GetRunNumber()));
1178   }
1179 }
1180
1181 //__________________________________________________________________
1182 //
1183 // EOF
1184 //