]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliQualAssChecker.cxx
Updated QA classes (Yves)
[u/mrichter/AliRoot.git] / STEER / AliQualAssChecker.cxx
CommitLineData
421ab0fb 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/* $Id: */
17
18///////////////////////////////////////////////////////////////////////////////
19// //
20// class for running the Quality Assurance Checker
21// //
22///////////////////////////////////////////////////////////////////////////////
23
24#include "AliLog.h"
25#include "AliModule.h"
421ab0fb 26#include "AliQualAss.h"
27#include "AliQualAssChecker.h"
a5fa6165 28#include "AliQualAssCheckerBase.h"
421ab0fb 29
a5fa6165 30#include <TKey.h>
421ab0fb 31#include <TObjArray.h>
a5fa6165 32#include <TPluginManager.h>
33#include <TROOT.h>
421ab0fb 34#include <TStopwatch.h>
35#include <TString.h>
a5fa6165 36#include <TSystem.h>
421ab0fb 37
38ClassImp(AliQualAssChecker)
a5fa6165 39 TFile * AliQualAssChecker::fgQAResultFile = 0x0 ;
40 TString AliQualAssChecker::fgQAResultDirName = "local://RUN/";
41 TString AliQualAssChecker::fgQAResultFileName = "QA.root" ;
421ab0fb 42
43//_____________________________________________________________________________
44AliQualAssChecker::AliQualAssChecker(const char* name, const char* title) :
45 TNamed(name, title),
a5fa6165 46 fDataFile(0x0),
47 fRefDirName("/QA/Ref/"),
48 fRefName("QA.root"),
49 fFoundDetectors(".")
421ab0fb 50{
a5fa6165 51 // ctor: initialise checkers and open the data file
52 for (Int_t det = 0 ; det < AliQualAss::kNDET ; det++)
53 fCheckers[det] = NULL ;
54
55 GetDataFile() ;
421ab0fb 56}
57
58//_____________________________________________________________________________
59AliQualAssChecker::AliQualAssChecker(const AliQualAssChecker& qac) :
60 TNamed(qac),
a5fa6165 61 fDataFile(qac.fDataFile),
62 fRefDirName(qac.fRefDirName),
63 fRefName(qac.fRefName),
64 fFoundDetectors(qac.fFoundDetectors)
421ab0fb 65{
a5fa6165 66 // copy constructor
67
68 for (Int_t det = 0 ; det < AliQualAss::kNDET ; det++)
69 fCheckers[det] = NULL ;
421ab0fb 70}
71
72//_____________________________________________________________________________
73AliQualAssChecker& AliQualAssChecker::operator = (const AliQualAssChecker& qac)
74{
75// assignment operator
76
77 this->~AliQualAssChecker();
78 new(this) AliQualAssChecker(qac);
79 return *this;
80}
81
82//_____________________________________________________________________________
83AliQualAssChecker::~AliQualAssChecker()
84{
85// clean up
421ab0fb 86}
87
88//_____________________________________________________________________________
a5fa6165 89TFile * AliQualAssChecker:: GetDataFile()
90{
91 // Open if necessary the Data file and return its pointer
92
93 if (!fDataFile)
94 fDataFile = TFile::Open(AliQualAss::GetDataName()) ;
95 if (!fDataFile)
96 AliFatal(Form("QA Data File %s does not exist", AliQualAss::GetDataName() )) ;
97 return fDataFile ;
98}
99
100//_____________________________________________________________________________
101TFile * AliQualAssChecker:: GetQAResultFile()
421ab0fb 102{
103 // Check if file to store QA exists, if not create it
104
a5fa6165 105 if (fgQAResultFile) {
106 if (fgQAResultFile->IsOpen()){
107 fgQAResultFile->Close() ;
108 fgQAResultFile = 0x0 ;
421ab0fb 109 }
110 }
a5fa6165 111 if ( fgQAResultFileName.Contains("local://"))
112 fgQAResultFileName.ReplaceAll("local:/", "") ;
113
421ab0fb 114 TString opt("") ;
a5fa6165 115 if ( !gSystem->AccessPathName(fgQAResultFileName) )
421ab0fb 116 opt = "UPDATE" ;
117 else
118 opt = "NEW" ;
a5fa6165 119 fgQAResultFile = TFile::Open(fgQAResultFileName, opt) ;
120
121 return fgQAResultFile ;
421ab0fb 122}
123
124//_____________________________________________________________________________
a5fa6165 125 AliQualAssCheckerBase * AliQualAssChecker::GetDetQualAssChecker(Int_t det)
421ab0fb 126{
a5fa6165 127 // Gets the Quality Assurance checker for the detector specified by its name
128
129 if (fCheckers[det])
130 return fCheckers[det];
131
132 TString detName(AliQualAss::GetDetName(det)) ;
133
134 AliInfo(Form("Retrieving QA checker for %s", detName.Data())) ;
135 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
136 TString qacName = "Ali" + detName + "QualAssChecker" ;
137
138 AliQualAssCheckerBase * qac = NULL ;
139 // first check if a plugin is defined for the quality assurance checker
140 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQualAssChecker", detName.Data());
141 // if not, add a plugin for it
142 if (!pluginHandler) {
143 //AliInfo(Form("defining plugin for %s", qacName.Data()));
144 TString libs = gSystem->GetLibraries();
145
146 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0))
147 pluginManager->AddHandler("AliQualAssChecker", detName, qacName, detName + "qac", qacName + "()");
148 else
149 pluginManager->AddHandler("AliQualAssChecker", detName, qacName, detName, qacName + "()");
150
151 pluginHandler = pluginManager->FindHandler("AliQualAssChecker", detName);
152
153 if (pluginHandler && (pluginHandler->LoadPlugin() == 0))
154 qac = (AliQualAssCheckerBase *) pluginHandler->ExecPlugin(0);
155
156 if (qac)
157 fCheckers[det] = qac ;
421ab0fb 158 }
a5fa6165 159
160 return qac ;
161}
162
163
164//_____________________________________________________________________________
165TDirectory * AliQualAssChecker::GetRefSubDir(const char * det, const char * task)
166{
167 // Opens and returns the file with the reference data
168 TFile * f = TFile::Open(fRefDirName, "READ") ;
169 if (!f)
170 AliFatal(Form("Cannot find reference file %s", fRefDirName.Data())) ;
171 TDirectory * rv = NULL ;
172 rv = f->GetDirectory(det) ;
173 if (!rv) {
174 AliWarning(Form("Directory %s not found in %d", det, fRefDirName.Data())) ;
175 } else {
176 rv = rv->GetDirectory(task) ;
177 if (!rv)
178 AliWarning(Form("Directory %s/%s not found in %s", det, task, fRefDirName.Data())) ;
179 }
180 return rv ;
421ab0fb 181}
182
183//_____________________________________________________________________________
184Bool_t AliQualAssChecker::Run()
185{
a5fa6165 186 // run the Quality Assurance Checker for all tasks Hits, SDigits, Digits, recpoints, tracksegments, recparticles and ESDs
421ab0fb 187
188 Bool_t rv = kFALSE ;
a5fa6165 189
421ab0fb 190 TStopwatch stopwatch;
191 stopwatch.Start();
192
a5fa6165 193 //search for all detectors QA directories
194 TList * detKeyList = GetDataFile()->GetListOfKeys() ;
195 TIter nextd(detKeyList) ;
196 TKey * detKey ;
197 while ( (detKey = dynamic_cast<TKey *>(nextd()) ) ) {
198 AliInfo(Form("Found %s", detKey->GetName())) ;
199 //Check which detector
200 TString detName ;
201 TString detNameQA(detKey->GetName()) ;
202 Int_t det ;
203 for ( det = 0; det < AliQualAss::kNDET ; det++) {
204 detName = AliQualAss::GetDetName(det) ;
205 if (detNameQA.Contains(detName)) {
206 fFoundDetectors+=detName ;
207 fFoundDetectors+="." ;
208 break ;
209 }
210 }
211 TDirectory * detDir = GetDataFile()->GetDirectory(detKey->GetName()) ;
212 TList * taskKeyList = detDir->GetListOfKeys() ;
213 TIter nextt(taskKeyList) ;
214 TKey * taskKey ;
215 // now search for the tasks dir
216 while ( (taskKey = static_cast<TKey *>(nextt()) ) ) {
217 TString taskName( taskKey->GetName() ) ;
218 AliInfo(Form("Found %s", taskName.Data())) ;
219 TDirectory * taskDir = detDir->GetDirectory(taskName.Data()) ;
220 taskDir->cd() ;
221 AliQualAssCheckerBase * qac = GetDetQualAssChecker(det) ;
222 if (qac)
223 AliInfo(Form("QA checker found for %s", detName.Data())) ;
224 if (!qac)
225 AliFatal(Form("QA checker not found for %s", detName.Data())) ;
226 AliQualAss::ALITASK index = AliQualAss::kNULLTASK ;
227 if ( taskName == AliQualAss::GetTaskName(AliQualAss::kHITS) )
228 index = AliQualAss::kSIM ;
229 if ( taskName == AliQualAss::GetTaskName(AliQualAss::kSDIGITS) )
230 index = AliQualAss::kSIM ;
231 if ( taskName == AliQualAss::GetTaskName(AliQualAss::kDIGITS) )
232 index = AliQualAss::kSIM ;
233 if ( taskName == AliQualAss::GetTaskName(AliQualAss::kRECPOINTS) )
234 index = AliQualAss::kREC ;
235 if ( taskName == AliQualAss::GetTaskName(AliQualAss::kTRACKSEGMENTS) )
236 index = AliQualAss::kREC ;
237 if ( taskName == AliQualAss::GetTaskName(AliQualAss::kRECPARTICLES) )
238 index = AliQualAss::kREC ;
239 if ( taskName == AliQualAss::GetTaskName(AliQualAss::kESDS) )
240 index = AliQualAss::kESD ;
241 qac->Init(AliQualAss::DETECTORINDEX(det)) ;
242 qac->SetRefandData(GetRefSubDir(detNameQA.Data(), taskName.Data()), taskDir) ;
243 qac->Run(index) ;
421ab0fb 244 }
a5fa6165 245 }
246 AliInfo("QA performed for following detectors:") ;
247 for ( Int_t det = 0; det < AliQualAss::kNDET; det++) {
248 if (fFoundDetectors.Contains(AliQualAss::GetDetName(det))) {
249 printf("%s, ",AliQualAss::GetDetName(det)) ;
250 fFoundDetectors.ReplaceAll(AliQualAss::GetDetName(det), "") ;
251 }
421ab0fb 252 }
a5fa6165 253 printf("\n") ;
254 rv = kTRUE ;
421ab0fb 255
a5fa6165 256 return rv ;
257
421ab0fb 258}
259
260//_____________________________________________________________________________
a5fa6165 261void AliQualAssChecker::SetQAResultDirName(const char * name)
421ab0fb 262{
263 // Set the root directory where to store the QA status object
264
a5fa6165 265 fgQAResultDirName.Prepend(name) ;
266 AliInfo(Form("QA results are in %s", fgQAResultDirName.Data())) ;
267 if ( fgQAResultDirName.Contains("local://"))
268 fgQAResultDirName.ReplaceAll("local:/", "") ;
269 fgQAResultFileName.Prepend(fgQAResultDirName) ;
421ab0fb 270}
271
272//_____________________________________________________________________________
a5fa6165 273void AliQualAssChecker::SetRefDirName(const char * name)
421ab0fb 274{
275 // Set the root directory of reference data
276
a5fa6165 277 fRefDirName.Prepend(name) ;
278 fRefDirName.Append(fRefName) ;
279 AliInfo(Form("Reference data are taken from %s", fRefDirName.Data())) ;
280 if ( fRefDirName.Contains("local://"))
281 fRefDirName.ReplaceAll("local:/", "") ;
421ab0fb 282}
283
284
285
286
287