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 ---
30 #include <TIterator.h>
34 // --- Standard library ---
36 // --- AliRoot header files ---
39 #include "AliQAChecker.h"
40 #include "AliHMPIDQAChecker.h"
41 #include "AliCDBEntry.h"
42 #include "AliQAManager.h"
44 ClassImp(AliHMPIDQAChecker)
45 //_________________________________________________________________
46 AliHMPIDQAChecker::AliHMPIDQAChecker() :
47 AliQACheckerBase("HMPID","HMPID Quality Assurance Data Checker"),
51 //ctor, fetches the reference data from OCDB
52 char * detOCDBDir = Form("HMPID/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
53 AliCDBEntry * QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
55 fQARefRec = dynamic_cast<TObjArray*> (QARefRec->GetObject()) ;
57 if (fQARefRec->GetEntries())
58 fNoReference = kFALSE ;
60 AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker!");
64 //_________________________________________________________________
65 AliHMPIDQAChecker::~AliHMPIDQAChecker()
67 if(fQARefRec) { fQARefRec->Delete() ; delete fQARefRec ; }
69 //_________________________________________________________________
70 void AliHMPIDQAChecker::Check(Double_t * check, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
73 // Main check function: Depending on the TASK, different checks are applied
74 // At the moment: check for empty histograms and checks for RecPoints
78 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
80 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
82 // checking for empy histograms
83 if(CheckEntries(list[specie]) == 0) {
84 AliWarning("histograms are empty");
85 check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
89 if(index == AliQAv1::kSIM) check[specie] = CheckSim(list[specie], fQARefRec);
91 // checking rec points
92 if(index == AliQAv1::kREC) check[specie] = CheckRec(list[specie], fQARefRec);
94 //default check response. It will be changed when reasonable checks will be considered
95 else check[specie] = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run
98 //_________________________________________________________________
99 Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
102 // check on the QA histograms on the input list:
105 Double_t test = 0.0 ;
108 if (list->GetEntries() == 0){
109 test = 1. ; // nothing to check
115 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
118 //Printf("hitogram %s has entries: %f ",hdata->GetName(),hdata->GetEntries());
119 if(hdata->GetEntries()>0)rv=1;
124 AliError("Data type cannot be processed") ;
130 AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
131 test = 0.; //upper limit value to set kWARNING flag for a task
139 //_________________________________________________________________
140 Double_t AliHMPIDQAChecker::CheckSim(TObjArray *listsim, TObjArray *listref) const
143 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
146 Float_t checkresponse = 0;
148 Float_t counter = 0 ;
149 TIter next(listsim) ;
151 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
152 //PH The histogram should have at least 10 bins with at least 5 entries
153 Int_t nbinsabove = 0;
154 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
155 if (histo->GetBinContent(ibin)>5) nbinsabove++;
158 if( nbinsabove < 10 ) counter++;
160 TString h = histo->GetTitle();
161 if(h.Contains("Zoom")){
162 histo->Fit("expo","LQ0","",5,50);
163 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
165 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
166 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
167 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
168 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
171 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listsim->GetEntries())
173 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
174 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
175 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
176 return checkresponse;
179 //___________________________________________________________________________________________________
180 Double_t AliHMPIDQAChecker::CheckRec(TObjArray *listrec, TObjArray *listref) const
183 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
186 Float_t checkresponse = 0;
188 Float_t counter = 0 ;
189 TIter next(listrec) ;
191 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
192 //PH The histogram should have at least 10 bins with at least 5 entries
193 Int_t nbinsabove = 0;
194 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
195 if (histo->GetBinContent(ibin)>5) nbinsabove++;
198 if( nbinsabove < 10 ) counter++;
200 TString h = histo->GetTitle();
201 if(h.Contains("Zoom")){
202 histo->Fit("expo","LQ0","",5,50);
203 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
205 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
206 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
207 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
208 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
211 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
213 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
214 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
215 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
216 return checkresponse;