]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDQAChecker.cxx
Related to the report #60402: protection in case of timestamp=0 with warnings added
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQAChecker.cxx
CommitLineData
e4e6ae61 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>
fdd7f404 29#include <TF1.h>
e4e6ae61 30#include <TIterator.h>
31#include <TKey.h>
32#include <TFile.h>
33
34// --- Standard library ---
35
36// --- AliRoot header files ---
37#include "AliLog.h"
4e25ac79 38#include "AliQAv1.h"
b990acb3 39#include "AliQAChecker.h"
40#include "AliHMPIDQAChecker.h"
fdd7f404 41#include "AliCDBEntry.h"
fc7e0df2 42#include "AliQAManager.h"
e4e6ae61 43
b990acb3 44ClassImp(AliHMPIDQAChecker)
e4e6ae61 45
fdd7f404 46//_________________________________________________________________
486788fc 47Double_t * AliHMPIDQAChecker::Check(AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
fdd7f404 48{
49//
50// Main check function: Depending on the TASK, different checks are applied
51// At the moment: check for empty histograms and checks for RecPoints
52
57acd2d2 53 Double_t * check = new Double_t[AliRecoParam::kNSpecies] ;
54
4e25ac79 55 char * detOCDBDir = Form("HMPID/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
fc7e0df2 56 AliCDBEntry *QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
2ac899f2 57 if(!QARefRec){
57acd2d2 58 AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker ...exiting");
59 return check;
fdd7f404 60 }
61
57acd2d2 62 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
63 check[specie] = 1.0;
4e25ac79 64 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
57acd2d2 65 continue ;
2ac899f2 66 // checking for empy histograms
57acd2d2 67 if(CheckEntries(list[specie]) == 0) {
68 AliWarning("histograms are empty");
69 check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
70 }
71
2ac899f2 72 //check sim
73 if(index == AliQAv1::kSIM) check[specie] = CheckSim(list[specie],(TObjArray *)QARefRec->GetObject());
74
57acd2d2 75 // checking rec points
2ac899f2 76 if(index == AliQAv1::kREC) check[specie] = CheckRec(list[specie],(TObjArray *)QARefRec->GetObject());
fdd7f404 77
57acd2d2 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
2ac899f2 80 } // species loop
81
fdd7f404 82 return check;
83
84}
85//_________________________________________________________________
86Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
87{
88 //
89 // check on the QA histograms on the input list:
90 //
91
92 Double_t test = 0.0 ;
93 Int_t count = 0 ;
94
95 if (list->GetEntries() == 0){
96 test = 1. ; // nothing to check
97 }
98 else {
99 TIter next(list) ;
100 TH1 * hdata ;
101 count = 0 ;
102 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
103 if (hdata) {
104 Double_t rv = 0.;
105 //Printf("hitogram %s has entries: %f ",hdata->GetName(),hdata->GetEntries());
106 if(hdata->GetEntries()>0)rv=1;
107 count++ ;
108 test += rv ;
109 }
110 else{
111 AliError("Data type cannot be processed") ;
112 }
113
114 }
115 if (count != 0) {
116 if (test==0) {
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
119 }
120 else test = 1 ;
121 }
122 }
123
124 return test ;
125}
126//_________________________________________________________________
2ac899f2 127Double_t AliHMPIDQAChecker::CheckSim(TObjArray *listsim, TObjArray *listref) const
128{
129 //
130 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
131 //
132
133 Float_t checkresponse = 0;
134
135 Float_t counter = 0 ;
136 TIter next(listsim) ;
137 TH1* histo;
138 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
5089b327 139 //PH The histogram should have at least 10 bins with at least 5 entries
2ac899f2 140 Int_t nbinsabove = 0;
5089b327 141 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
142 if (histo->GetBinContent(ibin)>5) nbinsabove++;
2ac899f2 143 }
144
5089b327 145 if( nbinsabove < 10 ) counter++;
2ac899f2 146 else {
147 TString h = histo->GetTitle();
148 if(h.Contains("Zoom")){
5089b327 149 histo->Fit("expo","LQ0","",5,50);
2ac899f2 150 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
151 }
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))));
156 }
157 }
158 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listsim->GetEntries())
159
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;
164}
fdd7f404 165
2ac899f2 166//___________________________________________________________________________________________________
167Double_t AliHMPIDQAChecker::CheckRec(TObjArray *listrec, TObjArray *listref) const
fdd7f404 168{
0f46609e 169 //
fbe0560f 170 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
171 //
fdd7f404 172
173 Float_t checkresponse = 0;
174
175 Float_t counter = 0 ;
176 TIter next(listrec) ;
177 TH1* histo;
178 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
5089b327 179 //PH The histogram should have at least 10 bins with at least 5 entries
3762c8fb 180 Int_t nbinsabove = 0;
5089b327 181 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
182 if (histo->GetBinContent(ibin)>5) nbinsabove++;
3762c8fb 183 }
184
5089b327 185 if( nbinsabove < 10 ) counter++;
fdd7f404 186 else {
187 TString h = histo->GetTitle();
188 if(h.Contains("Zoom")){
5089b327 189 histo->Fit("expo","LQ0","",5,50);
fdd7f404 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++;
5379c4a3 195 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
fdd7f404 196 }
197 }
198 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->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}
e4e6ae61 205