1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
20 // Checks the quality assurance.
21 // By comparing with reference data
25 // --- ROOT system ---
31 #include <TIterator.h>
35 #include <TParameter.h>
36 #include <TPaveText.h>
37 // --- Standard library ---
39 // --- AliRoot header files ---
42 #include "AliQAChecker.h"
43 #include "AliHMPIDQAChecker.h"
44 #include "AliCDBEntry.h"
45 #include "AliQAManager.h"
46 #include "AliQAThresholds.h"
48 ClassImp(AliHMPIDQAChecker)
49 //_________________________________________________________________
50 AliHMPIDQAChecker::AliHMPIDQAChecker() :
51 AliQACheckerBase("HMPID","HMPID Quality Assurance Data Checker"),
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),
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);
73 fQARefRec = dynamic_cast<TObjArray*> (QARefRec->GetObject()) ;
75 if (fQARefRec->GetEntries())
76 fNoReference = kFALSE ;
78 AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker!");
83 //_________________________________________________________________
84 AliHMPIDQAChecker::~AliHMPIDQAChecker()
86 if(fQARefRec) { fQARefRec->Delete() ; delete fQARefRec ; }
88 //_________________________________________________________________
89 void AliHMPIDQAChecker::Check(Double_t * check, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
92 // Main check function: Depending on the TASK, different checks are applied
93 // At the moment: check for empty histograms and checks for RecPoints
95 InitOnlineThresholds();
99 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
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
109 check[specie] = AliQAv1::kINFO ;
113 if(index == AliQAv1::kSIM) check[specie] = CheckSim(list[specie], fQARefRec);
115 // checking rec points
116 if(index == AliQAv1::kREC) check[specie] = CheckRec(list[specie], fQARefRec);
119 if(index == AliQAv1::kRAW) { check[specie] = CheckRaw(specie,list[specie]); }
124 //_________________________________________________________________
125 Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
128 // check on the QA histograms on the input list:
131 Double_t test = 0.0 ;
134 if (list->GetEntries() == 0){
135 test = 1. ; // nothing to check
141 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
144 if(hdata->GetEntries()>0)rv=1;
149 AliError("Data type cannot be processed") ;
155 AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
156 test = 0.; //upper limit value to set kWARNING flag for a task
164 //_________________________________________________________________
165 Double_t AliHMPIDQAChecker::CheckSim(TObjArray *listsim, TObjArray *listref) const
168 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
171 Float_t checkresponse = 0;
173 Float_t counter = 0 ;
174 TIter next(listsim) ;
176 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
177 //PH The histogram should have at least 10 bins with at least 5 entries
178 Int_t nbinsabove = 0;
179 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
180 if (histo->GetBinContent(ibin)>5) nbinsabove++;
183 if( nbinsabove < 10 ) counter++;
185 TString h = histo->GetTitle();
186 if(h.Contains("Zoom")){
187 histo->Fit("expo","LQ0","",5,50);
188 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
190 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
191 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
192 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
193 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
196 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listsim->GetEntries())
198 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
199 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
200 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
201 return checkresponse;
204 //___________________________________________________________________________________________________
205 Double_t AliHMPIDQAChecker::CheckRec(TObjArray *listrec, TObjArray *listref) const
208 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
211 Float_t checkresponse = 0;
213 Float_t counter = 0 ;
214 TIter next(listrec) ;
216 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
217 //PH The histogram should have at least 10 bins with at least 5 entries
218 Int_t nbinsabove = 0;
219 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
220 if (histo->GetBinContent(ibin)>5) nbinsabove++;
223 if( nbinsabove < 10 ) counter++;
225 TString h = histo->GetTitle();
226 if(h.Contains("Zoom")){
227 histo->Fit("expo","LQ0","",5,50);
228 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
230 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
231 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
232 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
233 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
236 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
238 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
239 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
240 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
241 return checkresponse;
243 //___________________________________________________________________________________________________
244 Double_t AliHMPIDQAChecker::CheckRaw(Int_t specie, TObjArray* list)
247 // Check the raw data for offline / online using default or updated thresholds from AMORE
248 // 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
249 // in AMORE. But we can pu undividual labels.
252 //Int_t raqQualFlag = AliQAv1::kNULLBit;
254 Int_t hmpQaFlags[4]={-1}; //init for the 4 shifter histos
256 TString histname ="null";
257 TPaveText text(0.65,0.8,0.9,0.99,"NDC");
258 TPaveText text1(0.65,0.8,0.9,0.99,"NDC");
261 Int_t entries = list->GetEntriesFast();
263 if ( entries == 0 ) {
264 AliWarning(Form("HMPID QA Checker RAWS: no object to analyse! Exiting..."));
265 return AliQAv1::kFATAL;
268 TLine* lineDdlDataSizeMin = new TLine(1536,fHmpQaThr_DataSizeUpperThreshold,1548,fHmpQaThr_DataSizeUpperThreshold);
269 TLine* lineDdlDataSizeMax = new TLine(1536,fHmpQaThr_DataSizeLowerThreshold,1548,fHmpQaThr_DataSizeLowerThreshold);
271 TLine* linePadOccupMin = new TLine(1536,fHmpQaThr_PadOccupancyUpperThreshold,1548,fHmpQaThr_PadOccupancyUpperThreshold);
272 TLine* linePadOccupMax = new TLine(1536,fHmpQaThr_PadOccupancyLowerThreshold,1548,fHmpQaThr_PadOccupancyLowerThreshold);
275 Int_t badDdlCnt = 0, badOccCnt = 0;
277 //___ check data size per ddl
278 if( list->FindObject(Form("%s_hHmpDdlDataSize",AliRecoParam::GetEventSpecieName(specie)) ))
280 TH1* h1 = dynamic_cast<TH1*>( list->FindObject(Form("%s_hHmpDdlDataSize",AliRecoParam::GetEventSpecieName(specie))));
282 if( h1->Integral() > 1 ) {
283 // no entres -> fatal
284 //if( h1->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
288 // clean up the text, stat, lines, ...
289 if( h1->GetListOfFunctions() ) h1->GetListOfFunctions()->Clear();
291 for ( Int_t iddl = 1 ; iddl <= 14; iddl++)
293 if( h1->GetBinContent(iddl) < fHmpQaThr_DataSizeLowerThreshold || h1->GetBinContent(iddl) > fHmpQaThr_DataSizeUpperThreshold ) badDdlCnt++;
295 //___ check if one or more DDLs are excluded
298 badDdlCnt -= fHmpQaThr_NumberOfExcludedDDL;
302 //___ check how many are bad
303 if ( badDdlCnt == 0 )
305 hmpQaFlags[0] = AliQAv1::kINFO;
307 text1.AddText(Form("OK (%d)",fIsOnlineThr));
308 text1.SetFillColor(kGreen);
309 h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());
310 lineDdlDataSizeMin->SetLineColor(kGreen);
311 lineDdlDataSizeMax->SetLineColor(kGreen);
312 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
313 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());
315 else if ( badDdlCnt == 1 )
317 hmpQaFlags[0] = AliQAv1::kWARNING;
319 text1.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
320 text1.SetFillColor(kOrange);
321 h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());
322 lineDdlDataSizeMin->SetLineColor(kOrange);
323 lineDdlDataSizeMax->SetLineColor(kOrange);
324 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
325 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());
327 else if ( badDdlCnt >= 2 )
329 hmpQaFlags[0] = AliQAv1::kERROR;
331 text1.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
332 text1.SetFillColor(kRed);
333 h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());
334 lineDdlDataSizeMin->SetLineColor(kRed);
335 lineDdlDataSizeMax->SetLineColor(kRed);
336 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
337 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());
341 hmpQaFlags[0] = AliQAv1::kFATAL;
343 text1.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
344 text1.SetFillColor(kRed);
345 h1->GetListOfFunctions()->Add((TPaveText*)text1.Clone());
346 lineDdlDataSizeMin->SetLineColor(kRed);
347 lineDdlDataSizeMax->SetLineColor(kRed);
348 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMin->Clone());
349 h1->GetListOfFunctions()->Add((TLine*)lineDdlDataSizeMax->Clone());
352 }//the histo is filled
353 }//___hHmpDdlDataSize
357 if( list->FindObject(Form("%s_fHmpPadOcc",AliRecoParam::GetEventSpecieName(specie)) ))
361 TH1* h1 = dynamic_cast<TH1*>( list->FindObject(Form("%s_fHmpPadOcc",AliRecoParam::GetEventSpecieName(specie))));
363 if( h1->Integral() > 1 ) {
364 // no entres -> fatal
365 //if( h1->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
367 // clean up the text, stat, lines, ...
368 if( h1->GetListOfFunctions() ) h1->GetListOfFunctions()->Clear();
370 for ( Int_t iddl = 1 ; iddl <= 14; iddl++)
372 if( h1->GetBinContent(iddl) < fHmpQaThr_PadOccupancyLowerThreshold || h1->GetBinContent(iddl) > fHmpQaThr_PadOccupancyUpperThreshold ) badOccCnt++;
375 badOccCnt -= fHmpQaThr_NumberOfExcludedDDL;
377 //___ check how many are bad
378 if ( badOccCnt == 0 )
380 hmpQaFlags[1] = AliQAv1::kINFO;
382 text.AddText(Form("OK (%d)",fIsOnlineThr));
383 text.SetFillColor(kGreen);
384 h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());
385 linePadOccupMin->SetLineColor(kGreen);
386 linePadOccupMax->SetLineColor(kGreen);
387 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
388 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());
390 else if ( badOccCnt == 1 )
392 hmpQaFlags[1] = AliQAv1::kWARNING;
394 text.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
395 text.SetFillColor(kOrange);
396 h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());
397 linePadOccupMin->SetLineColor(kGreen);
398 linePadOccupMax->SetLineColor(kGreen);
399 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
400 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());
402 else if ( badOccCnt == 2 )
404 hmpQaFlags[1] = AliQAv1::kERROR;
406 text.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
407 text.SetFillColor(kRed);
408 h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());
409 linePadOccupMin->SetLineColor(kGreen);
410 linePadOccupMax->SetLineColor(kGreen);
411 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
412 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());
416 hmpQaFlags[1] = AliQAv1::kFATAL;
418 text.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
419 text.SetFillColor(kRed);
420 h1->GetListOfFunctions()->Add((TPaveText*)text.Clone());
421 linePadOccupMin->SetLineColor(kGreen);
422 linePadOccupMax->SetLineColor(kGreen);
423 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMin->Clone());
424 h1->GetListOfFunctions()->Add((TLine*)linePadOccupMax->Clone());
430 Int_t sumPadMapChROR[7]={0};
431 Int_t sumPadMapChROL[7]={0};
432 Int_t bigMapFlag = AliQAv1::kINFO;
433 Int_t errCntBigMap = 0;
435 if( list->FindObject(Form("%s_hHmpBigMap",AliRecoParam::GetEventSpecieName(specie)) ))
437 TH2* h2 = dynamic_cast<TH2*>( list->FindObject(Form("%s_hHmpBigMap",AliRecoParam::GetEventSpecieName(specie))));
439 if( h2->Integral() > 1 ) {
440 // no entres -> fatal
441 // if( h2->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
443 // clean up the text, stat, lines, ...
444 if( h2->GetListOfFunctions() ) h2->GetListOfFunctions()->Clear();
446 //calculate missing pad fraction
447 for(Int_t ich = 0; ich < 7; ich++)
449 for(Int_t iy=1+ich*144;iy<=144+ich*144;iy++) {
450 for(Int_t ix=1;ix<=80;ix++) if(h2->GetBinContent(ix,iy) > 0) sumPadMapChROL[ich]++;
451 for(Int_t ix=81;ix<=160;ix++) if(h2->GetBinContent(ix,iy) > 0) sumPadMapChROR[ich]++;
454 //check the calculated missing pad fraction
455 for(Int_t ich = 0; ich < 7; ich++)
458 bigMapFlag = AliQAv1::kINFO;
459 if( (1-sumPadMapChROL[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionWarningThreshold ||
460 (1-sumPadMapChROR[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionWarningThreshold )
462 bigMapFlag = AliQAv1::kWARNING;
464 if( (1-sumPadMapChROL[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionErrorThreshold ||
465 (1-sumPadMapChROR[ich]/1.0/11520) > fHmpQaThr_MissingPadFractionErrorThreshold ) {
466 bigMapFlag = AliQAv1::kERROR;
471 if( errCntBigMap > 0 ) bigMapFlag = AliQAv1::kERROR;
474 if ( bigMapFlag == AliQAv1::kINFO )
476 hmpQaFlags[2] = AliQAv1::kINFO;
478 text.AddText(Form("OK (%d)",fIsOnlineThr));
479 text.SetFillColor(kGreen);
480 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
482 else if ( bigMapFlag == AliQAv1::kWARNING )
484 hmpQaFlags[2] = AliQAv1::kWARNING;
486 text.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
487 text.SetFillColor(kOrange);
488 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
490 else if ( bigMapFlag == AliQAv1::kERROR )
492 hmpQaFlags[2] = AliQAv1::kERROR;
494 text.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
495 text.SetFillColor(kRed);
496 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
500 hmpQaFlags[2] = AliQAv1::kFATAL;
502 text.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
503 text.SetFillColor(kRed);
504 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
511 Int_t numSectorsMissing = 0, numSectorsGainLoss = 0;
512 Double_t hSumQPerSector[42]={0};
514 if( list->FindObject(Form("%s_fHmpHvSectorQ",AliRecoParam::GetEventSpecieName(specie))))
517 TH2* h2 = dynamic_cast<TH2*>( list->FindObject(Form("%s_fHmpHvSectorQ",AliRecoParam::GetEventSpecieName(specie))));
519 if(h2->Integral() > 0 ) {
520 // no entres -> fatal
521 // if( h2->GetEntries() == 0) {raqQualFlag = AliQAv1::kFATAL;}
523 // clean up the text, stat, lines, ...
524 if( h2->GetListOfFunctions() ) h2->GetListOfFunctions()->Clear();
527 for(Int_t isec = 1 ; isec <= 42; isec++)
529 for(Int_t ibiny=100;ibiny<410;ibiny++) {hSumQPerSector[isec-1] += h2->GetBinContent(ibiny,isec); }
530 if(hSumQPerSector[isec-1] < 0.001) {numSectorsGainLoss++;} // there is no photon and mip peak, gain loss
531 if(h2->GetBinContent(1,isec) < 0.01 ) {numSectorsMissing++; } //practically there is no charge , the sector is missing
533 Int_t sectorErrors = numSectorsGainLoss+numSectorsMissing;
536 if ( sectorErrors <= 3)
538 hmpQaFlags[3] = AliQAv1::kINFO;
540 text.AddText(Form("OK (%d)",fIsOnlineThr));
541 text.SetFillColor(kGreen);
542 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
544 else if ( sectorErrors > fHmpQaThr_SectorGainLossWarningThreshold)
546 hmpQaFlags[3] = AliQAv1::kWARNING;
548 text.AddText(Form("WARNING CHECK TWIKI (%d)",fIsOnlineThr));
549 if(numSectorsMissing > 0 ) text.AddText(Form("MISSING SECTOR?"));
550 text.SetFillColor(kOrange);
551 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
553 else if ( sectorErrors > fHmpQaThr_SectorGainLossErrorThreshold)
555 hmpQaFlags[3] = AliQAv1::kERROR;
557 text.AddText(Form("ERROR CALL ONCALL (%d)",fIsOnlineThr));
558 if(numSectorsMissing > 0 ) text.AddText(Form("MISSING SECTOR?"));
559 text.SetFillColor(kRed);
560 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
564 hmpQaFlags[3] = AliQAv1::kFATAL;
566 text.AddText(Form("FATAL CALL ONCALL (%d)",fIsOnlineThr));
567 text.SetFillColor(kRed);
568 h2->GetListOfFunctions()->Add((TPaveText*)text.Clone());
576 lineDdlDataSizeMin->Delete();
577 lineDdlDataSizeMax->Delete();
578 linePadOccupMin->Delete();
579 linePadOccupMax->Delete();
582 switch ( TMath::MaxElement(4,hmpQaFlags))
588 case AliQAv1::kWARNING:
592 case AliQAv1::kERROR:
596 case AliQAv1::kFATAL:
601 dflag = AliQAv1::kNULLBit;
611 //___________________________________________________________________________________________________
612 void AliHMPIDQAChecker::InitOnlineThresholds()
615 // Init the online thresholds from GRP generated by AMORE
618 AliCDBManager* man = AliCDBManager::Instance();
619 if(!man) { fIsOnlineThr = kFALSE; return; }
622 if(!man->Get("GRP/Calib/QAThresholds")) { fIsOnlineThr = kFALSE; return; }
624 AliCDBEntry* entry = man->Get("GRP/Calib/QAThresholds");
625 if(!entry) { fIsOnlineThr = kFALSE; return; }
627 TObjArray* branch = (TObjArray*) entry->GetObject();
628 if(!branch ) { fIsOnlineThr = kFALSE; return; }
630 AliQAThresholds* thresholds = (AliQAThresholds*) branch->FindObject("HMP");
631 if(!thresholds) { fIsOnlineThr = kFALSE; return; }
633 fIsOnlineThr = kTRUE;
637 TString parName = "zero";
638 while ( thresholds->GetThreshold(teCnt))
640 if(!((thresholds->GetThreshold(teCnt))->GetName())) return;
642 parName = (thresholds->GetThreshold(teCnt))->GetName();
644 if( parName.Contains("HmpNumberOfExcludedDDLthreshold") )
646 TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt);
647 fHmpQaThr_NumberOfExcludedDDL = myParam->GetVal();
651 if( parName.Contains("HmpDataSizeLowerThreshold") )
653 TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt);
654 fHmpQaThr_DataSizeLowerThreshold = myParam->GetVal();
657 if( parName.Contains("HmpDataSizeUpperThreshold") )
659 TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt);
660 fHmpQaThr_DataSizeUpperThreshold = myParam->GetVal();
663 if( parName.Contains("HmpPadOccupancyLowerThreshold") )
665 TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt);
666 fHmpQaThr_PadOccupancyLowerThreshold = myParam->GetVal();
669 if( parName.Contains("HmpPadOccupancyUpperThreshold") )
671 TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt);
672 fHmpQaThr_PadOccupancyUpperThreshold = myParam->GetVal();
675 if( parName.Contains("HmpSectorGainLossWarningThreshold") )
677 TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt);
678 fHmpQaThr_SectorGainLossWarningThreshold = myParam->GetVal();
681 if( parName.Contains("HmpSectorGainLossErrorThreshold") )
683 TParameter<int>* myParam = (TParameter<int>*) thresholds->GetThreshold(teCnt);
684 fHmpQaThr_SectorGainLossErrorThreshold = myParam->GetVal();
687 if( parName.Contains("HmpMissingPadFractionWarningThreshold") )
689 TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt);
690 fHmpQaThr_MissingPadFractionWarningThreshold = myParam->GetVal();
693 if( parName.Contains("HmpMissingPadFractionErrorThreshold") )
695 TParameter<float>* myParam = (TParameter<float>*) thresholds->GetThreshold(teCnt);
696 fHmpQaThr_MissingPadFractionErrorThreshold = myParam->GetVal();
703 // PrintThresholds();
705 //___________________________________________________________________________________________________
706 void AliHMPIDQAChecker::PrintThresholds()
708 Printf("--- Printing thresholds ---");
709 Printf("--- Default or online: %i ---",fIsOnlineThr);
710 Printf("--- fHmpQaThr_NumberOfExcludedDDL: %i ---",fHmpQaThr_NumberOfExcludedDDL);
711 Printf("--- fHmpQaThr_DataSizeLowerThreshold: %i ---",fHmpQaThr_DataSizeLowerThreshold);
712 Printf("--- fHmpQaThr_DataSizeUpperThreshold: %i ---",fHmpQaThr_DataSizeUpperThreshold);
713 Printf("--- fHmpQaThr_PadOccupancyLowerThreshold: %f ---",fHmpQaThr_PadOccupancyLowerThreshold);
714 Printf("--- fHmpQaThr_PadOccupancyUpperThreshold: %f ---",fHmpQaThr_PadOccupancyUpperThreshold);
715 Printf("--- fHmpQaThr_SectorGainLossWarningThreshold: %i ---",fHmpQaThr_SectorGainLossWarningThreshold);
716 Printf("--- fHmpQaThr_SectorGainLossErrorThreshold: %i ---",fHmpQaThr_SectorGainLossErrorThreshold);
717 Printf("--- fHmpQaThr_MissingPadFractionWarningThreshold: %f ---",fHmpQaThr_MissingPadFractionWarningThreshold);
718 Printf("--- fHmpQaThr_MissingPadFractionErrorThreshold: %f ---",fHmpQaThr_MissingPadFractionErrorThreshold);
719 Printf("--- Printing thresholds done ---");
725 //___________________________________________________________________________________________________