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