AliRoot Core  edcc906 (edcc906)
AliFMDQAChecker.cxx
Go to the documentation of this file.
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 //__________________________________________________________________
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  fMinMPV(.2),
267  fMaxXi(.3),
268  fMaxSigma(.6),
269  fNoFits(false)
270 {
271 }
272 
273 //__________________________________________________________________
274 void
276 {
277  ProcessExternalParam("ELossLowCut", fELossLowCut);
278  ProcessExternalParam("ELossNRMS", fELossNRMS);
279  ProcessExternalParam("ELossBadChi2Nu", fELossBadChi2Nu);
280  ProcessExternalParam("ELossFkupChi2Nu", fELossFkupChi2Nu);
281  ProcessExternalParam("ELossGoodParError", fELossGoodParError);
282  ProcessExternalParam("ROErrorsBad", fROErrorsBad);
283  ProcessExternalParam("ROErrorsFkup", fROErrorsFkup);
284  ProcessExternalParam("ELossMinSharing", fELossMinSharing);
285  Double_t tmp = 0;
286  ProcessExternalParam("CommonScale", tmp);
287  fDoScale = tmp > 0; tmp = fELossMinEntries;
288  ProcessExternalParam("ELossMinEntries", tmp);
289  fELossMinEntries = tmp; tmp = fELossMaxEntries;
290  ProcessExternalParam("ELossMaxEntries", tmp);
291  fELossMaxEntries = tmp; tmp = fMaxNProblem;
292  ProcessExternalParam("MaxNProblem", tmp);
293  fMaxNProblem = tmp; tmp = 0;
294  fELossMaxEntries = tmp; tmp = fMaxNBad;
295  ProcessExternalParam("MaxNBad", tmp);
296  fMaxNBad = tmp; tmp = 0;
297  ProcessExternalParam("NoFits", tmp);
298  fNoFits = tmp > 0; tmp = 0;
299 
300  ProcessExternalParam("ELossMinMPV", fMinMPV);
301  ProcessExternalParam("ELossMaxXi", fMaxXi);
302  ProcessExternalParam("ELossMaxSigma", fMaxSigma);
303 
304  GetThresholds();
305 
306  fDidExternal = true;
307 }
308 //__________________________________________________________________
309 void
310 AliFMDQAChecker::ProcessExternalParam(const char* name, Double_t& v)
311 {
312  TObject* o = fExternParamList->FindObject(name);
313  if (!o) return;
314  TParameter<double>* p = static_cast<TParameter<double>*>(o);
315  v = p->GetVal();
316  AliDebugF(3, "External parameter: %-20s=%lf", name, v);
317 }
318 
319 //__________________________________________________________________
320 void
322 {
323  const char* path = "GRP/Calib/QAThresholds";
325  AliCDBEntry* cdbEnt = cdbMan->Get(path);
326  if (!cdbEnt) {
327  AliWarningF("Failed to get CDB entry at %s", path);
328  return;
329  }
330 
331  TObjArray* cdbObj = static_cast<TObjArray*>(cdbEnt->GetObject());
332  if (!cdbObj) {
333  AliWarningF("Failed to get CDB object at %s", path);
334  return;
335  }
336 
337  TObject* fmdObj = cdbObj->FindObject("FMD");
338  if (!fmdObj) {
339  AliWarningF("Failed to get FMD object at from CDB %s", path);
340  return;
341  }
342 
343  AliQAThresholds* qaThr = static_cast<AliQAThresholds*>(fmdObj);
344  Int_t nThr = qaThr->GetSize();
345  for (Int_t i = 0; i < nThr; i++) {
346  TObject* thr = qaThr->GetThreshold(i);
347  if (!thr) continue;
348 
349  TParameter<double>* d = dynamic_cast<TParameter<double>*>(thr);
350  if (!d) {
351  AliWarningF("Parameter %s not of type double", thr->GetName());
352  continue;
353  }
354  Double_t val = d->GetVal();
355  TString name(thr->GetName());
356  if (name.EqualTo("ELossBadChi2Nu")) fELossBadChi2Nu = val;
357  else if (name.EqualTo("ELossFkupChi2Nu")) fELossFkupChi2Nu = val;
358  else if (name.EqualTo("ELossGoodParError")) fELossGoodParError = val;
359  else if (name.EqualTo("ROErrorsBad")) fROErrorsBad = val;
360  else if (name.EqualTo("ROErrorsFkup")) fROErrorsFkup = val;
361  else if (name.EqualTo("MaxNProblem")) fMaxNProblem = val;
362  else if (name.EqualTo("MaxNBad")) fMaxNBad = val;
363  else if (name.EqualTo("ELossMinMPV")) fMinMPV = val;
364  else if (name.EqualTo("ELossMaxXi")) fMaxXi = val;
365  else if (name.EqualTo("ELossMaxSigma")) fMaxSigma = val;
366 
367  AliDebugF(3, "Threshold %s=%f", name.Data(), val);
368  }
369 }
370 
371 //__________________________________________________________________
373 AliFMDQAChecker::Quality2Bit(UShort_t value) const
374 {
375  AliQAv1::QABIT_t ret = AliQAv1::kINFO; // Assume success
376  if (value >= kWhatTheFk) ret = AliQAv1::kFATAL;
377  else if (value >= kBad) ret = AliQAv1::kERROR;
378  else if (value >= kProblem) ret = AliQAv1::kWARNING;
379 
380  return ret;
381 }
382 
383 //__________________________________________________________________
384 void
385 AliFMDQAChecker::SetQA(AliQAv1::ALITASK_t index, Double_t* values) const
386 {
387  AliQAv1 * qa = AliQAv1::Instance(index) ;
388 
389  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
390  // Check if specie is defined
392  continue ;
393 
394  // No checker is implemented, set all QA to Fatal
395  if (!values) {
396  qa->Set(AliQAv1::kFATAL, specie) ;
397  continue;
398  }
399 
400  UShort_t value = values[specie];
401  AliQAv1::QABIT_t ret = Quality2Bit(value);
402 
403  qa->Set(ret, AliRecoParam::ConvertIndex(specie));
404  AliDebugF(3, "Quality of %s: %d -> %d",
405  AliRecoParam::GetEventSpecieName(specie), value, ret);
406  }
407 }
408 
409 //__________________________________________________________________
410 UShort_t
412 {
413  if (hist->GetEntries() <= 0) return kOK;
414  return (hist->GetMean() > 0 ? kOK : kProblem);
415 }
416 
417 //__________________________________________________________________
418 UShort_t
421  TH1* hist) const
422 {
423  if(what == AliQAv1::kESD) return CheckESD(specie, hist);
424  if(what == AliQAv1::kRAW) return CheckRaw(specie, hist);
425  if(what == AliQAv1::kSIM) return CheckSim(specie, hist);
426  if(what == AliQAv1::kREC) return CheckRec(specie, hist);
427  return 0;
428 }
429 //__________________________________________________________________
430 UShort_t
432  TH1* hist) const
433 {
434  return BasicCheck(hist);
435 }
436 namespace {
437  Double_t Chi2Scale(TH1* h, Double_t base=10000)
438  {
439  return 1. / TMath::Max(1., h->GetEntries() / base);
440  }
441  void AddLine(TObjArray* lines,
442  Double_t x1, Double_t x2, Double_t x3,
443  Double_t y, Double_t dy,
444  const char* name, Double_t val, Double_t lim,
445  Bool_t ok, Int_t color)
446  {
447  TString n; n.Form("%s:", name);
448  TLatex* ltx = new TLatex(x1, y, n);
449  ltx->SetNDC(true);
450  ltx->SetTextSize(dy-0.01);
451  ltx->SetTextColor(color);
452  lines->Add(ltx);
453 
454  n.Form("%7.3f", val);
455  ltx = new TLatex(x2, y, n);
456  ltx->SetNDC(true);
457  ltx->SetTextSize(dy-0.01);
458  ltx->SetTextColor(color);
459  lines->Add(ltx);
460 
461  if (lim < 0) n = "(ignored)";
462  else n.Form("%c %4.2f", ok ? '<' : '>', lim);
463  ltx = new TLatex(x3, y, n);
464  ltx->SetNDC(true);
465  ltx->SetTextSize(dy-0.01);
466  ltx->SetTextColor(color);
467  lines->Add(ltx);
468  }
469 }
470 
471 //__________________________________________________________________
472 UShort_t
473 AliFMDQAChecker::CheckFit(TH1* hist, const TFitResultPtr& res,
474  Double_t low, Double_t high, Int_t& color) const
475 {
476  color = kGreen+4;
477 
478  // Check if there's indeed a result - if not, flag as OK
479  if (!res.Get()) return 0;
480 
481  UShort_t ret = 0;
482  Int_t nPar = res->NPar();
483  Double_t dy = .06;
484  Double_t x = .2;
485  Double_t x2 = .3;
486  Double_t x3 = .4;
487  Double_t y = .9-dy;
488  Double_t chi2 = res->Chi2();
489  Int_t nu = res->Ndf();
490  Double_t s = Chi2Scale(hist,fELossMinEntries);
491  Double_t red = (nu == 0 ? fELossFkupChi2Nu : chi2 / nu);
492  TObjArray* lines = 0;
493  // TLatex* lRed = 0;
494  TLatex* ltx = 0;
495  Int_t chi2Check = 0;
496  Double_t chi2Lim = fELossBadChi2Nu;
497  if (AliDebugLevel() > 0)
498  printf("FIT: %s, 1, %d, %f, %f\n", hist->GetName(),
499  Int_t(hist->GetEntries()), red, s * red);
500  red *= s;
501  if (red > fELossBadChi2Nu) { // || res->Prob() < .01) {
502  // AliWarningF("Fit gave chi^2/nu=%f/%d=%f>%f (%f)",
503  // res->Chi2(), res->Ndf(), red, fELossBadChi2Nu,
504  // fELossFkupChi2Nu);
505  // res->Print();
506  chi2Check++;
507  if (red > fELossFkupChi2Nu) {
508  chi2Check++;
509  chi2Lim = fELossFkupChi2Nu;
510  }
511  }
512  ret += chi2Check;
513 
514  if (fShowFitResults) {
515  lines = new TObjArray(nPar+3);
516  lines->SetName("lines");
517  lines->SetOwner(true);
518 
519  AddLine(lines, x, x2, x3, y, dy, "#chi^{2}/#nu", red, chi2Lim,
520  chi2Check < 1, chi2Check < 1 ? color :
521  chi2Check < 2 ? kOrange+2 : kRed+2);
522 
523  Double_t x1 = .85;
524  Double_t y1 = .5;
525 
526  // y1 -= dy;
527  ltx = new TLatex(x1, y1, Form("Fit range: [%6.2f,%6.2f]", low, high));
528  ltx->SetTextColor(kGray+3);
529  ltx->SetTextSize(dy-.01);
530  ltx->SetTextAlign(31);
531  ltx->SetNDC(true);
532  lines->Add(ltx);
533 
534  y1 -= dy;
535  ltx = new TLatex(x1, y1, Form("Entries: %d (%d)",
536  Int_t(hist->GetEffectiveEntries()),
538  ltx->SetTextColor(kGray+3);
539  ltx->SetTextSize(dy-.01);
540  ltx->SetTextAlign(31);
541  ltx->SetNDC(true);
542  lines->Add(ltx);
543 
544  y1 -= dy;
545  ltx = new TLatex(x1, y1, Form("%s: %f #pm %f",
546  res->ParName(1).c_str(),
547  res->Parameter(1),
548  res->ParError(1)));
549  ltx->SetTextColor(kGray+3);
550  ltx->SetTextSize(dy-.01);
551  ltx->SetTextAlign(31);
552  ltx->SetNDC(true);
553  lines->Add(ltx);
554  }
555 
556  // Now check the relative error on the fit parameters
557  Int_t parsOk = 0;
558  for (Int_t i = 0; i < nPar; i++) {
559  if (res->IsParameterFixed(i)) continue;
560  Double_t thr = fELossGoodParError;
561  Double_t pv = res->Parameter(i);
562  Double_t pe = res->ParError(i);
563  Double_t rel = (pv == 0 ? 100 : pe / pv);
564  Bool_t ok = (i == 3) || (rel < thr);
565  if (i == 1 && pv < fMinMPV) { // Low peak
566  ok = false; ret++; } // Double penalty
567  if (i == 2 && pv > fMaxXi) { // Large xi
568  ok = false; ret++; } // Double penalty
569  if (i == 3 && pv > fMaxSigma) { // Large sigma
570  ok = false; ret++; } // Double penalty
571  if (lines) {
572  y -= dy;
573  AddLine(lines, x, x2, x3, y, dy,Form("#delta%s/%s",
574  res->ParName(i).c_str(),
575  res->ParName(i).c_str()),
576  rel, (i == 3 ? -1 : thr), ok, ok ? color : kOrange+2);
577  }
578  if (i == 3) continue; // Skip sigma
579  if (ok) parsOk++;
580  }
581  if (parsOk > 0)
582  ret = TMath::Max(ret-(parsOk-1),0);
583  if (ret > 1) color = kRed+2;
584  if (ret > 0) color = kOrange+2;
585 
586  if (lines) {
587  TList* lf = hist->GetListOfFunctions();
588  TObject* old = lf->FindObject(lines->GetName());
589  if (old) {
590  lf->Remove(old);
591  delete old;
592  }
593  lf->Add(lines);
594  }
595  hist->SetStats(false);
596 
597  return ret;
598 }
599 
600 //__________________________________________________________________
601 UShort_t
603  TH1* hist) const
604 {
605  Int_t ret = BasicCheck(hist);
606  TString name(hist->GetName());
607  if (name.Contains("readouterrors", TString::kIgnoreCase)) {
608  // Check the mean number of errors per event
609  TH2* roErrors = static_cast<TH2*>(hist);
610  Int_t nY = roErrors->GetNbinsY();
611 
612  TLatex* ltx = new TLatex(.15, .9, Form("Thresholds: %5.2f,%5.2f",
614  ltx->SetName("thresholds");
615  ltx->SetTextColor(kGray+3);
616  ltx->SetNDC();
617 
618  TList* ll = hist->GetListOfFunctions();
619  TObject* old = ll->FindObject(ltx->GetName());
620  if (old) {
621  ll->Remove(old);
622  delete old;
623  }
624  ll->Add(ltx);
625 
626  for (Int_t i = 1; i <= 3; i++) {
627  Double_t sum = 0;
628  Int_t cnt = 0;
629  for (Int_t j = 1; j <= nY; j++) {
630  Int_t n = roErrors->GetBinContent(i, j);
631  sum += n * roErrors->GetYaxis()->GetBinCenter(j);
632  cnt += n;
633  }
634  Double_t mean = (cnt <= 0 ? 0 : sum / cnt);
635  Double_t x = ((i-.5) * (1-0.1-0.1) / 3 + 0.1);
636 
637  ltx = new TLatex(x, kROErrorsLabelY, Form("Mean: %6.3f", mean));
638  ltx->SetName(Form("FMD%d", i));
639  ltx->SetNDC();
640  ltx->SetTextAngle(90);
641  ltx->SetTextColor(kGreen+4);
642  old = ll->FindObject(ltx->GetName());
643  if (old) {
644  ll->Remove(old);
645  delete old;
646  }
647  ll->Add(ltx);
648 
649  if (mean > fROErrorsBad) {
650  AliWarningF("Mean of readout errors for FMD%d = %f > %f (%f)",
651  i, mean, fROErrorsBad, fROErrorsFkup);
652  ret++;
653  ltx->SetTextColor(kOrange+2);
654  if (mean > fROErrorsFkup) {
655  ret++;
656  ltx->SetTextColor(kRed+2);
657  }
658  }
659  }
660  }
661  else if (name.Contains("eloss",TString::kIgnoreCase)) {
662  // If we' asked to not fit the data, return immediately
663  if (fNoFits) return ret;
664  // Do not fit cosmic or calibration data
665  if (specie == AliRecoParam::kCosmic ||
666  specie == AliRecoParam::kCalib) return ret;
667  // Do not fit `expert' histograms
668  if (hist->TestBit(AliQAv1::GetExpertBit())) return ret;
669  // Do not fit histograms with too little data
670  if (hist->GetEntries() < fELossMinEntries) return ret;
671 
672  // Try to fit a function to the histogram
673  Double_t xMin = hist->GetXaxis()->GetXmin();
674  Double_t xMax = hist->GetXaxis()->GetXmax();
675 
676  hist->GetXaxis()->SetRangeUser(fELossLowCut, xMax);
677  Int_t bMaxY = hist->GetMaximumBin();
678  Double_t xMaxY = hist->GetXaxis()->GetBinCenter(bMaxY);
679  Double_t rms = hist->GetRMS();
680  Double_t low = hist->GetXaxis()->GetBinCenter(bMaxY-4);
681  hist->GetXaxis()->SetRangeUser(0.2, xMaxY+(fELossNRMS+1)*rms);
682  rms = hist->GetRMS();
683  hist->GetXaxis()->SetRange(0,-1);
684  TF1* func = makeLandauGaus(name);
685  func->SetParameter(1, xMaxY);
686  func->SetLineColor(kGreen+4);
687  // func->SetLineStyle(2);
688  Double_t high = xMax; // xMaxY+fELossNRMS*rms;
689  if (fELossNRMS > 0) high = xMaxY+fELossNRMS*rms;
690 
691  // Check we don't have an empty fit range
692  if (low >= high) return ret;
693 
694  // Check that we have enough counts in the fit range
695  Int_t bLow = hist->FindBin(low);
696  Int_t bHigh = hist->FindBin(high);
697  if (bLow >= bHigh || hist->Integral(bLow, bHigh) < fELossMinEntries)
698  return ret;
699 
700  // Set our fit function
701  TString fitOpt("QS");
702  TFitResultPtr res = hist->Fit(func, fitOpt, "", low, high);
703  Int_t color = func->GetLineColor();
704  UShort_t qual = CheckFit(hist, res, low, high, color);
705 
706  // Make sure we save the function in the full range of the histogram
707  func = hist->GetFunction("landauGaus");
708  if (fELossNRMS <= 0) func->SetRange(xMin, xMax);
709  // func->SetParent(hist);
710  func->Save(xMin, xMax, 0, 0, 0, 0);
711  func->SetLineColor(color);
712 
713  fitOpt.Append("+");
714  res = hist->Fit("pol2", fitOpt, "", fELossMinSharing, low-0.05);
715  func = hist->GetFunction("pol2");
716  Double_t s = Chi2Scale(hist,fELossMinEntries*100);
717  Double_t chi2 = (!res.Get() ? 0 : res->Chi2());
718  Int_t nu = (!res.Get() ? 1 : res->Ndf());
719  Double_t red = s * (nu == 0 ? fELossFkupChi2Nu : chi2 / nu);
720  if (AliDebugLevel())
721  printf("FIT: %s, 2, %d, %f, %f\n", hist->GetName(),
722  Int_t(hist->GetEntries()), red, s * red);
723  red *= s;
724  if (red > fELossFkupChi2Nu) func->SetLineColor(kRed);
725  else func->SetLineColor(kGreen+4);
726 
727  // Now check if this histogram should be cleared or not
728  if (fELossMaxEntries > 0 && hist->GetEntries() > fELossMaxEntries)
729  hist->SetBit(AliFMDQADataMakerRec::kResetBit);
730  if (qual > 0) {
731  func->SetLineWidth(3);
732  func->SetLineStyle(1);
733  if (qual > 1)
734  func->SetLineWidth(4);
735  }
736  ret += qual;
737  }
738 
739  return ret;
740 }
741 //__________________________________________________________________
742 UShort_t
744  TH1* hist) const
745 {
746  //
747  // Check simulated hits
748  //
749  return BasicCheck(hist);
750 }
751 //__________________________________________________________________
752 UShort_t
754  TH1* hist) const
755 {
756  //
757  // Check reconstructed data
758  //
759  return BasicCheck(hist);
760 }
761 
762 //__________________________________________________________________
763 void AliFMDQAChecker::AddStatusPave(TH1* hist, Int_t qual,
764  Double_t xl, Double_t yl,
765  Double_t xh, Double_t yh) const
766 {
767  //
768  // Add a status pave to a plot
769  //
770  if (xh < 0) xh = gStyle->GetStatX();
771  if (xl < 0) xl = xh - gStyle->GetStatW();
772  if (yh < 0) yh = gStyle->GetStatY();
773  if (yl < 0) yl = xl - gStyle->GetStatH();
774 
775  TPaveText* text = new TPaveText(xl, yl, xh, yh, "brNDC");
776  Int_t bg = kGreen-10;
777  Int_t fg = kBlack;
778  TString msg = "OK";
779  if (qual >= kWhatTheFk) { bg = kRed+1; fg = kWhite; msg = "Argh!"; }
780  else if (qual >= kBad) { bg = kRed-3; fg = kWhite; msg = "Bad"; }
781  else if (qual >= kProblem) { bg = kOrange-4; msg = "Warning"; }
782  text->AddText(msg);
783  text->SetTextFont(62);
784  text->SetTextColor(fg);
785  text->SetFillColor(bg);
786 
787  TList* ll = hist->GetListOfFunctions();
788  TObject* old = ll->FindObject(text->GetName());
789  if (old) {
790  ll->Remove(old);
791  delete old;
792  }
793  ll->Add(text);
794 }
795 
796 //__________________________________________________________________
797 void AliFMDQAChecker::Check(Double_t* rv,
798  AliQAv1::ALITASK_t what,
799  TObjArray** list,
800  const AliDetectorRecoParam* /*t*/)
801 {
802  //
803  // Member function called to do the actual checking
804  //
805  // Parameters:
806  // rv Array of return values.
807  // what What to check
808  // list Array of arrays of histograms. There's one arrat for
809  // each 'specie'
810  // t Reconstruction parameters - not used.
811  //
812  // The bounds defined for RV are
813  //
814  // FATAL: [-1, 0.0]
815  // ERROR: (0.0, 0.002]
816  // WARNING: (0.002,0.5]
817  // INFO: (0.5, 1.0]
818  //
819  // Double_t* rv = new Double_t[AliRecoParam::kNSpecies] ;
820 
822 
823  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
824  // Int_t count = 0;
825  rv[specie] = 0.;
826 
827  if (!AliQAv1::Instance()->IsEventSpecieSet(specie) )
828  continue ;
829 
830  if(!list[specie]) continue;
831 
832  TH1* hist = 0;
833  Int_t nHist = list[specie]->GetEntriesFast();
834 
835  // Find the status histogram if any
836  TH2* status = 0;
837  Int_t istatus = AliFMDQADataMakerRec::GetHalfringIndex(4, 'i', 0, 0);
838  if (istatus < nHist)
839  status = dynamic_cast<TH2*>(list[specie]->At(istatus));
840 
841  UShort_t ret = 0;
842  for(Int_t i= 0; i< nHist; i++) {
843  if (!(hist = static_cast<TH1*>(list[specie]->At(i)))) continue;
844  if (hist == status) continue;
845 
846  Int_t qual = CheckOne(what, AliRecoParam::ConvertIndex(specie), hist);
847  hist->SetUniqueID(Quality2Bit(qual));
848  hist->SetStats(0);
849  AddStatusPave(hist, qual);
850  ret += qual;
851 
852  if (!status) continue;
853 
854  // Parse out the detector and ring, calculate the bin, and fill
855  // status histogram.
856  TString nme(hist->GetName());
857  Char_t cD = nme[nme.Length()-2];
858  Char_t cR = nme[nme.Length()-1];
859  Int_t xbin = 0;
860  switch (cD) {
861  case '1': xbin = 1; break;
862  case '2': xbin = 2 + ((cR == 'i' || cR == 'I') ? 0 : 1); break;
863  case '3': xbin = 4 + ((cR == 'i' || cR == 'I') ? 0 : 1); break;
864  }
865  if (xbin == 0) continue;
866  status->Fill(xbin, qual);
867 
868  } // for (int i ...)
869  rv[specie] = ret;
870  // if (ret > kWhatTheFk) rv[specie] = fLowTestValue[AliQAv1::kFATAL];
871  // else if (ret > kBad) rv[specie] = fUpTestValue[AliQAv1::kERROR];
872  // else if (ret > kProblem) rv[specie] = fUpTestValue[AliQAv1::kWARNING];
873  // else rv[specie] = fUpTestValue[AliQAv1::kINFO];
874  AliDebugF(3, "Combined sum is %d -> %f", ret, rv[specie]);
875 
876  if (status) {
877  Int_t nProblem = 0;
878  Int_t nBad = 0;
879  for (Int_t i = 1; i < status->GetXaxis()->GetNbins(); i++) {
880  nProblem += status->GetBinContent(i, 3);
881  nBad += status->GetBinContent(i, 4);
882  }
883  Int_t qual = 0;
884  if (nProblem > fMaxNProblem) qual++;
885  if (nBad > fMaxNBad) qual += 2;
886  status->SetUniqueID(Quality2Bit(qual));
887  AddStatusPave(status, qual);
888  }
889  // if (count != 0) rv[specie] /= count;
890  }
891  // return rv;
892 }
893 
894 namespace {
895  Int_t CheckForLog(TAxis* axis,
896  TVirtualPad* pad,
897  Int_t xyz)
898  {
899  Int_t ret = 0;
900  TString t(axis->GetTitle());
901  if (!t.Contains("[log]", TString::kIgnoreCase)) return 0;
902  t.ReplaceAll("[log]", "");
903  switch (xyz) {
904  case 1: pad->SetLogx(); ret |= 0x1; break;
905  case 2: pad->SetLogy(); ret |= 0x2; break;
906  case 3: pad->SetLogz(); ret |= 0x4; break;
907  }
908  axis->SetTitle(t);
909  return ret;
910  }
911  void RestoreLog(TAxis* axis, Bool_t log)
912  {
913  if (!log) return;
914  TString t(axis->GetTitle());
915  t.Append("[log]");
916  axis->SetTitle(t);
917  }
918 }
919 
920 namespace {
921  void FindMinMax(TH1* h, Double_t& min, Double_t& max)
922  {
923  Double_t tmin = 1e9;
924  Double_t tmax = 0;
925  for (Int_t i = 1; i <= h->GetNbinsX(); i++) {
926  Double_t c = h->GetBinContent(i);
927  if (c < 1e-8) continue;
928  tmin = TMath::Min(tmin, c);
929  tmax = TMath::Max(tmax, c);
930  }
931  min = tmin;
932  max = tmax;
933  }
934 }
935 
936 namespace {
937  Int_t GetHalfringPad(TH1* h) {
938  TString nme(h->GetName());
939  Char_t cD = nme[nme.Length()-2];
940  Char_t cR = nme[nme.Length()-1];
941  Int_t xbin = 0;
942  switch (cD) {
943  case '1': xbin = 1; break;
944  case '2': xbin = ((cR == 'i' || cR == 'I') ? 2 : 5); break;
945  case '3': xbin = ((cR == 'i' || cR == 'I') ? 3 : 6); break;
946  }
947  return xbin;
948  }
949 }
950 
951 //____________________________________________________________________________
952 void
954  AliQAv1::TASKINDEX_t task,
955  AliQAv1::MODE_t mode)
956 {
957  // makes the QA image for sim and rec
958  //
959  // Parameters:
960  // task What to check
961  // list Array of arrays of histograms. There's one array for
962  // each 'specie'
963  // t Reconstruction parameters - not used.
964  //
965  Int_t nImages = 0 ;
966  Double_t max = 0;
967  Double_t min = 10000;
968 
969  // Loop over all species
970  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
972  if (!AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))
973  ->IsEventSpecieSet(spe))
974  continue;
975 
976  // If nothing is defined for this specie, go on.
977  if(!list[specie] || list[specie]->GetEntriesFast() == 0) continue;
978 
979  // Loop over the histograms and figure out how many histograms we
980  // have and the min/max
981  TH1* hist = 0;
982  Int_t nHist = list[specie]->GetEntriesFast();
983  for(Int_t i= 0; i< nHist; i++) {
984  hist = static_cast<TH1F*>(list[specie]->At(i));
985  if (hist && hist->TestBit(AliQAv1::GetImageBit())) {
986  nImages++;
987  TString name(hist->GetName());
988  if (name.Contains("readouterrors", TString::kIgnoreCase) ||
989  name.Contains("status", TString::kIgnoreCase)) continue;
990 
991  // Double_t hMax = hist->GetMaximum();
992  // hist->GetBinContent(hist->GetMaximumBin());
993  // Double_t hMin = hist->GetMinimum();
994  // hist->GetBinContent(hist->GetMinimumBin());
995  Double_t hMax, hMin;
996  FindMinMax(hist, hMin, hMax);
997  max = TMath::Max(max, hMax);
998  min = TMath::Min(min, hMin);
999  // AliInfoF("Min/max of %40s: %f/%f, global -> %f/%f",
1000  // hist->GetName(), hMin, hMax, min, max);
1001  }
1002  }
1003  break ;
1004  }
1005  // AliInfoF("Global min/max=%f/%f", min, max);
1006  min = TMath::Max(1e-1, min);
1007  max = TMath::Max(1e5, max);
1008 
1009  // IF no images, go on.
1010  if (nImages == 0) {
1012  Form("No histogram will be plotted for %s %s\n", GetName(),
1013  AliQAv1::GetTaskName(task).Data()));
1014  return;
1015  }
1016 
1018  Form("%d histograms will be plotted for %s %s\n",
1019  nImages, GetName(), AliQAv1::GetTaskName(task).Data()));
1020  gStyle->SetOptStat(0);
1021 
1022  // Again loop over species and draw a canvas
1023  for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
1024  if (!AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))
1025  ->IsEventSpecieSet(specie)) continue;
1026 
1027  // if Nothing here, go on
1028  if(!list[specie] || list[specie]->GetEntries() <= 0 ||
1029  nImages <= 0) continue;
1030 
1031  // Form the title
1032  const Char_t * title = Form("QA_%s_%s_%s", GetName(),
1033  AliQAv1::GetTaskName(task).Data(),
1035  if (!fImage[specie]) fImage[specie] = new TCanvas(title, title) ;
1036  fImage[specie]->Clear() ;
1037  fImage[specie]->SetTitle(title) ;
1038  fImage[specie]->cd() ;
1039 
1040  // Put something in the canvas - even if empty
1041  TPaveText someText(0.015, 0.015, 0.98, 0.98) ;
1042  someText.AddText(title) ;
1043  someText.SetFillColor(0);
1044  someText.SetFillStyle(0);
1045  someText.SetBorderSize(0);
1046  someText.SetTextColor(kRed+1);
1047  someText.Draw() ;
1048  TString outName(Form("%s%s%d.%s", AliQAv1::GetImageFileName(),
1049  AliQAv1::GetModeName(mode),
1050  AliQAChecker::Instance()->GetRunNumber(),
1052  fImage[specie]->Print(outName, "ps") ;
1053 
1054  // Now set some parameters on the canvas
1055  fImage[specie]->Clear();
1056  fImage[specie]->SetTopMargin(0.10);
1057  fImage[specie]->SetBottomMargin(0.15);
1058  fImage[specie]->SetLeftMargin(0.15);
1059  fImage[specie]->SetRightMargin(0.05);
1060 
1061  // Put title on top
1062  const char* topT = Form("Mode: %s, Task: %s, Specie: %s, Run: %d",
1063  AliQAv1::GetModeName(mode),
1064  AliQAv1::GetTaskName(task).Data(),
1066  AliQAChecker::Instance()->GetRunNumber());
1067  TLatex* topText = new TLatex(.5, .99, topT);
1068  topText->SetTextAlign(23);
1069  topText->SetTextSize(.038);
1070  topText->SetTextFont(42);
1071  topText->SetTextColor(kBlue+3);
1072  topText->SetNDC();
1073  topText->Draw();
1074 
1075  // Find the status histogram if any
1076  TH2* status = 0;
1077  Int_t istatus = AliFMDQADataMakerRec::GetHalfringIndex(4, 'i', 0, 0);
1078  if (istatus < list[specie]->GetEntriesFast())
1079  status = dynamic_cast<TH2*>(list[specie]->At(istatus));
1080 
1081  // Divide canvas
1082  // if (fDoScale)
1083  TVirtualPad* plots = fImage[specie];
1084  TVirtualPad* stat = 0;
1085  if (status) {
1086  // AliWarning("Drawing plots sub-pad");
1087  TPad* pM = new TPad("plots", "Plots Pad", 0, .2, 1., .9, 0, 0);
1088  fImage[specie]->cd();
1089  pM->Draw();
1090  plots = pM;
1091  // AliWarning("Drawing status sub-pad");
1092  TPad* pS = new TPad("status", "Status Pad", 0, 0, 1., .2, 0, 0);
1093  fImage[specie]->cd();
1094  pS->Draw();
1095  pS->SetLogz();
1096  stat = pS;
1097  // status->DrawCopy("colz");
1098  }
1099  // AliWarningF("fImage[specie]=%p, plots=%p", fImage[specie], plots);
1100  // plots->cd();
1101  Int_t nx = 3;
1102  Int_t ny = (nImages + .5) / nx;
1103  plots->Divide(nx, ny, 0, 0);
1104  // else fImage[specie]->Divide(nx, ny);
1105 
1106 
1107  // Loop over histograms
1108  TH1* hist = 0;
1109  Int_t nHist = list[specie]->GetEntriesFast();
1110  for (Int_t i = 0; i < nHist; i++) {
1111  hist = static_cast<TH1*>(list[specie]->At(i));
1112  if (!hist || !hist->TestBit(AliQAv1::GetImageBit())) continue;
1113  if (hist == status) continue;
1114  TString name(hist->GetName());
1115  Bool_t isROE = name.Contains("readouterrors", TString::kIgnoreCase);
1116 
1117  // Go to sub-pad
1118  TVirtualPad* pad = 0;
1119  if (isROE) pad = plots->cd(4);
1120  else pad = plots->cd(GetHalfringPad(hist));
1121 
1122  pad->SetRightMargin(0.01);
1123  if (!fDoScale) {
1124  pad->SetLeftMargin(0.10);
1125  pad->SetBottomMargin(0.10);
1126  }
1127 
1128  // Check for log scale
1129  Int_t logOpts = 0;
1130  logOpts |= CheckForLog(hist->GetXaxis(), pad, 1);
1131  logOpts |= CheckForLog(hist->GetYaxis(), pad, 2);
1132  logOpts |= CheckForLog(hist->GetZaxis(), pad, 3);
1133 
1134  // Figure out special cases
1135  TString opt("");
1136  if (isROE) {
1137  pad->SetRightMargin(0.15);
1138  pad->SetBottomMargin(0.10);
1139  // pad->SetTopMargin(0.02);
1140  opt="COLZ";
1141  }
1142  else {
1143  pad->SetGridx();
1144  pad->SetGridy();
1145  if (fDoScale) {
1146  hist->SetMinimum(min);
1147  hist->SetMaximum(max);
1148  }
1149  else {
1150  hist->SetMinimum();
1151  hist->SetMaximum();
1152  }
1153  }
1154  // Draw (As a copy)
1155  hist->DrawCopy(opt);
1156 
1157  // Special cases
1158  if (!name.Contains("readouterrors", TString::kIgnoreCase)) {
1159  gStyle->SetOptTitle(0);
1160  TPad* insert = new TPad("insert", "Zoom",
1161  .5,.5, .99, .95, 0, 0, 0);
1162  insert->SetTopMargin(0.01);
1163  insert->SetRightMargin(0.01);
1164  insert->SetFillColor(0);
1165  insert->SetBorderSize(1);
1166  insert->SetBorderMode(0);
1167  insert->Draw();
1168  insert->cd();
1169  if (logOpts & 0x1) insert->SetLogx();
1170  if (logOpts & 0x2) insert->SetLogy();
1171  if (logOpts & 0x4) insert->SetLogz();
1172  hist->GetXaxis()->SetRange(1, hist->GetNbinsX()/8);
1173  TH1* copy = hist->DrawCopy(opt);
1174  copy->GetXaxis()->SetNdivisions(408, false);
1175  // Restore full range
1176  hist->GetXaxis()->SetRange(0, 0);
1177  gStyle->SetOptTitle(1);
1178  }
1179  pad->cd();
1180  // Possibly restore the log options
1181  RestoreLog(hist->GetXaxis(), logOpts & 0x1);
1182  RestoreLog(hist->GetYaxis(), logOpts & 0x2);
1183  RestoreLog(hist->GetZaxis(), logOpts & 0x4);
1184  }
1185  if (status && stat) {
1186  stat->cd();
1187  status->DrawCopy("BOX TEXT");
1188  }
1189  // Print to a post-script file
1190  fImage[specie]->Print(outName, "ps");
1191  if (AliDebugLevel() > 0)
1192  fImage[specie]->Print(Form("%s_%d.png",
1194  AliQAChecker::Instance()->GetRunNumber()));
1195  }
1196 }
1197 
1198 //__________________________________________________________________
1199 //
1200 // EOF
1201 //
static const char * GetEventSpecieName(EventSpecie_t es)
static UInt_t GetImageBit()
Definition: AliQAv1.h:55
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
UShort_t CheckRaw(AliRecoParam::EventSpecie_t specie, TH1 *hist) const
AliQAv1::QABIT_t Quality2Bit(UShort_t qual) const
TStyle * gStyle
static UInt_t GetExpertBit()
Definition: AliQAv1.h:54
#define AliDebugLevel()
#define TObjArray
UShort_t BasicCheck(TH1 *hist) const
Double_t landau(Double_t *xp, Double_t *pp)
Definition: DrawESD.C:38
UShort_t CheckSim(AliRecoParam::EventSpecie_t specie, TH1 *hist) const
const char * path
UShort_t CheckESD(AliRecoParam::EventSpecie_t specie, TH1 *hist) const
ALITASK_t
Definition: AliQAv1.h:26
QABIT_t
Definition: AliQAv1.h:28
Float_t p[]
Definition: kNNTest.C:133
static const char * GetImageFileFormat()
Definition: AliQAv1.h:57
MODE_t
Definition: AliQAv1.h:32
void MakeImage(TObjArray **list, AliQAv1::TASKINDEX_t task, AliQAv1::MODE_t mode)
Quality assurance checker for the FMD.
static const char * GetImageFileName()
Definition: AliQAv1.h:56
Double_t fELossFkupChi2Nu
void Check(Double_t *rv, AliQAv1::ALITASK_t what, TObjArray **list, const AliDetectorRecoParam *t)
TList * fExternParamList
flag to print the images or not
TStatToolkit stat
Definition: AnalyzeLaser.C:5
Double_t fELossBadChi2Nu
AliCDBEntry * Get(const AliCDBId &query, Bool_t forceCaching=kFALSE)
void sum()
Bool_t bg
Definition: RunAnaESD.C:6
UShort_t CheckRec(AliRecoParam::EventSpecie_t specie, TH1 *hist) const
TObject * GetObject()
Definition: AliCDBEntry.h:56
static DETECTORINDEX_t GetDetIndex(const char *name)
Definition: AliQAv1.cxx:348
Double_t chi2
Definition: AnalyzeLaser.C:7
static EventSpecie_t ConvertIndex(Int_t index)
Double_t fROErrorsFkup
Double_t fELossMinSharing
#define AliWarningF(message,...)
Definition: AliLog.h:556
Definition: AliCDBEntry.h:18
void ProcessExternalParam(const char *name, Double_t &v)
#define AliDebug(logLevel, message)
Definition: AliLog.h:300
static const char * GetModeName(MODE_t mode)
Definition: AliQAv1.h:91
void Set(QABIT_t bit, AliRecoParam::EventSpecie_t es)
Definition: AliQAv1.cxx:763
UShort_t CheckOne(AliQAv1::ALITASK_t what, AliRecoParam::EventSpecie_t specie, TH1 *hist) const
void SetQA(AliQAv1::ALITASK_t index, Double_t *values) const
TASKINDEX_t
Definition: AliQAv1.h:30
Bool_t IsEventSpecieSet(AliRecoParam::EventSpecie_t es) const
Definition: AliQAv1.h:92
void AddStatusPave(TH1 *hist, Int_t qual, Double_t xl=-1, Double_t yl=-1, Double_t xh=-1, Double_t yh=-1) const
static AliCDBManager * Instance(TMap *entryCache=NULL, Int_t run=-1)
void res(Char_t i)
Definition: Resolution.C:2
static Int_t GetQADebugLevel()
Definition: AliQAv1.h:72
#define AliDebugF(logLevel, format,...)
Definition: AliLog.h:338
static AliQAChecker * Instance()
Double_t fELossGoodParError
static Int_t GetHalfringIndex(UShort_t det, Char_t ring, UShort_t board, UShort_t monitor=0)
static TString GetTaskName(UInt_t tsk)
Definition: AliQAv1.h:90
UShort_t CheckFit(TH1 *hist, const TFitResultPtr &res, Double_t low, Double_t high, Int_t &color) const
static AliQAv1 * Instance()
Definition: AliQAv1.cxx:585
TObject * GetThreshold(Int_t i)