]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQAChecker.cxx
The present commit corresponds to an important change in the way the
[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 "AliCDBEntry.h"
25 #include "AliCDBManager.h"
26 #include "AliCDBStorage.h"
27 #include "AliRunInfo.h" 
28 #include "AliLog.h"
29 #include "AliModule.h" 
30 #include "AliQA.h"
31 #include "AliQAChecker.h"
32 #include "AliQACheckerBase.h"
33 #include "AliCorrQAChecker.h"
34 #include "AliGlobalQAChecker.h"
35 #include "AliGRPObject.h"
36
37 #include <TKey.h>
38 #include <TObjArray.h>
39 #include <TObjString.h>
40 #include <TPluginManager.h> 
41 #include <TROOT.h>
42 #include <TStopwatch.h> 
43 #include <TString.h> 
44 #include <TSystem.h> 
45 #include <TList.h>
46 #include <TNtupleD.h>
47
48 ClassImp(AliQAChecker)
49   AliQAChecker * AliQAChecker::fgQAChecker = 0x0 ;
50
51 //_____________________________________________________________________________
52 AliQAChecker::AliQAChecker(const char* name, const char* title) :
53   TNamed(name, title),
54   fDataFile(0x0), 
55   fRunInfo(0x0), 
56   fRunInfoOwner(kFALSE), 
57   fRefFile(0x0), 
58   fFoundDetectors("."), 
59   fEventSpecie(AliRecoParam::kDefault) 
60 {
61   // ctor: initialise checkers and open the data file   
62   for (Int_t det = 0 ; det < AliQA::kNDET ; det++) 
63     fCheckers[det] = NULL ; 
64 }
65
66 //_____________________________________________________________________________
67 AliQAChecker::AliQAChecker(const AliQAChecker& qac) :
68   TNamed(qac),
69   fDataFile(qac.fDataFile), 
70   fRunInfo(qac.fRunInfo), 
71   fRunInfoOwner(kFALSE),   
72   fRefFile(qac.fRefFile), 
73   fFoundDetectors(qac.fFoundDetectors),
74   fEventSpecie(qac.fEventSpecie)
75 {
76   // copy constructor
77   
78   for (Int_t det = 0 ; det < AliQA::kNDET ; det++) 
79     fCheckers[det] = NULL ; 
80 }
81
82 //_____________________________________________________________________________
83 AliQAChecker& AliQAChecker::operator = (const AliQAChecker& qac)
84 {
85 // assignment operator
86
87   this->~AliQAChecker();
88   new(this) AliQAChecker(qac);
89   return *this;
90 }
91
92 //_____________________________________________________________________________
93 AliQAChecker::~AliQAChecker()
94 {
95 // clean up
96   if (fRunInfo)
97     delete fRunInfo ; 
98   delete [] fCheckers ; 
99   AliQA::Close() ; 
100 }
101
102 //_____________________________________________________________________________
103   AliQACheckerBase * AliQAChecker::GetDetQAChecker(Int_t det)
104 {
105         // Gets the Quality Assurance checker for the detector specified by its name
106         
107         if (fCheckers[det]) 
108     return fCheckers[det];
109
110         AliQACheckerBase * qac = NULL ;
111
112         TString detName(AliQA::GetDetName(det)) ; 
113         
114         if (det == AliQA::kGLOBAL) {
115                 qac = new AliGlobalQAChecker() ; 
116         } else if (det == AliQA::kCORR) {
117                 qac = new AliCorrQAChecker() ; 
118         } else {
119                 AliDebug(1, Form("Retrieving QA checker for %s", detName.Data())) ; 
120                 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
121                 TString qacName = "Ali" + detName + "QAChecker" ;
122
123                 // first check if a plugin is defined for the quality assurance checker
124                 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQAChecker", detName.Data());
125                 // if not, add a plugin for it
126                 if (!pluginHandler) {
127                         //AliInfo(Form("defining plugin for %s", qacName.Data()));
128                         TString libs = gSystem->GetLibraries();
129                 
130                         if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0))
131                                 pluginManager->AddHandler("AliQAChecker", detName, qacName, detName + "qac", qacName + "()");
132                         else 
133                                 pluginManager->AddHandler("AliQAChecker", detName, qacName, detName, qacName + "()");
134
135                         pluginHandler = pluginManager->FindHandler("AliQAChecker", detName);    
136
137                         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) 
138                                 qac = (AliQACheckerBase *) pluginHandler->ExecPlugin(0);
139   
140                 }
141         }
142         if (qac) 
143                 fCheckers[det] = qac ;
144         
145         return qac ; 
146 }
147  
148 //_____________________________________________________________________________
149 void AliQAChecker::GetRefSubDir(const char * det, const char * task, TDirectory *& dirFile, TObjArray **& dirOCDB)     
150
151   // Opens and returns the file with the reference data 
152   dirFile = NULL ; 
153   TString refStorage(AliQA::GetQARefStorage()) ;
154 //  if (refStorage.Contains(AliQA::GetLabLocalFile())) {        
155 //    refStorage.ReplaceAll(AliQA::GetLabLocalFile(), "") ; 
156 //    refStorage += AliQA::GetQARefFileName() ;
157 //    if ( fRefFile ) 
158 //      if ( fRefFile->IsOpen() ) 
159 //                                      fRefFile->Close() ; 
160 //    fRefFile = TFile::Open(refStorage.Data()) ; 
161 //    if (!fRefFile) { 
162 //      AliError(Form("Cannot find reference file %s", refStorage.Data())) ; 
163 //      dirFile = NULL ; 
164 //    }
165 //    dirFile = fRefFile->GetDirectory(det) ; 
166 //    if (!dirFile) {
167 //      AliWarning(Form("Directory %s not found in %d", det, refStorage.Data())) ; 
168 //    } else {
169 //                      dirFile = dirFile->GetDirectory(task) ; 
170 //      if (!dirFile) 
171 //                              AliWarning(Form("Directory %s/%s not found in %s", det, task, refStorage.Data())) ; 
172 //    }  
173 //  } else 
174   if (!refStorage.Contains(AliQA::GetLabLocalOCDB()) && !refStorage.Contains(AliQA::GetLabAliEnOCDB())) {
175     AliError(Form("%s is not a valid location for reference data", refStorage.Data())) ; 
176     return ; 
177   } else {
178     AliCDBManager* man = AliCDBManager::Instance() ;
179     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
180       if ( !AliQA::Instance()->IsEventSpecieSet(specie) ) 
181         continue ; 
182       //if ( strcmp(AliQA::GetRefDataDirName(), "") == 0 ) { // the name of the last level of the directory is not set (EventSpecie)
183         // Get it from RunInfo
184         //if (!fRunInfo)  // not yet set, get the info from GRP
185         //  LoadRunInfoFromGRP() ; 
186       AliQA::SetQARefDataDirName(specie) ;
187       //}
188       if ( ! man->GetLock() ) { 
189         man->SetDefaultStorage(AliQA::GetQARefStorage()) ; 
190         man->SetSpecificStorage("*", AliQA::GetQARefStorage()) ;
191       }
192       char * detOCDBDir = Form("%s/%s/%s", det, AliQA::GetRefOCDBDirName(), AliQA::GetRefDataDirName()) ; 
193       AliCDBEntry * entry = man->Get(detOCDBDir, man->GetRun()) ;
194       if (entry) {
195         dirOCDB = new TObjArray*[AliRecoParam::kNSpecies] ;     
196         TList * listDetQAD = dynamic_cast<TList *>(entry->GetObject()) ;
197         TIter next(listDetQAD) ;
198         TObjArray * ar ; 
199         while ( (ar = (TObjArray*)next()) )
200           if ( listDetQAD ) 
201           dirOCDB[specie] = dynamic_cast<TObjArray *>(listDetQAD->FindObject(Form("%s/%s", task, AliRecoParam::GetEventSpecieName(specie)))) ; 
202       }
203     }
204   }
205 }
206
207 //_____________________________________________________________________________
208 AliQAChecker * AliQAChecker::Instance()
209 {
210         // returns unique instance of the checker
211   if ( ! fgQAChecker ) 
212    fgQAChecker = new AliQAChecker() ; 
213   return fgQAChecker ;  
214 }
215
216 //_____________________________________________________________________________
217 void AliQAChecker::LoadRunInfoFromGRP() 
218 {
219   AliCDBManager* man = AliCDBManager::Instance() ;
220   AliCDBEntry* entry = man->Get(AliQA::GetGRPPath().Data());
221   AliGRPObject* grpObject = 0x0;
222   if (entry) {
223
224           TMap* m = dynamic_cast<TMap*>(entry->GetObject());  // old GRP entry
225
226           if (m) {
227             AliInfo("It is a map");
228             //m->Print();
229             grpObject = new AliGRPObject();
230                  grpObject->ReadValuesFromMap(m);
231     }
232
233     else {
234             AliInfo("It is a new GRP object");
235         grpObject = dynamic_cast<AliGRPObject*>(entry->GetObject());  // new GRP entry
236     }
237
238     entry->SetOwner(0);
239     AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
240   }
241
242   if (!grpObject) {
243      AliFatal("No GRP entry found in OCDB!");
244   }
245
246   TString lhcState = grpObject->GetLHCState();
247   if (lhcState==AliGRPObject::GetInvalidString()) {
248     AliError("GRP/GRP/Data entry:  missing value for the LHC state ! Using UNKNOWN");
249     lhcState = "UNKNOWN";
250   }
251
252   TString beamType = grpObject->GetBeamType();
253   if (beamType==AliGRPObject::GetInvalidString()) {
254     AliError("GRP/GRP/Data entry:  missing value for the beam type ! Using UNKNOWN");
255     beamType = "UNKNOWN";
256   }
257
258   Float_t beamEnergy = grpObject->GetBeamEnergy();
259   if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
260     AliError("GRP/GRP/Data entry:  missing value for the beam energy ! Using 0");
261     beamEnergy = 0;
262   }
263
264   TString runType = grpObject->GetRunType();
265   if (runType==AliGRPObject::GetInvalidString()) {
266     AliError("GRP/GRP/Data entry:  missing value for the run type ! Using UNKNOWN");
267     runType = "UNKNOWN";
268   }
269
270   Int_t activeDetectors = grpObject->GetDetectorMask();
271   if (activeDetectors==AliGRPObject::GetInvalidInt()) {
272     AliError("GRP/GRP/Data entry:  missing value for the detector mask ! Using 1074790399");
273     activeDetectors = 1074790399;
274   }
275
276   fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
277
278   fRunInfoOwner = kTRUE ; 
279
280   // set the event specie
281   fEventSpecie = AliRecoParam::kDefault ;
282   if (strcmp(runType,"PHYSICS")) {
283     // Not a physics run, the event specie is set to kCalib
284     fEventSpecie = AliRecoParam::kCalib ;
285     return;
286   }
287   if (strcmp(lhcState,"STABLE_BEAMS") == 0) {
288     // Heavy ion run (any beam tha is not pp, the event specie is set to kHighMult
289     fEventSpecie = AliRecoParam::kHighMult ;
290     if ((strcmp(beamType,"p-p") == 0) ||
291         (strcmp(beamType,"p-")  == 0) ||
292         (strcmp(beamType,"-p")  == 0) ||
293         (strcmp(beamType,"P-P") == 0) ||
294         (strcmp(beamType,"P-")  == 0) ||
295         (strcmp(beamType,"-P")  == 0)) {
296       // Proton run, the event specie is set to kLowMult
297       fEventSpecie = AliRecoParam::kLowMult ;
298     }
299     else if (strcmp(beamType,"-") == 0) {
300       // No beams, we assume cosmic data
301       fEventSpecie = AliRecoParam::kCosmic ;
302     }
303     else if (strcmp(beamType,"UNKNOWN") == 0) {
304       // No LHC beam information is available, we use the default event specie
305       fEventSpecie = AliRecoParam::kDefault ;
306     }
307   }
308 }
309
310 //_____________________________________________________________________________
311 Bool_t AliQAChecker::Run(const char * fileName)
312 {
313   // run the Quality Assurance Checker for all tasks Hits, SDigits, Digits, recpoints, tracksegments, recparticles and ESDs
314   // starting from data in file  
315   TStopwatch stopwatch;
316   stopwatch.Start();
317
318   //search for all detectors QA directories
319   TList * detKeyList = AliQA::GetQADataFile(fileName)->GetListOfKeys() ; 
320   TIter nextd(detKeyList) ; 
321   TKey * detKey ; 
322   while ( (detKey = dynamic_cast<TKey *>(nextd()) ) ) {
323     AliDebug(1, Form("Found %s", detKey->GetName())) ;
324     //Check which detector
325     TString detName ; 
326     TString detNameQA(detKey->GetName()) ; 
327     Int_t det ; 
328     for ( det = 0; det < AliQA::kNDET ; det++) {
329       detName = AliQA::GetDetName(det) ; 
330       if (detNameQA.Contains(detName)) {
331         fFoundDetectors+=detName ; 
332         fFoundDetectors+="." ;          
333         break ; 
334       }
335     } 
336     TDirectory * detDir = AliQA::GetQADataFile(fileName)->GetDirectory(detKey->GetName()) ; 
337     TList * taskKeyList = detDir->GetListOfKeys() ;
338     TIter nextt(taskKeyList) ; 
339     TKey * taskKey ; 
340     // now search for the tasks dir
341     while ( (taskKey = static_cast<TKey *>(nextt()) ) ) {
342       TString taskName( taskKey->GetName() ) ; 
343       AliInfo(Form("Found %s", taskName.Data())) ;
344       TDirectory * taskDir = detDir->GetDirectory(taskName.Data()) ; 
345       taskDir->cd() ; 
346       AliQACheckerBase * qac = GetDetQAChecker(det) ; 
347       if (qac)
348         AliInfo(Form("QA checker found for %s", detName.Data())) ; 
349       if (!qac)
350         AliFatal(Form("QA checker not found for %s", detName.Data())) ; 
351       AliQA::ALITASK_t index = AliQA::kNULLTASK ; 
352       if ( taskName == AliQA::GetTaskName(AliQA::kHITS) ) 
353         index = AliQA::kSIM ; 
354       if ( taskName == AliQA::GetTaskName(AliQA::kSDIGITS) ) 
355         index = AliQA::kSIM ; 
356       if ( taskName == AliQA::GetTaskName(AliQA::kDIGITS) ) 
357         index = AliQA::kSIM ; 
358       if ( taskName == AliQA::GetTaskName(AliQA::kRECPOINTS) ) 
359         index = AliQA::kREC ; 
360       if ( taskName == AliQA::GetTaskName(AliQA::kTRACKSEGMENTS) ) 
361         index = AliQA::kREC ; 
362       if ( taskName == AliQA::GetTaskName(AliQA::kRECPARTICLES) ) 
363         index = AliQA::kREC ; 
364       if ( taskName == AliQA::GetTaskName(AliQA::kESDS) ) 
365         index = AliQA::kESD ; 
366       qac->Init(AliQA::DETECTORINDEX_t(det)) ; 
367       
368       TDirectory * refDir     = NULL ; 
369       TObjArray ** refOCDBDir = NULL ;  
370       GetRefSubDir(detNameQA.Data(), taskName.Data(), refDir, refOCDBDir) ;
371                   qac->SetRefandData(refDir, refOCDBDir, taskDir) ;
372                   qac->Run(index) ; 
373     }
374   }
375   AliInfo("QA performed for following detectors:") ; 
376   for ( Int_t det = 0; det < AliQA::kNDET; det++) {
377     if (fFoundDetectors.Contains(AliQA::GetDetName(det))) {
378       printf("%s, ",AliQA::GetDetName(det)) ; 
379       fFoundDetectors.ReplaceAll(AliQA::GetDetName(det), "") ; 
380     }   
381   }
382   printf("\n") ; 
383   return kTRUE ; 
384 }
385
386 //_____________________________________________________________________________
387 Bool_t AliQAChecker::Run(AliQA::DETECTORINDEX_t det, AliQA::TASKINDEX_t task, TObjArray ** list)
388 {
389         // run the Quality Assurance Checker for detector det, for task task starting from data in list
390
391         AliQACheckerBase * qac = GetDetQAChecker(det) ; 
392         if (qac)
393                 AliDebug(1, Form("QA checker found for %s", AliQA::GetDetName(det).Data())) ;
394         if (!qac)
395                 AliError(Form("QA checker not found for %s", AliQA::GetDetName(det).Data())) ; 
396   
397         AliQA::ALITASK_t index = AliQA::kNULLTASK ; 
398         if ( task == AliQA::kRAWS ) 
399                 index = AliQA::kRAW ; 
400         else if ( task == AliQA::kHITS ) 
401                 index = AliQA::kSIM ; 
402         else if ( task == AliQA::kSDIGITS ) 
403                 index = AliQA::kSIM ; 
404         else if ( task == AliQA::kDIGITS ) 
405                 index = AliQA::kSIM ; 
406         else if ( task == AliQA::kRECPOINTS ) 
407                 index = AliQA::kREC ; 
408         else if ( task == AliQA::kTRACKSEGMENTS ) 
409                 index = AliQA::kREC ; 
410         else if ( task == AliQA::kRECPARTICLES ) 
411                 index = AliQA::kREC ; 
412         else if ( task == AliQA::kESDS ) 
413                 index = AliQA::kESD ; 
414
415         TDirectory * refDir     = NULL ; 
416         TObjArray ** refOCDBDir = NULL  ;       
417   qac->Init(det) ; 
418   GetRefSubDir(AliQA::GetDetName(det), AliQA::GetTaskName(task), refDir, refOCDBDir) ;
419   qac->SetRefandData(refDir, refOCDBDir) ; 
420   qac->Run(index, list) ; 
421         return kTRUE ; 
422 }
423
424 //_____________________________________________________________________________
425 Bool_t AliQAChecker::Run(AliQA::DETECTORINDEX_t det, AliQA::TASKINDEX_t task, TNtupleD ** list)
426 {
427         // run the Quality Assurance Checker for detector det, for task task starting from data in list
428   
429         AliQACheckerBase * qac = GetDetQAChecker(det) ; 
430         if (qac)
431                 AliDebug(1, Form("QA checker found for %s", AliQA::GetDetName(det).Data())) ;
432         if (!qac)
433                 AliError(Form("QA checker not found for %s", AliQA::GetDetName(det).Data())) ; 
434   
435         AliQA::ALITASK_t index = AliQA::kNULLTASK ; 
436         if ( task == AliQA::kRAWS ) 
437                 index = AliQA::kRAW ; 
438         else if ( task == AliQA::kHITS ) 
439                 index = AliQA::kSIM ; 
440         else if ( task == AliQA::kSDIGITS ) 
441                 index = AliQA::kSIM ; 
442         else if ( task == AliQA::kDIGITS ) 
443                 index = AliQA::kSIM ; 
444         else if ( task == AliQA::kRECPOINTS ) 
445                 index = AliQA::kREC ; 
446         else if ( task == AliQA::kTRACKSEGMENTS ) 
447                 index = AliQA::kREC ; 
448         else if ( task == AliQA::kRECPARTICLES ) 
449                 index = AliQA::kREC ; 
450         else if ( task == AliQA::kESDS ) 
451                 index = AliQA::kESD ; 
452   
453         TDirectory * refDir     = NULL ; 
454         TObjArray ** refOCDBDir = NULL ;        
455   qac->Init(det) ; 
456   GetRefSubDir(AliQA::GetDetName(det), AliQA::GetTaskName(task), refDir, refOCDBDir) ;
457   qac->SetRefandData(refDir, refOCDBDir) ; 
458   qac->Run(index, list) ; 
459   return kTRUE ; 
460 }