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