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)
46 //_________________________________________________________________
47 Double_t * AliHMPIDQAChecker::Check(AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
50 // Main check function: Depending on the TASK, different checks are applied
51 // At the moment: check for empty histograms and checks for RecPoints
53 Double_t * check = new Double_t[AliRecoParam::kNSpecies] ;
55 char * detOCDBDir = Form("HMPID/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
56 AliCDBEntry *QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
58 AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker ...exiting");
62 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
64 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
66 // checking for empy histograms
67 if(CheckEntries(list[specie]) == 0) {
68 AliWarning("histograms are empty");
69 check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
73 if(index == AliQAv1::kSIM) check[specie] = CheckSim(list[specie],(TObjArray *)QARefRec->GetObject());
75 // checking rec points
76 if(index == AliQAv1::kREC) check[specie] = CheckRec(list[specie],(TObjArray *)QARefRec->GetObject());
78 //default check response. It will be changed when reasonable checks will be considered
79 else check[specie] = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run
85 //_________________________________________________________________
86 Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
89 // check on the QA histograms on the input list:
95 if (list->GetEntries() == 0){
96 test = 1. ; // nothing to check
102 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
105 //Printf("hitogram %s has entries: %f ",hdata->GetName(),hdata->GetEntries());
106 if(hdata->GetEntries()>0)rv=1;
111 AliError("Data type cannot be processed") ;
117 AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
118 test = 0.; //upper limit value to set kWARNING flag for a task
126 //_________________________________________________________________
127 Double_t AliHMPIDQAChecker::CheckSim(TObjArray *listsim, TObjArray *listref) const
130 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
133 Float_t checkresponse = 0;
135 Float_t counter = 0 ;
136 TIter next(listsim) ;
138 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
139 //PH The histogram should have at least 3 bins with entries
140 Int_t nbinsabove = 0;
141 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) { //1,50 is the fit region, see histo->Fit("expo","Q0","",1,50);
142 if (histo->GetBinContent(ibin)>0) nbinsabove++;
145 if( nbinsabove < 3 ) counter++;
147 TString h = histo->GetTitle();
148 if(h.Contains("Zoom")){
149 histo->Fit("expo","LQ0","",1,50);
150 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
152 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
153 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
154 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
155 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
158 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listsim->GetEntries())
160 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
161 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
162 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
163 return checkresponse;
166 //___________________________________________________________________________________________________
167 Double_t AliHMPIDQAChecker::CheckRec(TObjArray *listrec, TObjArray *listref) const
170 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
173 Float_t checkresponse = 0;
175 Float_t counter = 0 ;
176 TIter next(listrec) ;
178 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
179 //PH The histogram should have at least 3 bins with entries
180 Int_t nbinsabove = 0;
181 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) { //1,50 is the fit region, see histo->Fit("expo","Q0","",1,50);
182 if (histo->GetBinContent(ibin)>0) nbinsabove++;
185 if( nbinsabove < 3 ) counter++;
187 TString h = histo->GetTitle();
188 if(h.Contains("Zoom")){
189 histo->Fit("expo","LQ0","",1,50);
190 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
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))));
198 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
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;