]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDQAChecker.cxx
Changes for #89427: Porting modification for material budget issues to the Release...
[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)
34fbdb5b 45 //_________________________________________________________________
46AliHMPIDQAChecker::AliHMPIDQAChecker() :
47AliQACheckerBase("HMPID","HMPID Quality Assurance Data Checker"),
09b0a300 48fNoReference(kTRUE),
49fQARefRec(NULL)
34fbdb5b 50{
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);
54 if(QARefRec) {
55 fQARefRec = dynamic_cast<TObjArray*> (QARefRec->GetObject()) ;
56 if (fQARefRec)
57 if (fQARefRec->GetEntries())
58 fNoReference = kFALSE ;
59 if (fNoReference)
60 AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker!");
61 }
62}
63
64//_________________________________________________________________
65AliHMPIDQAChecker::AliHMPIDQAChecker(const AliHMPIDQAChecker& qac) :
66AliQACheckerBase(qac.GetName(), qac.GetTitle()),
67fNoReference(qac.fNoReference),
68fQARefRec(NULL)
69{
70 fNoReference = qac.fNoReference ;
71 if (qac.fQARefRec) {
72 fQARefRec = new TObjArray(qac.fQARefRec->GetEntries()) ;
73 for (Int_t index=0; index < qac.fQARefRec->GetEntries(); index++)
74 fQARefRec->Add(qac.fQARefRec->At(index)) ;
75 }
76}
e4e6ae61 77
34fbdb5b 78//_________________________________________________________________
79AliHMPIDQAChecker::~AliHMPIDQAChecker()
80{
71edfefd 81 if(fQARefRec) { fQARefRec->Delete() ; delete fQARefRec ; }
34fbdb5b 82}
fdd7f404 83//_________________________________________________________________
86017bd8 84void AliHMPIDQAChecker::Check(Double_t * check, AliQAv1::ALITASK_t index, TObjArray ** list, const AliDetectorRecoParam * /*recoParam*/)
fdd7f404 85{
86//
87// Main check function: Depending on the TASK, different checks are applied
88// At the moment: check for empty histograms and checks for RecPoints
89
34fbdb5b 90 if(fNoReference)
fdd7f404 91
57acd2d2 92 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
93 check[specie] = 1.0;
4e25ac79 94 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
57acd2d2 95 continue ;
2ac899f2 96 // checking for empy histograms
57acd2d2 97 if(CheckEntries(list[specie]) == 0) {
98 AliWarning("histograms are empty");
99 check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
100 }
101
2ac899f2 102 //check sim
34fbdb5b 103 if(index == AliQAv1::kSIM) check[specie] = CheckSim(list[specie], fQARefRec);
2ac899f2 104
57acd2d2 105 // checking rec points
34fbdb5b 106 if(index == AliQAv1::kREC) check[specie] = CheckRec(list[specie], fQARefRec);
107
57acd2d2 108 //default check response. It will be changed when reasonable checks will be considered
109 else check[specie] = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run
2ac899f2 110 } // species loop
fdd7f404 111}
112//_________________________________________________________________
113Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
114{
115 //
116 // check on the QA histograms on the input list:
117 //
118
119 Double_t test = 0.0 ;
120 Int_t count = 0 ;
121
122 if (list->GetEntries() == 0){
123 test = 1. ; // nothing to check
124 }
125 else {
126 TIter next(list) ;
127 TH1 * hdata ;
128 count = 0 ;
129 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
130 if (hdata) {
131 Double_t rv = 0.;
132 //Printf("hitogram %s has entries: %f ",hdata->GetName(),hdata->GetEntries());
133 if(hdata->GetEntries()>0)rv=1;
134 count++ ;
135 test += rv ;
136 }
137 else{
138 AliError("Data type cannot be processed") ;
139 }
140
141 }
142 if (count != 0) {
143 if (test==0) {
144 AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
145 test = 0.; //upper limit value to set kWARNING flag for a task
146 }
147 else test = 1 ;
148 }
149 }
150
151 return test ;
152}
153//_________________________________________________________________
2ac899f2 154Double_t AliHMPIDQAChecker::CheckSim(TObjArray *listsim, TObjArray *listref) const
155{
156 //
157 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
158 //
159
160 Float_t checkresponse = 0;
161
162 Float_t counter = 0 ;
163 TIter next(listsim) ;
164 TH1* histo;
165 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
5089b327 166 //PH The histogram should have at least 10 bins with at least 5 entries
2ac899f2 167 Int_t nbinsabove = 0;
5089b327 168 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
169 if (histo->GetBinContent(ibin)>5) nbinsabove++;
2ac899f2 170 }
171
5089b327 172 if( nbinsabove < 10 ) counter++;
2ac899f2 173 else {
174 TString h = histo->GetTitle();
175 if(h.Contains("Zoom")){
5089b327 176 histo->Fit("expo","LQ0","",5,50);
2ac899f2 177 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
178 }
179 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
180 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
181 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
182 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
183 }
184 }
185 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listsim->GetEntries())
186
187 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
188 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
189 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
190 return checkresponse;
191}
fdd7f404 192
2ac899f2 193//___________________________________________________________________________________________________
194Double_t AliHMPIDQAChecker::CheckRec(TObjArray *listrec, TObjArray *listref) const
fdd7f404 195{
0f46609e 196 //
fbe0560f 197 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
198 //
fdd7f404 199
200 Float_t checkresponse = 0;
201
202 Float_t counter = 0 ;
203 TIter next(listrec) ;
204 TH1* histo;
205 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
5089b327 206 //PH The histogram should have at least 10 bins with at least 5 entries
3762c8fb 207 Int_t nbinsabove = 0;
5089b327 208 for (Int_t ibin=histo->FindBin(1); ibin<=histo->FindBin(50); ibin++) {
209 if (histo->GetBinContent(ibin)>5) nbinsabove++;
3762c8fb 210 }
211
5089b327 212 if( nbinsabove < 10 ) counter++;
fdd7f404 213 else {
214 TString h = histo->GetTitle();
215 if(h.Contains("Zoom")){
5089b327 216 histo->Fit("expo","LQ0","",5,50);
fdd7f404 217 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
218 }
219 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
220 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
221 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
5379c4a3 222 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
fdd7f404 223 }
224 }
225 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
226
227 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
228 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
229 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
230 return checkresponse;
231}
e4e6ae61 232