]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0QAChecker.cxx
Upadte from Taku
[u/mrichter/AliRoot.git] / T0 / AliT0QAChecker.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 //  Checks the quality assurance. 
17 //  By comparing with reference data
18 //  Skeleton for T0
19 //---------------------------------------------
20 //checkig without reference data:
21 //for RAW QA all histograms should have approximatly the same 
22 //number of entries as RefPoint
23 //for Rec Points checks 
24 //  - amplitude measured by 2 methos
25 //  - online and offline T0 measurements
26 // for ESD quality of reconstruction ( and measurements):
27 // RMS of vertex and T0 less than 75ps
28 //
29 // Alla.Maevskaya@cern.ch   
30 //...
31
32 // --- ROOT system ---
33 #include <Riostream.h>
34 #include <TClass.h>
35 #include <TH1F.h> 
36 #include <TF1.h> 
37 #include <TFitResultPtr.h>
38 #include <TH2.h> 
39 #include <TIterator.h> 
40 #include <TKey.h> 
41 #include <TFile.h> 
42 #include <TMath.h>
43 #include <TString.h>
44 #include <TPaveText.h>
45
46 // --- Standard library ---
47
48 // --- AliRoot header files ---
49 #include "AliLog.h"
50 #include "AliQAv1.h"
51 #include "AliQAChecker.h"
52 #include "AliCDBEntry.h"
53 #include "AliQAManager.h"
54 #include "AliT0QAChecker.h"
55
56 ClassImp(AliT0QAChecker)
57 //____________________________________________________________________________
58 AliT0QAChecker::AliT0QAChecker() :
59 AliQACheckerBase("T0","T0 Quality Assurance Checker")
60
61 {
62   // Standard constructor
63
64 }
65
66 //____________________________________________________________________________
67 AliT0QAChecker::AliT0QAChecker(const AliT0QAChecker& qac):
68   AliQACheckerBase(qac.GetName(), qac.GetTitle()) 
69 {
70   // copy constructor
71   AliError("Copy should not be used with this class\n");
72 }
73 //____________________________________________________________________________
74 AliT0QAChecker& AliT0QAChecker::operator=(const AliT0QAChecker& qac){
75   // assignment operator
76   this->~AliT0QAChecker();
77   new(this)AliT0QAChecker(qac);
78   return *this;
79 }
80
81
82 //____________________________________________________________________________
83 AliT0QAChecker::~AliT0QAChecker(){
84   // destructor
85
86 }
87
88 //__________________________________________________________________
89 void AliT0QAChecker::Check(Double_t *  test, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
90 {
91   // Super-basic check on the QA histograms on the input list:
92   // look whether they are empty!
93     
94  char * detOCDBDir = Form("T0/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
95
96   AliCDBEntry *QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
97   //  QARefRec->Dump();
98   if( !QARefRec){
99     AliInfo("QA reference data NOT retrieved for Reconstruction check. No T0 reference distribution");
100   }
101
102     
103   for(Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++){ 
104     test[specie]    = 1.0; //FK//  initiate qa flag for the whole set of histograms as good 
105   }
106
107
108   for(Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
109     if(!(AliQAv1::Instance()->IsEventSpecieSet(specie) && list[specie]) || list[specie]->GetEntries() == 0) {
110       continue;
111     }
112
113     if(index == AliQAv1::kRAW){
114
115       if(AliRecoParam::ConvertIndex(specie) == AliRecoParam::kCalib){//      if (index == AliQAv1::kRAW )
116         //check laser data efficiencies   
117         Double_t qaFlag = CheckLaser(list[specie]);
118         if(qaFlag < test[specie]) test[specie] = qaFlag;
119       }
120
121       if(AliRecoParam::ConvertIndex(specie) == AliRecoParam::kCalib   ||
122          AliRecoParam::ConvertIndex(specie) == AliRecoParam::kHighMult ||
123          AliRecoParam::ConvertIndex(specie) == AliRecoParam::kLowMult){ 
124          //AliRecoParam::ConvertIndex(specie) == AliRecoParam::kDefault ||
125          
126         //check BCID   
127         Double_t qaFlag = CheckBCID(list[specie]);
128         if(qaFlag < test[specie]) test[specie] = qaFlag;
129       }
130
131       if(AliRecoParam::ConvertIndex(specie) == AliRecoParam::kHighMult ||
132         AliRecoParam::ConvertIndex(specie) == AliRecoParam::kLowMult){ 
133         //AliRecoParam::ConvertIndex(specie) == AliRecoParam::kDefault ||
134         //check physics 
135         Double_t qaFlag = CheckRaw(list[specie]);
136         if(qaFlag < test[specie]) test[specie] = qaFlag;
137       }
138     }
139
140     if(index == AliQAv1::kESD && AliRecoParam::Convert(specie) != AliRecoParam::kCalib){
141       test[specie] = CheckESD(list[specie]);
142     } 
143   }
144 }
145
146 //--------------------------------------------------------------------------
147 Double_t AliT0QAChecker::CheckLaser(TObjArray *listrec) const {
148    
149   TH1 *hdata; 
150   TH1 *fhRawEff[10];
151   Int_t nEffHistos=0;
152
153   //thresholds for warning and error on efficiencies
154   Float_t thrWarning = 0.5; //FK//  warning level
155   Float_t thrError   = 0.2; //FK//  error level
156
157   const int kNumberOfHistos = 3; 
158   Int_t consecutiveHistoNumber[kNumberOfHistos] = { 207, 208, 209};  //Checked histos   fhCDFeff, hEffLED, hEffQTC
159   Int_t qualityFlag[kNumberOfHistos]; //quality flag for a given histogram
160
161   for(Int_t ir=0; ir<kNumberOfHistos; ir++){
162     qualityFlag[ir] = kT0Info; //init quality flag for a given histogram
163
164     hdata = (TH1*) listrec->UncheckedAt(consecutiveHistoNumber[ir]);
165     if(hdata){
166       fhRawEff[nEffHistos] = hdata;
167       nEffHistos++;
168     }
169   }
170      
171   TLine linelowyellow(0, thrWarning, 24, thrWarning);    
172   linelowyellow.SetLineColor(5);
173   linelowyellow.SetLineStyle(3);
174   linelowyellow.SetLineWidth(4);
175   TLine linelowred(0, thrError, 24, thrError);    
176   linelowred.SetLineColor(2);
177   linelowred.SetLineStyle(3);
178   linelowred.SetLineWidth(4);
179
180   Bool_t bEffHistosNotEmpty = kFALSE; //check if all histograms have some counts
181  
182   for(Int_t ih = 0; ih < nEffHistos; ih++){
183      
184     EraseOldMessages((TH1*) fhRawEff[ih]);// clean objects added at previous checks
185  
186     fhRawEff[ih]->SetLineWidth(2);
187     fhRawEff[ih]->SetMaximum(2.);
188     fhRawEff[ih]->SetMinimum(0.);
189     fhRawEff[ih]->GetListOfFunctions()->Add((TLine*)linelowyellow.Clone());
190     fhRawEff[ih]->GetListOfFunctions()->Add((TLine*)linelowred.Clone());
191
192     if(fhRawEff[ih]->Integral()>0) bEffHistosNotEmpty = kTRUE; //this histo does have some counts in it
193
194     Int_t nbins= fhRawEff[ih]->GetNbinsX();
195     for(Int_t ib=1; ib<=nbins; ib++){ //loop over bins and check if the efficiency is above level
196         
197       Float_t chcont = fhRawEff[ih]->GetBinContent(ib);
198       if(chcont < thrWarning && qualityFlag[ih] > kT0Error  ) qualityFlag[ih] = kT0Warning;//Warning level
199       if(chcont < thrError)  qualityFlag[ih] = kT0Error;//Error level
200     }
201    
202     if(qualityFlag[ih] == kT0Info ){
203       AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 efficiency  %s  is good", fhRawEff[ih]->GetName() ));
204     }else if(qualityFlag[ih] == kT0Warning){ 
205       AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 efficiency  %s  is not so good", fhRawEff[ih]->GetName() ));
206     }else if(qualityFlag[ih] == kT0Error){
207       AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 efficiency  %s  is not good", fhRawEff[ih]->GetName() ));
208     }
209   }
210
211   //executive summary
212   int lowestQualityFlag = (int) kT0Info;
213   for(Int_t ih = 0; ih < nEffHistos; ih++){
214
215     if(!bEffHistosNotEmpty){ //all laser efficiency plots are empty
216
217       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
218       text.AddText(Form("1) T0 is in BEAMTUNIG: empty plots are ok")); 
219       text.AddText(Form("2) T0 is in READY: check calibriation trigger")); 
220       text.AddText(Form("if also physics data are empty report"));
221       text.AddText(Form("readout problem to the T0 on-call expert")); 
222       text.SetBorderSize(0);
223       text.SetFillStyle(0);
224       fhRawEff[ih]->GetListOfFunctions()->Add((TPaveText*)text.Clone());              
225     }
226  
227     if( qualityFlag[ih] <lowestQualityFlag )  lowestQualityFlag = qualityFlag[ih];
228   }
229    
230   return ConvertQualityFlagToDouble(lowestQualityFlag); 
231 }
232 //--------------------------------------------------------------------------
233 Double_t AliT0QAChecker::CheckBCID(TObjArray *listrec) const {
234    
235   Int_t qualityFlagBCID = kT0Info; //init quality flag for a given histogram; 
236
237   TH2F *hBCID = (TH2F*) listrec->UncheckedAt(224); //BCID versus TRM  BCID
238
239      
240   // clean objects added at previous checks
241   EraseOldMessages((TH1*)hBCID);
242
243   if(hBCID->Integral()>0){
244     //BCID does have some counts in it
245
246     Int_t nbinsX = hBCID->GetNbinsX();
247     Int_t nbinsY = hBCID->GetNbinsY();
248     double entriesOnDiagonal  = 0; //count diagonal and off diagonal entries
249     double entriesOffDiagonal = 0;
250     Int_t offset = 37; 
251     Int_t period = 3564; 
252
253     for(Int_t itrm=1; itrm<=nbinsX; itrm++){ //BCID TRM
254       for(Int_t ibcid=1; ibcid<=nbinsY; ibcid++){ //BCID
255         if(TMath::Abs(itrm-ibcid)<6 || TMath::Abs( (itrm+offset)%period - ibcid)<2 ){ 
256            entriesOnDiagonal  += hBCID->GetBinContent(itrm,ibcid); //On  Diagonal
257         }else{
258          entriesOffDiagonal += hBCID->GetBinContent(itrm,ibcid); //Off Diagonal
259         }
260       }
261     }
262     if(entriesOnDiagonal<1 || entriesOffDiagonal>20){
263       qualityFlagBCID = kT0Warning; //no entries on diagonal
264       AliDebug(AliQAv1::GetQADebugLevel(), Form("T0   %s is not diagonal", hBCID->GetName() ));
265
266       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
267       text.AddText(Form("Check if entries are on a diagonal.")); 
268       text.AddText(Form("Report readout problem to the T0 on-call expert")); 
269       text.SetBorderSize(0);
270       text.SetFillStyle(0);
271       hBCID->GetListOfFunctions()->Add((TPaveText*)text.Clone());             
272     }
273   }else{ //BCID empty
274
275     qualityFlagBCID = kT0Error;
276     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 :  %s has NO entries", hBCID->GetName() ));
277
278     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
279     text.AddText(Form("NO ENTRIES!!!")); 
280     text.AddText(Form("If T0 is READY report")); 
281     text.AddText(Form("readout problem to the T0 on-call expert"));
282     text.SetBorderSize(0);
283     text.SetFillStyle(0);
284     hBCID->GetListOfFunctions()->Add((TPaveText*)text.Clone());       
285   }
286
287   //executive summary
288   int lowestQualityFlag = (int) qualityFlagBCID;
289
290   return ConvertQualityFlagToDouble(lowestQualityFlag); 
291   
292 }
293
294 //--------------------------------------------------------------------------
295 Double_t AliT0QAChecker::CheckRaw(TObjArray *listrec) const {
296    
297
298   Int_t qualityFlagTrigger = kT0Info; //init quality flag for a given histogram; 
299
300   TH1F *hTrigger = (TH1F*) listrec->UncheckedAt(169);//hRawTrigger 
301
302      
303   // clean objects added at previous checks
304   EraseOldMessages((TH1*) hTrigger); 
305
306   if(hTrigger->Integral()>0){
307     //trigger plot does have some counts in it
308     //are Mean, ORA and ORC not empty?  
309     if( hTrigger->GetBinContent(1)<0.001 || hTrigger->GetBinContent(3)<0.001 || hTrigger->GetBinContent(4)<0.001){
310       qualityFlagTrigger = kT0Error; //no entries on diagonal
311       AliDebug(AliQAv1::GetQADebugLevel(), Form("T0: too little ORA and ORC in  %s", hTrigger->GetName() ));
312
313       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
314       text.AddText(Form("Check ORA and ORC")); 
315       text.AddText(Form("Report problem to the T0 on-call expert")); 
316       text.SetBorderSize(0);
317       text.SetFillStyle(0);
318       hTrigger->GetListOfFunctions()->Add((TPaveText*)text.Clone());          
319     }
320
321   }else{ //Trigger histo empty
322
323     qualityFlagTrigger = kT0Error;
324     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hTrigger->GetName() ));
325
326     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
327     text.AddText(Form("NO ENTRIES!!!")); 
328     text.AddText(Form("If T0 is READY report")); 
329     text.AddText(Form("readout problem to the T0 on-call expert")); 
330     text.SetBorderSize(0);
331     text.SetFillStyle(0);
332     hTrigger->GetListOfFunctions()->Add((TPaveText*)text.Clone());            
333
334   }
335
336   //--------- timing plots ---------
337   Int_t qualityFlagMeanBest = kT0Info; //init quality flag for a given histogram;
338   TH1F *hMeanBest    = (TH1F*) listrec->UncheckedAt(223);  //fhMeanBest  time
339   EraseOldMessages((TH1*) hMeanBest); 
340
341   if(hMeanBest->Integral()<1){
342     qualityFlagMeanBest = kT0Error;
343     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hMeanBest->GetName() ));
344
345     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
346     text.AddText(Form("NO ENTRIES!!!")); 
347     text.AddText(Form("If T0 is READY and beam is on report")); 
348     text.AddText(Form("the problem to the T0 on-call expert"));
349     text.SetBorderSize(0);
350     text.SetFillStyle(0);
351     hMeanBest->GetListOfFunctions()->Add((TPaveText*)text.Clone());           
352   }
353
354
355
356
357   //----------- vetext plots ---------
358
359   Int_t qualityFlagVetrexFirst = kT0Info; //init quality flag for a given histogram;
360   Int_t qualityFlagVetrexBest  = kT0Info; //init quality flag for a given histogram;
361  
362   TH1F *hVertexFirst     = (TH1F*) listrec->UncheckedAt(225);//fhVertex1st
363   TH1F *hVertexBest      = (TH1F*) listrec->UncheckedAt(226);//fhVertexBest
364   TH1F *hOrCminOrATvdcOn = (TH1F*) listrec->UncheckedAt(217);//vertex with TVDC on
365
366   // clean objects added at previous checks
367   EraseOldMessages((TH1*) hVertexFirst); 
368   EraseOldMessages((TH1*) hVertexBest); 
369   float errorLevelDifference   = 400; 
370   float warningLevelDifference = 200; 
371
372   if(hVertexFirst->Integral()>0 && hOrCminOrATvdcOn->Integral()>0){
373     Float_t thr = 0.05 * hVertexFirst->GetBinContent(hVertexFirst->GetMaximumBin());
374      
375     float meanVertexFirst = GetMeanAboveThreshold( hVertexFirst, thr);
376     float meanVertexTVDC  = GetMeanAboveThreshold( hOrCminOrATvdcOn,thr);
377     float diff = TMath::Abs(meanVertexFirst - meanVertexTVDC);
378     if(diff > errorLevelDifference){
379       qualityFlagVetrexFirst = kT0Warning; //init quality flag for a given histogram;
380
381       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
382       text.AddText(Form("large diff. between TVDC vertex")); 
383       text.AddText(Form("and the first vertex."));
384       text.SetBorderSize(0);
385       text.SetFillStyle(0);
386       hVertexFirst->GetListOfFunctions()->Add((TPaveText*)text.Clone());              
387
388     }else if(diff > warningLevelDifference){
389       qualityFlagVetrexFirst = kT0Warning; //init quality flag for a given histogram;
390
391       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
392       text.AddText(Form("large diff. between TVDC vertex")); 
393       text.AddText(Form("and the first vertex."));
394       text.SetBorderSize(0);
395       text.SetFillStyle(0);
396       hVertexFirst->GetListOfFunctions()->Add((TPaveText*)text.Clone());              
397
398     }
399   }else{
400     qualityFlagVetrexFirst = kT0Error;
401     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hVertexFirst->GetName() ));
402
403     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
404     text.AddText(Form("NO ENTRIES!!!")); 
405     text.AddText(Form("If T0 is READY and beam is on report")); 
406     text.AddText(Form("the problem to the T0 on-call expert"));
407     text.SetBorderSize(0);
408     text.SetFillStyle(0);
409     hVertexFirst->GetListOfFunctions()->Add((TPaveText*)text.Clone());        
410   }
411
412   if(hVertexBest->Integral()>0 && hOrCminOrATvdcOn->Integral()>0){
413
414     Float_t thr = 0.05 * hVertexBest->GetBinContent(hVertexBest->GetMaximumBin());
415     float meanVertexBest  = GetMeanAboveThreshold(hVertexBest,thr);
416     float meanVertexTVDC =  GetMeanAboveThreshold(hOrCminOrATvdcOn,thr);
417     float diff = TMath::Abs(meanVertexBest - meanVertexTVDC);
418     if(diff > errorLevelDifference){
419       qualityFlagVetrexBest = kT0Warning; //init quality flag for a given histogram;
420
421       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
422       text.AddText(Form("Large. diff. between TVDC vertex")); 
423       text.AddText(Form("and the best vertex."));
424       text.SetBorderSize(0);
425       text.SetFillStyle(0);
426       hVertexBest->GetListOfFunctions()->Add((TPaveText*)text.Clone());       
427
428     }else if(diff > warningLevelDifference){
429       qualityFlagVetrexBest = kT0Warning; //init quality flag for a given histogram;
430
431       TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
432       text.AddText(Form("Large. diff. between TVDC vertex")); 
433       text.AddText(Form("and the best vertex."));
434       text.SetBorderSize(0);
435       text.SetFillStyle(0);
436       hVertexBest->GetListOfFunctions()->Add((TPaveText*)text.Clone());       
437     }
438   }else{
439     qualityFlagVetrexBest = kT0Error;
440     AliDebug(AliQAv1::GetQADebugLevel(), Form("T0 histogram  %s has NO entries", hVertexBest->GetName() ));
441
442     TPaveText text(0.20,0.50,0.99,0.99,"NDC");   
443     text.AddText(Form("NO ENTRIES!!!")); 
444     text.AddText(Form("If T0 is READY and beam is on report")); 
445     text.AddText(Form("the problem to the T0 on-call expert"));
446     text.SetBorderSize(0);
447     text.SetFillStyle(0);
448     hVertexBest->GetListOfFunctions()->Add((TPaveText*)text.Clone());         
449   }
450
451
452   //executive summary
453   int lowestQualityFlag = (int) qualityFlagTrigger;
454   if(qualityFlagVetrexBest<lowestQualityFlag)  lowestQualityFlag = qualityFlagVetrexBest;
455   if(qualityFlagVetrexFirst<lowestQualityFlag) lowestQualityFlag = qualityFlagVetrexFirst;
456   if(qualityFlagMeanBest <lowestQualityFlag)   lowestQualityFlag = qualityFlagMeanBest; 
457   
458
459   return ConvertQualityFlagToDouble(lowestQualityFlag); 
460   
461 }
462
463 //--------------------------------------------------------------------------
464 Double_t AliT0QAChecker::CheckESD(TObjArray *listrec ) const
465 {
466   Float_t checkr = 0;
467   TH1 *fhESD;
468  
469   fhESD  = (TH1*) listrec->UncheckedAt(2);
470   if(fhESD){
471     AliDebug(AliQAv1::GetQADebugLevel(), Form("count %s ", fhESD->GetName()) );
472      TF1 *f1 = new TF1("f1","gaus",-1,1);
473     fhESD->Fit("f1","R","Q", -1,1);
474     Double_t par[3];
475     f1->GetParameters(&par[0]);
476     
477     TPaveText text(0.30,0.50,0.99,0.99,"NDC");
478     
479     text.AddText(Form("T0 RUN %d ",AliCDBManager::Instance()->GetRun()));
480     
481     AliDebug(AliQAv1::GetQADebugLevel(), Form("numentries %d mean %f  #sigma %f", (int)fhESD->GetEntries(),par[1], par[2]));
482     
483     
484     if (par[2] > 0.07 && par[2] < 1.) {
485       checkr=0.5;
486        text.AddText(Form("not good resolution :\n %f ns\n", par[2] ));
487        text.SetFillColor(5);
488        printf("T0 detector resolution is not good enouph: %f ns\n",par[2] );
489     }
490     if(TMath::Abs(par[1])>0.05) {
491       checkr = 0.5;
492       text.AddText(Form(" Check clock shift on %f ns", par[1]));
493       text.SetFillColor(5);
494     }
495     if (par[2] >  1. || TMath::Abs(par[1])>0.1) {
496       checkr = 0.25;
497       text.AddText(Form(" Bad resolution:\n mean %f ns sigma %f ns", par[1], par[2]));
498       text.SetFillColor(2);
499       { // RS Clean previous additions
500         TList* lstF = fhESD->GetListOfFunctions();
501         if (lstF) {
502           TObject *stats = lstF->FindObject("stats");
503           lstF->Remove(stats);
504           TObject *obj;
505           while ((obj = lstF->First())) {
506             while(lstF->Remove(obj)) { }
507             delete obj;
508           }
509           if (stats) lstF->Add(stats);
510         } 
511       }
512       fhESD->GetListOfFunctions()->Add((TPaveText*)text.Clone());       
513       AliDebug(AliQAv1::GetQADebugLevel(),
514                Form("Please, check calibration: shift= %f resolution %f test=%f\n",
515                     par[1], par[2], checkr) ) ; 
516     }
517   }
518   else
519     {
520       AliDebug(AliQAv1::GetQADebugLevel(),
521                Form("No ESD QA histogram found, nothing to check"));
522         checkr=0;
523     }
524   
525   
526   return checkr;
527 }
528
529
530 //--------------------------------------------------------------------------
531 void AliT0QAChecker::EraseOldMessages(TH1* h) const 
532 {
533   //erase the old captions 
534   TList* lstF = h->GetListOfFunctions();
535   if(lstF){  
536      TObject *stats = lstF->FindObject("stats");
537      lstF->Remove(stats);
538      TObject *obj;
539      while ((obj = lstF->First())) {
540        while(lstF->Remove(obj)) { }
541          delete obj;
542     }
543     if (stats) lstF->Add(stats);
544   }
545 }
546 //--------------------------------------------------------------------------
547 Double_t AliT0QAChecker::ConvertQualityFlagToDouble(int qualityFlag) const 
548 {
549   //covert quality flag to double
550   Double_t checkr=1.0;
551
552   switch ( qualityFlag ){
553     case kT0Info:
554         checkr = 1.0; break;
555     case kT0Warning:
556         checkr = 0.75; break;
557     case kT0Error:
558           checkr = 0.25; break;
559     case kT0Fatal:
560         checkr = -1.0; break;
561     default:
562          AliError("Invalid ecc value. FIXME !");
563          checkr = 0.25; break;
564   };
565
566   return checkr; 
567 }
568  
569 //--------------------------------------------------------------------------
570 Float_t AliT0QAChecker::GetMeanAboveThreshold(TH1F* hV, Float_t thr) const{
571   //caculate mean value of histo bins above threshold
572   Int_t nBins = hV->GetNbinsX();
573   Int_t nBinsAboveThr = 0; 
574   Float_t sum = 0;
575
576   for(Int_t ib=1;ib<=nBins;ib++){
577     Float_t val = hV->GetBinContent(ib);
578     if(val<thr){
579       sum+=val;
580       nBinsAboveThr++; 
581     }
582   }
583
584   if(nBinsAboveThr>0) return sum/nBinsAboveThr;
585   else return hV->GetMean();
586 }