]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDQAChecker.cxx
In clustering code the method name is changed, Bug fixed in AliPMDrechit, Initializat...
[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
57acd2d2 46//____________________________________________________________________________
4e25ac79 47Double_t * AliHMPIDQAChecker::Check(AliQAv1::ALITASK_t /*index*/)
57acd2d2 48{
49 Double_t * rv = new Double_t[AliRecoParam::kNSpecies] ;
50 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++)
51 rv[specie] = 0.0 ;
52 return rv ;
53}
54
fdd7f404 55//_________________________________________________________________
4e25ac79 56Double_t * AliHMPIDQAChecker::Check(AliQAv1::ALITASK_t index, TObjArray ** list)
fdd7f404 57{
58//
59// Main check function: Depending on the TASK, different checks are applied
60// At the moment: check for empty histograms and checks for RecPoints
61
57acd2d2 62 Double_t * check = new Double_t[AliRecoParam::kNSpecies] ;
63
57acd2d2 64 AliInfo(Form("Fix needed ....."));
4e25ac79 65 char * detOCDBDir = Form("HMPID/%s/%s", AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName()) ;
fc7e0df2 66 AliCDBEntry *QARefRec = AliQAManager::QAManager()->Get(detOCDBDir);
fdd7f404 67 if( !QARefRec){
57acd2d2 68 AliInfo("QA reference data NOT retrieved for Reconstruction check. No HMPIDChecker ...exiting");
69 return check;
fdd7f404 70 }
71
72// checking for empy histograms
57acd2d2 73 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
74 check[specie] = 1.0;
4e25ac79 75 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
57acd2d2 76 continue ;
77 if(CheckEntries(list[specie]) == 0) {
78 AliWarning("histograms are empty");
79 check[specie] = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
80 }
81
82 // checking rec points
4e25ac79 83 if(index == AliQAv1::kREC) check[specie] = CheckRecPoints(list[specie],(TObjArray *)QARefRec->GetObject());
fdd7f404 84
57acd2d2 85 //default check response. It will be changed when reasonable checks will be considered
86 else check[specie] = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run
87 }
fdd7f404 88 return check;
89
90}
91//_________________________________________________________________
92Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
93{
94 //
95 // check on the QA histograms on the input list:
96 //
97
98 Double_t test = 0.0 ;
99 Int_t count = 0 ;
100
101 if (list->GetEntries() == 0){
102 test = 1. ; // nothing to check
103 }
104 else {
105 TIter next(list) ;
106 TH1 * hdata ;
107 count = 0 ;
108 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
109 if (hdata) {
110 Double_t rv = 0.;
111 //Printf("hitogram %s has entries: %f ",hdata->GetName(),hdata->GetEntries());
112 if(hdata->GetEntries()>0)rv=1;
113 count++ ;
114 test += rv ;
115 }
116 else{
117 AliError("Data type cannot be processed") ;
118 }
119
120 }
121 if (count != 0) {
122 if (test==0) {
123 AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
124 test = 0.; //upper limit value to set kWARNING flag for a task
125 }
126 else test = 1 ;
127 }
128 }
129
130 return test ;
131}
132//_________________________________________________________________
133
134Double_t AliHMPIDQAChecker::CheckRecPoints(TObjArray *listrec, TObjArray *listref) const
135{
0f46609e 136 //
fbe0560f 137 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
138 //
fdd7f404 139
140 Float_t checkresponse = 0;
141
142 Float_t counter = 0 ;
143 TIter next(listrec) ;
144 TH1* histo;
145 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
3762c8fb 146 //PH The histogram should have at least 3 bins with entries
147 Int_t nbinsabove = 0;
326d0a78 148 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);
3762c8fb 149 if (histo->GetBinContent(ibin)>0) nbinsabove++;
150 }
151
152 if( nbinsabove < 3 ) counter++;
fdd7f404 153 else {
154 TString h = histo->GetTitle();
155 if(h.Contains("Zoom")){
326d0a78 156 histo->Fit("expo","LQ0","",1,50);
fdd7f404 157 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
158 }
159 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
160 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
161 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
5379c4a3 162 AliDebug(AliQAv1::GetQADebugLevel(),Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
fdd7f404 163 }
164 }
165 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
166
167 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
168 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
169 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
170 return checkresponse;
171}
172//________________________________________________________________
e4e6ae61 173