]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDQAChecker.cxx
AliEveTrackCounter, AliEveTrackCounterEditor
[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"
b990acb3 38#include "AliQA.h"
39#include "AliQAChecker.h"
40#include "AliHMPIDQAChecker.h"
fdd7f404 41#include "AliCDBEntry.h"
42#include "AliCDBManager.h"
e4e6ae61 43
b990acb3 44ClassImp(AliHMPIDQAChecker)
e4e6ae61 45
fdd7f404 46//_________________________________________________________________
47const Double_t AliHMPIDQAChecker::Check(AliQA::ALITASK_t index, TObjArray * list)
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
53 AliDebug(1,Form("AliHMPIDChecker"));
54 AliCDBEntry *QARefRec = AliCDBManager::Instance()->Get("HMPID/QARef/Rec");
55 if( !QARefRec){
56 AliInfo("QA reference data NOT retrieved for Recostruction check. No HMPIDChecker ...exiting");
57 return 1.;
58 }
59
60// checking for empy histograms
61 Double_t check =0;
62 if(CheckEntries(list) == 0) {
63 AliWarning("histograms are empty");
64 check = 0.4;//-> Corresponds to kWARNING see AliQACheckerBase::Run
65 return check;
66 }
67
68// checking rec points
69 if(index == AliQA::kREC) check = CheckRecPoints(list,(TObjArray *)QARefRec->GetObject());
70
71//default check response. It will be changed when reasonable checks will be considered
72 else check = 0.7 ; // /-> Corresponds to kINFO see AliQACheckerBase::Run
73
74 return check;
75
76}
77//_________________________________________________________________
78Double_t AliHMPIDQAChecker::CheckEntries(TObjArray * list) const
79{
80 //
81 // check on the QA histograms on the input list:
82 //
83
84 Double_t test = 0.0 ;
85 Int_t count = 0 ;
86
87 if (list->GetEntries() == 0){
88 test = 1. ; // nothing to check
89 }
90 else {
91 TIter next(list) ;
92 TH1 * hdata ;
93 count = 0 ;
94 while ( (hdata = dynamic_cast<TH1 *>(next())) ) {
95 if (hdata) {
96 Double_t rv = 0.;
97 //Printf("hitogram %s has entries: %f ",hdata->GetName(),hdata->GetEntries());
98 if(hdata->GetEntries()>0)rv=1;
99 count++ ;
100 test += rv ;
101 }
102 else{
103 AliError("Data type cannot be processed") ;
104 }
105
106 }
107 if (count != 0) {
108 if (test==0) {
109 AliWarning("Histograms are booked for THIS specific Task, but they are all empty: setting flag to kWARNING");
110 test = 0.; //upper limit value to set kWARNING flag for a task
111 }
112 else test = 1 ;
113 }
114 }
115
116 return test ;
117}
118//_________________________________________________________________
119
120Double_t AliHMPIDQAChecker::CheckRecPoints(TObjArray *listrec, TObjArray *listref) const
121{
0f46609e 122 //
123 // check on the HMPID RecPoints by using expo fit and Kolmogorov Test:
124 //
fdd7f404 125
126 Float_t checkresponse = 0;
127
128 Float_t counter = 0 ;
129 TIter next(listrec) ;
130 TH1* histo;
131 while ( (histo = dynamic_cast<TH1 *>(next())) ) {
0f46609e 132 if( histo->GetEntries() < 3 ) counter++;
fdd7f404 133 else {
134 TString h = histo->GetTitle();
135 if(h.Contains("Zoom")){
136 histo->Fit("expo","Q0","",1,50);
137 if(histo->GetFunction("expo")->GetParameter(1) !=0 ) if(TMath::Abs((-1./(histo->GetFunction("expo"))->GetParameter(1)) - 35 ) > 5) counter++;
138 }
139 if(h.Contains("size MIP")) if(TMath::Abs(histo->GetMean()-5) > 2) counter++;
140 if(h.Contains("size Phots")) if(TMath::Abs(histo->GetMean()-2) > 2) counter++;
141 if(h.Contains("distribution")) if(histo->KolmogorovTest((TH1F *)listref->At(0))<0.8) counter++;
142 AliDebug(1,Form(" Kolm. test : %f ",histo->KolmogorovTest((TH1F *)listref->At(0))));
143 }
144 }
145 Float_t response = counter/(7.+7.+42.+42.); // 7.+7.+42 +42 = N checked histograms (-> To be replaced by listrec->GetEntries())
146
147 if(response < 0.1) checkresponse = 0.7; // <10% of the check histograms show a failing check -> Corresponds to kINFO see AliQACheckerBase::Run
148 else if(response < 0.5) checkresponse = 0.4; // 50% of the check histograms show a failing check -> Corresponds to kWARNING see AliQACheckerBase::Run
149 else checkresponse = 0.001; // > 50% of the check histograms show a failing check -> Corresponds to kERROR see AliQACheckerBase::Run
150 return checkresponse;
151}
152//________________________________________________________________
e4e6ae61 153