]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQAChecker.cxx
Setting of run number for the CDB manager deleted (fixing bug)
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQAChecker.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
17 /* $Id$ */
18
19 //...
20 //  Checks the quality assurance. 
21 //  By comparing with reference data
22 //  Skeleton for HMPID
23 //...
24
25 // --- ROOT system ---
26 #include <TClass.h>
27 #include <TH1F.h> 
28 #include <TH1I.h>
29 #include <TH2.h> 
30 #include <TF1.h> 
31 #include <TIterator.h> 
32 #include <TKey.h> 
33 #include <TFile.h>
34 #include <TLine.h>
35 #include <TParameter.h> 
36 #include <TPaveText.h>
37 // --- Standard library ---
38
39 // --- AliRoot header files ---
40 #include "AliLog.h"
41 #include "AliQAv1.h"
42 #include "AliQAChecker.h"
43 #include "AliHMPIDQAChecker.h"
44 #include "AliCDBEntry.h"
45 #include "AliQAManager.h"
46 #include "AliQAThresholds.h"
47
48 ClassImp(AliHMPIDQAChecker)
49  //_________________________________________________________________
50 AliHMPIDQAChecker::AliHMPIDQAChecker() : 
51 AliQACheckerBase("HMPID","HMPID Quality Assurance Data Checker"), 
52 fNoReference(kTRUE),
53 fQARefRec(NULL),
54
55 fHmpQaThr_NumberOfExcludedDDL(0),
56 fHmpQaThr_DataSizeLowerThreshold(900),
57 fHmpQaThr_DataSizeUpperThreshold(1500),
58 fHmpQaThr_PadOccupancyLowerThreshold(0.005),
59 fHmpQaThr_PadOccupancyUpperThreshold(0.8),
60 fHmpQaThr_SectorGainLossWarningThreshold(3),
61 fHmpQaThr_SectorGainLossErrorThreshold(6),
62 fHmpQaThr_MissingPadFractionWarningThreshold(0.3),
63 fHmpQaThr_MissingPadFractionErrorThreshold(0.5),
64 fIsOnlineThr(kFALSE)                                 
65  
66
67
68 {
69     //ctor, fetches the reference data from OCDB 
70   char * detOCDBDir = Form("HMPID/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ; 
71   AliCDBEntry * QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
72   if(QARefRec) {
73     fQARefRec = dynamic_cast<TObjArray*> (QARefRec->GetObject()) ; 
74     if (fQARefRec)
75       if (fQARefRec->GetEntries()) 
76         fNoReference = kFALSE ;            
77     if (fNoReference) 
78       AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker!");
79   }
80   
81   InitOnlineThresholds();
82   
83 }
84
85 //_________________________________________________________________
86 AliHMPIDQAChecker::~AliHMPIDQAChecker() 
87 {
88   if(fQARefRec) { fQARefRec->Delete() ;   delete fQARefRec ; }
89 }
90 //_________________________________________________________________
91 void AliHMPIDQAChecker::Check(Double_t *  check, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/) 
92 {
93 //
94 // Main check function: Depending on the TASK, different checks are applied
95 // At the moment:       check for empty histograms and checks for RecPoints
96
97   if(fNoReference)  
98
99   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
100     check[specie] = 1.0;    
101     //printf("+++++++++++++++++++++ specie %d name: %s \n",specie,AliRecoParam::GetEventSpecieName(specie));
102     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) continue ;     
103     // checking for empy histograms
104     if(CheckEntries(list[specie]) == 0)  {
105       AliWarning("histograms are empty");
106       check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
107     }
108   
109     check[specie] = AliQAv1::kINFO ;
110     
111     
112     //check sim
113     if(index == AliQAv1::kSIM) check[specie] = CheckSim(list[specie], fQARefRec);
114
115     // checking rec points
116     if(index == AliQAv1::kREC) check[specie] = CheckRec(list[specie], fQARefRec);
117    
118     //checking raw data
119     if(index == AliQAv1::kRAW) {       check[specie] = CheckRaw(specie,list[specie]);       }
120                                
121     
122
123
124   } // species loop
125 }
126 //_________________________________________________________________
127 Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
128 {
129   //
130   //  check on the QA histograms on the input list: 
131   // 
132
133   Double_t test = 0.0  ;
134   Int_t count = 0 ; 
135   
136   if (list->GetEntries() == 0){  
137     test = 1. ; // nothing to check
138   }
139   else {
140     TIter next(list) ; 
141     TH1 * hdata ;
142     count = 0 ; 
143     while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
144       if (hdata) { 
145         Double_t rv = 0.;
146         if(hdata->GetEntries()>0)rv=1; 
147         count++ ; 
148         test += rv ; 
149       }
150       else{
151         AliError("Data type cannot be processed") ;
152       }
153       
154     }
155     if (count != 0) { 
156       if (test==0) {
157         AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
158         test = 0.;  //upper limit value to set kWARNING flag for a task
159       }
160       else test = 1 ;
161     }
162   }
163
164   return test ; 
165 }  
166 //_________________________________________________________________
167 Double_t AliHMPIDQAChecker::CheckSim(TObjArray *listsim, TObjArray *listref) const
168 {
169   //
170   //  check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
171   //
172
173    Float_t checkresponse = 0;
174
175    Float_t counter = 0 ;
176    TIter next(listsim) ;
177    TH1* histo;
178    while ( (histo = dynamic_cast<TH1 *>(next())) ) {
179      //PH The histogram should have at least 10 bins with at least 5 entries
180      Int_t nbinsabove = 0;
181      for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) { 
182        if (histo->GetBinContent(ibin)>5) nbinsabove++;
183      }
184
185    if( nbinsabove < 10 ) counter++;
186    else {
187     TString h = histo->GetTitle();
188     if(h.Contains("Zoom")){
189     histo->Fit("expo","LQ0","",5,50);
190     if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
191     }
192     if(h.Contains("size  MIP"))   if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
193     if(h.Contains("size  Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
194     if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
195     AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));  
196    }
197   }
198  Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listsim->GetEntries())
199  
200  if(response < 0.1) checkresponse = 0.7;      // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
201  else if(response < 0.5) checkresponse = 0.4; //  50%  of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
202  else checkresponse = 0.001;                  // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
203  return checkresponse;
204 }
205
206 //___________________________________________________________________________________________________
207 Double_t AliHMPIDQAChecker::CheckRec(TObjArray *listrec, TObjArray *listref) const
208 {
209   //
210   //  check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
211   //
212
213    Float_t checkresponse = 0;
214
215    Float_t counter = 0 ;
216    TIter next(listrec) ;
217    TH1* histo;
218    while ( (histo = dynamic_cast<TH1 *>(next())) ) {
219      //PH The histogram should have at least 10 bins with at least 5 entries
220      Int_t nbinsabove = 0;
221      for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) { 
222        if (histo->GetBinContent(ibin)>5) nbinsabove++;
223      }
224
225    if( nbinsabove < 10 ) counter++;
226    else {
227     TString h = histo->GetTitle();
228     if(h.Contains("Zoom")){
229     histo->Fit("expo","LQ0","",5,50);
230     if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
231     }
232     if(h.Contains("size  MIP"))   if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
233     if(h.Contains("size  Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
234     if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
235     AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));  
236    }
237   }
238  Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
239  
240  if(response < 0.1) checkresponse = 0.7;      // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
241  else if(response < 0.5) checkresponse = 0.4; //  50%  of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
242  else checkresponse = 0.001;                  // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
243  return checkresponse;
244 }
245 //___________________________________________________________________________________________________
246 Double_t AliHMPIDQAChecker::CheckRaw(Int_t specie, TObjArray* list)
247 {
248   //
249   // Check the raw data for offline / online using default or updated thresholds from AMORE
250   // As of now (06/07/2012) the quality flag of all histos in raw will be est to the result of teh CheckRaw and displayed
251   // in AMORE. But we can pu undividual labels.
252   //
253   Int_t raqQualFlag = AliQAv1::kNULLBit;
254   
255   Int_t hmpQaFlags[4]={-1}; //init for the 4 shifter histos
256   
257   TString histname ="null";
258   TPaveText text(0.65,0.8,0.9,0.99,"NDC"); 
259   TPaveText text1(0.65,0.8,0.9,0.99,"NDC"); 
260    
261  
262   Int_t entries = list->GetEntriesFast();
263   
264   if ( entries == 0 ) {
265      AliWarning(Form("HMPID QA Checker RAWS: no object to analyse! Exiting..."));
266      return AliQAv1::kFATAL;
267   }
268   
269   TLine* lineDdlDataSizeMin = new TLine(1536,fHmpQaThr_DataSizeUpperThreshold,1548,fHmpQaThr_DataSizeUpperThreshold);
270   TLine* lineDdlDataSizeMax = new TLine(1536,fHmpQaThr_DataSizeLowerThreshold,1548,fHmpQaThr_DataSizeLowerThreshold);
271   
272   TLine* linePadOccupMin = new TLine(1536,fHmpQaThr_PadOccupancyUpperThreshold,1548,fHmpQaThr_PadOccupancyUpperThreshold);
273   TLine* linePadOccupMax = new TLine(1536,fHmpQaThr_PadOccupancyLowerThreshold,1548,fHmpQaThr_PadOccupancyLowerThreshold);
274   
275   
276   Int_t badDdlCnt = 0, badOccCnt = 0;
277   
278   //___ check data size per ddl
279   if( list->FindObject(Form("%s_hHmpDdlDataSize",AliRecoParam::GetEventSpecieName(specie)) )) 
280     {
281       TH1* h1 = dynamic_cast<TH1*>( list->FindObject(Form("%s_hHmpDdlDataSize",AliRecoParam::GetEventSpecieName(specie)))); 
282       if(h1) {
283       if( h1->Integral() > 1 ) {
284       // no entres -> fatal
285       if( h1->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
286       h1->SetStats(0);
287       
288       
289       // clean up the text, stat, lines, ...
290       if( h1->GetListOfFunctions() ) h1->GetListOfFunctions()->Clear();
291       
292       for ( Int_t iddl = 1 ; iddl <= 14; iddl++)
293         {
294           if( h1->GetBinContent(iddl) < fHmpQaThr_DataSizeLowerThreshold || h1->GetBinContent(iddl) > fHmpQaThr_DataSizeUpperThreshold ) badDdlCnt++; 
295         }
296       //___ check if one or more DDLs are excluded
297         
298       
299       badDdlCnt -= fHmpQaThr_NumberOfExcludedDDL;
300       
301       
302       
303       //___ check how many are bad
304       if ( badDdlCnt == 0 )  
305         {         
306           hmpQaFlags[0] = AliQAv1::kINFO;
307           text1.Clear();
308           text1.AddText(Form("OK (%d)",fIsOnlineThr));
309           text1.SetFillColor(kGreen);
310           h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());      
311           lineDdlDataSizeMin->SetLineColor(kGreen);
312           lineDdlDataSizeMax->SetLineColor(kGreen);
313           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
314           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());
315         }
316       else if ( badDdlCnt == 1 )  
317         {         
318           hmpQaFlags[0]  = AliQAv1::kWARNING;
319           text1.Clear();
320           text1.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
321           text1.SetFillColor(kOrange);
322           h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());     
323           lineDdlDataSizeMin->SetLineColor(kOrange);
324           lineDdlDataSizeMax->SetLineColor(kOrange);
325           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
326           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());   
327         }
328       else if (  badDdlCnt >= 2 )  
329         {
330           hmpQaFlags[0]  = AliQAv1::kERROR;         
331           text1.Clear();
332           text1.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
333           text1.SetFillColor(kRed);
334           h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());           
335           lineDdlDataSizeMin->SetLineColor(kRed);
336           lineDdlDataSizeMax->SetLineColor(kRed);
337           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
338           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());   
339         }
340       else 
341         {
342           hmpQaFlags[0]  = AliQAv1::kFATAL;
343           text1.Clear();
344           text1.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
345           text1.SetFillColor(kRed);
346           h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());           
347           lineDdlDataSizeMin->SetLineColor(kRed);
348           lineDdlDataSizeMax->SetLineColor(kRed);
349           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
350           h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());   
351         }
352       }
353       }//the histo is filled
354     }//___hHmpDdlDataSize
355   
356   
357   
358   if( list->FindObject(Form("%s_fHmpPadOcc",AliRecoParam::GetEventSpecieName(specie)) )) 
359     {
360      
361       
362       TH1* h1 = dynamic_cast<TH1*>( list->FindObject(Form("%s_fHmpPadOcc",AliRecoParam::GetEventSpecieName(specie)))); 
363        if(h1) {
364        if( h1->Integral() > 1 ) {
365       // no entres -> fatal
366       if( h1->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
367        h1->SetStats(0);
368       // clean up the text, stat, lines, ...
369       if( h1->GetListOfFunctions() ) h1->GetListOfFunctions()->Clear();
370       
371       for ( Int_t iddl = 1 ; iddl <= 14; iddl++)
372         {
373           if( h1->GetBinContent(iddl) < fHmpQaThr_PadOccupancyLowerThreshold || h1->GetBinContent(iddl) > fHmpQaThr_PadOccupancyUpperThreshold ) badOccCnt++; 
374         }
375        
376       badOccCnt -= fHmpQaThr_NumberOfExcludedDDL;
377          
378       //___ check how many are bad
379       if ( badOccCnt == 0 )  
380         {
381           hmpQaFlags[1]  = AliQAv1::kINFO;
382           text.Clear();
383           text.AddText(Form("OK (%d)",fIsOnlineThr));
384           text.SetFillColor(kGreen);
385           h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());      
386           linePadOccupMin->SetLineColor(kGreen);
387           linePadOccupMax->SetLineColor(kGreen);
388           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
389           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());         
390         }
391       else if ( badOccCnt == 1 )  
392         {
393           hmpQaFlags[1]  = AliQAv1::kWARNING;
394           text.Clear();
395           text.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
396           text.SetFillColor(kOrange);
397           h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());      
398           linePadOccupMin->SetLineColor(kGreen);
399           linePadOccupMax->SetLineColor(kGreen);
400           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
401           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());
402         }
403       else if (  badOccCnt == 2 )  
404         {
405           hmpQaFlags[1]  = AliQAv1::kERROR;
406           text.Clear();
407           text.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
408           text.SetFillColor(kRed);
409           h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());            
410           linePadOccupMin->SetLineColor(kGreen);
411           linePadOccupMax->SetLineColor(kGreen);
412           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
413           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());
414         }
415       else 
416         {
417           hmpQaFlags[1] = AliQAv1::kFATAL;
418           text.Clear();
419           text.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
420           text.SetFillColor(kRed);
421           h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());            
422           linePadOccupMin->SetLineColor(kGreen);
423           linePadOccupMax->SetLineColor(kGreen);
424           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
425           h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());   
426         }
427       }
428     }
429     }//___HmpPadOcc
430   
431   Int_t sumPadMapChROR[7]={0};
432   Int_t sumPadMapChROL[7]={0};
433   Int_t bigMapFlag =  AliQAv1::kINFO;
434   Int_t errCntBigMap = 0;
435
436   if( list->FindObject(Form("%s_hHmpBigMap",AliRecoParam::GetEventSpecieName(specie)) )) 
437     {
438       TH2* h2 = dynamic_cast<TH2*>( list->FindObject(Form("%s_hHmpBigMap",AliRecoParam::GetEventSpecieName(specie)))); 
439       if(h2) {
440         if( h2->Integral() > 1 ) {
441       // no entres -> fatal
442       if( h2->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
443        h2->SetStats(0);
444       // clean up the text, stat, lines, ...
445       if( h2->GetListOfFunctions() ) h2->GetListOfFunctions()->Clear();
446
447       //calculate missing pad fraction
448       for(Int_t ich = 0; ich < 7; ich++)
449         {
450           for(Int_t iy=1+ich*144;iy<=144+ich*144;iy++) {
451             for(Int_t ix=1;ix<=80;ix++)    if(h2->GetBinContent(ix,iy) > 0) sumPadMapChROL[ich]++;
452             for(Int_t ix=81;ix<=160;ix++)  if(h2->GetBinContent(ix,iy) > 0) sumPadMapChROR[ich]++;        
453           }
454         }//ch loop
455        //check the calculated missing pad fraction
456       for(Int_t ich = 0; ich < 7; ich++)
457         {
458           
459           bigMapFlag =  AliQAv1::kINFO; 
460           if( (1-sumPadMapChROL[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionWarningThreshold ||
461               (1-sumPadMapChROR[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionWarningThreshold ) 
462           {
463             bigMapFlag =  AliQAv1::kWARNING;
464             }
465           if( (1-sumPadMapChROL[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionErrorThreshold ||
466               (1-sumPadMapChROR[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionErrorThreshold ) {
467              bigMapFlag =  AliQAv1::kERROR; 
468              errCntBigMap++;
469            }      
470           
471         }//ch loop
472       if( errCntBigMap > 0 ) bigMapFlag =  AliQAv1::kERROR;
473
474       //update labels
475        if (  bigMapFlag == AliQAv1::kINFO )  
476         {
477           hmpQaFlags[2]  = AliQAv1::kINFO;
478           text.Clear();
479           text.AddText(Form("OK (%d)",fIsOnlineThr));
480           text.SetFillColor(kGreen);
481           h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());      
482         }
483        else if (  bigMapFlag == AliQAv1::kWARNING )  
484         {
485           hmpQaFlags[2]  = AliQAv1::kWARNING;
486           text.Clear();
487           text.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
488           text.SetFillColor(kOrange);
489           h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());      
490         }
491        else if (  bigMapFlag == AliQAv1::kERROR )  
492         {
493           hmpQaFlags[2]  = AliQAv1::kERROR;
494           text.Clear();
495           text.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
496           text.SetFillColor(kRed);
497           h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());      
498         }
499        else
500          {
501            hmpQaFlags[2]  = AliQAv1::kFATAL;
502            text.Clear();
503            text.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
504            text.SetFillColor(kRed);
505            h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());     
506          }      
507         }       
508      }
509     }//___HmpBigMap
510   
511
512   Int_t numSectorsMissing = 0, numSectorsGainLoss = 0;
513   Double_t  hSumQPerSector[42]={0};
514
515   if( list->FindObject(Form("%s_fHmpHvSectorQ",AliRecoParam::GetEventSpecieName(specie))))
516     {
517       
518       TH2* h2 = dynamic_cast<TH2*>( list->FindObject(Form("%s_fHmpHvSectorQ",AliRecoParam::GetEventSpecieName(specie)))); 
519       if(h2) {
520         if(h2->Integral() > 0 ) {
521       // no entres -> fatal
522       if( h2->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
523              h2->SetStats(0);
524       // clean up the text, stat, lines, ...
525       if( h2->GetListOfFunctions() ) h2->GetListOfFunctions()->Clear();
526      
527       //___ check sectors 
528       for(Int_t isec = 1 ; isec <= 42; isec++)
529         {
530           for(Int_t ibiny=100;ibiny<410;ibiny++) {hSumQPerSector[isec-1] += h2->GetBinContent(ibiny,isec); }
531           if(hSumQPerSector[isec-1]  < 0.001) {numSectorsGainLoss++;}  // there is no photon and mip peak, gain loss
532           if(h2->GetBinContent(1,isec) < 0.01 ) {numSectorsMissing++; } //practically there is no charge , the sector is missing
533         }
534       Int_t  sectorErrors = numSectorsGainLoss+numSectorsMissing;
535      
536       
537         if ( sectorErrors <=  3)
538           {
539             hmpQaFlags[3]  = AliQAv1::kINFO;
540             text.Clear();
541             text.AddText(Form("OK (%d)",fIsOnlineThr));
542             text.SetFillColor(kGreen);
543             h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());       
544           }
545         else if ( sectorErrors > fHmpQaThr_SectorGainLossWarningThreshold)
546           {
547             hmpQaFlags[3]  = AliQAv1::kWARNING;
548             text.Clear();
549             text.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
550             if(numSectorsMissing > 0 ) text.AddText(Form("MISSING SECTOR?"));              
551             text.SetFillColor(kOrange);
552             h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
553           }
554         else if ( sectorErrors > fHmpQaThr_SectorGainLossErrorThreshold)
555           {
556             hmpQaFlags[3]  = AliQAv1::kERROR;
557             text.Clear();
558             text.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
559             if(numSectorsMissing > 0 ) text.AddText(Form("MISSING SECTOR?"));
560             text.SetFillColor(kRed);
561             h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());    
562           }
563         else
564           {
565             hmpQaFlags[3]  = AliQAv1::kFATAL;
566             text.Clear();
567             text.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
568             text.SetFillColor(kRed);
569             h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());    
570           }
571         }
572     }
573   }
574   
575   
576   //del lines, ...
577   lineDdlDataSizeMin->Delete();
578   lineDdlDataSizeMax->Delete();
579   linePadOccupMin->Delete();
580   linePadOccupMax->Delete();
581   
582   Double_t dflag = -1;
583   switch ( TMath::MaxElement(4,hmpQaFlags))
584     {
585     case  AliQAv1::kINFO:
586       dflag = 1.0;
587       
588       break;
589     case AliQAv1::kWARNING:
590       dflag = 0.75;
591       
592       break;
593     case AliQAv1::kERROR:
594       dflag = 0.25;
595       
596       break;
597     case AliQAv1::kFATAL:
598       dflag = -1.0;
599      
600       break;
601     default:
602       dflag = AliQAv1::kNULLBit;
603      
604       break;
605     }     
606      
607   return dflag;
608   
609   
610   
611 }
612 //___________________________________________________________________________________________________
613 void AliHMPIDQAChecker::InitOnlineThresholds()
614 {
615   //
616   // Init the online thresholds from GRP generated by AMORE
617   //
618     
619   AliCDBManager* man = AliCDBManager::Instance(); 
620   if(!man)  {     fIsOnlineThr = kFALSE;     return;   }
621   
622 //  man->SetRun(AliQAChecker::Instance()->GetRunNumber()); //this should be the normal code
623  //  man->SetSpecificStorage("GRP/Calib/QAThresholds","local://$ALICE_ROOT/OCDB/");// need for development machine testing
624
625   if(!man->Get("GRP/Calib/QAThresholds"))  {     fIsOnlineThr = kFALSE;    return;    }
626  
627   AliCDBEntry* entry = man->Get("GRP/Calib/QAThresholds");
628   if(!entry)    {     fIsOnlineThr = kFALSE;    return;    }
629   
630   TObjArray* branch = (TObjArray*) entry->GetObject();
631   if(!branch ) {     fIsOnlineThr = kFALSE;     return;    }
632
633   AliQAThresholds* thresholds = (AliQAThresholds*) branch->FindObject("HMP");  
634   if(!thresholds) { fIsOnlineThr = kFALSE; return;   }
635   else 
636      fIsOnlineThr = kTRUE; 
637   
638   Int_t teCnt = 0;
639   TString parName = "zero";
640   while ( thresholds->GetThreshold(teCnt)) 
641   {
642      
643      parName = (thresholds->GetThreshold(teCnt))->GetName();
644
645           
646      if( parName.Contains("HmpNumberOfExcludedDDLthreshold") )
647      {
648        TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt); 
649        fHmpQaThr_NumberOfExcludedDDL = myParam->GetVal();
650        myParam->Delete();
651      }
652    
653      
654      if( parName.Contains("HmpDataSizeLowerThreshold") ) 
655      {
656        TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt); 
657        fHmpQaThr_DataSizeLowerThreshold = myParam->GetVal();
658        myParam->Delete();        
659      }
660              
661      if( parName.Contains("HmpDataSizeUpperThreshold") ) 
662      {
663        TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt); 
664        fHmpQaThr_DataSizeUpperThreshold = myParam->GetVal();
665        myParam->Delete();        
666      }
667      
668      
669      if( parName.Contains("HmpPadOccupancyLowerThreshold") ) 
670      {
671        TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt); 
672        fHmpQaThr_PadOccupancyLowerThreshold = myParam->GetVal();
673        myParam->Delete();        
674      }
675            
676      if( parName.Contains("HmpPadOccupancyUpperThreshold") ) 
677      {
678        TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt); 
679        fHmpQaThr_PadOccupancyUpperThreshold = myParam->GetVal();
680        myParam->Delete();        
681      }
682        
683              
684      if( parName.Contains("HmpSectorGainLossWarningThreshold") ) 
685       {
686        TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt); 
687        fHmpQaThr_SectorGainLossWarningThreshold = myParam->GetVal();
688        myParam->Delete();        
689      } 
690        
691      if( parName.Contains("HmpSectorGainLossErrorThreshold") ) 
692        {
693        TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt); 
694        fHmpQaThr_SectorGainLossErrorThreshold = myParam->GetVal();
695        myParam->Delete();        
696      } 
697        
698      if( parName.Contains("HmpMissingPadFractionWarningThreshold") ) 
699        {
700        TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt); 
701        fHmpQaThr_MissingPadFractionWarningThreshold = myParam->GetVal();
702        myParam->Delete();        
703      }
704      
705      if( parName.Contains("HmpMissingPadFractionErrorThreshold") ) 
706          {
707        TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt); 
708        fHmpQaThr_MissingPadFractionErrorThreshold = myParam->GetVal();
709        myParam->Delete();        
710      }
711       
712      teCnt++;    
713   }//while
714   
715   
716   // PrintThresholds();
717 }
718 //___________________________________________________________________________________________________
719 void AliHMPIDQAChecker::PrintThresholds()
720 {
721   Printf("--- Printing thresholds ---");  
722   Printf(Form("--- Default or online: %d ---",fIsOnlineThr));  
723   Printf(Form("--- fHmpQaThr_NumberOfExcludedDDL %d ---",fHmpQaThr_NumberOfExcludedDDL));
724   Printf(Form("--- fHmpQaThr_DataSizeLowerThreshold %d ---",fHmpQaThr_DataSizeLowerThreshold));
725   Printf(Form("--- fHmpQaThr_DataSizeUpperThreshold %d ---",fHmpQaThr_DataSizeUpperThreshold));
726   Printf(Form("--- fHmpQaThr_PadOccupancyLowerThreshold %f ---",fHmpQaThr_PadOccupancyLowerThreshold));
727   Printf(Form("--- fHmpQaThr_PadOccupancyUpperThreshold %f ---",fHmpQaThr_PadOccupancyUpperThreshold));
728   Printf(Form("--- fHmpQaThr_SectorGainLossWarningThreshold %d ---",fHmpQaThr_SectorGainLossWarningThreshold));
729   Printf(Form("--- fHmpQaThr_SectorGainLossErrorThreshold %d ---",fHmpQaThr_SectorGainLossErrorThreshold));
730   Printf(Form("--- fHmpQaThr_MissingPadFractionWarningThreshold %f ---",fHmpQaThr_MissingPadFractionWarningThreshold));
731   Printf(Form("--- fHmpQaThr_MissingPadFractionErrorThreshold %f ---",fHmpQaThr_MissingPadFractionErrorThreshold));
732   Printf("--- Printing thresholds done ---");  
733
734   
735   
736   
737 }
738 //___________________________________________________________________________________________________
739
740