]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMakerSteer.cxx
Possibility to run the checker starting from a file
[u/mrichter/AliRoot.git] / STEER / AliQADataMakerSteer.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 #include <TFile.h>
19 #include <TPluginManager.h>
20 #include <TROOT.h>
21 #include <TString.h>
22 #include <TSystem.h>
23
24 #include "AliESDEvent.h"
25 #include "AliLog.h"
26 #include "AliModule.h"
27 #include "AliQA.h"
28 #include "AliQADataMaker.h"
29 #include "AliQADataMakerSteer.h" 
30 #include "AliRawReaderRoot.h"
31 #include "AliRun.h"
32 #include "AliRunLoader.h"
33
34 ClassImp(AliQADataMakerSteer) 
35
36 //_____________________________________________________________________________
37 AliQADataMakerSteer::AliQADataMakerSteer(const char* gAliceFilename, const char * name, const char * title) :
38         TNamed(name, title), 
39         fCycleSame(kFALSE),
40         fESD(NULL), 
41         fESDTree(NULL),
42         fFirst(kTRUE),  
43         fGAliceFileName(gAliceFilename), 
44         fRunNumber(0), 
45         fNumberOfEvents(0), 
46         fRawReader(NULL), 
47         fRunLoader(NULL)  
48 {
49         // default ctor
50         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
51                 fLoader[iDet]      = NULL ;
52                 fQADataMaker[iDet] = NULL ;
53                 fQACycles[iDet]    = 999999 ;   
54         }       
55 }
56
57 //_____________________________________________________________________________
58 AliQADataMakerSteer::AliQADataMakerSteer(const AliQADataMakerSteer & qas) : 
59         TNamed(qas), 
60         fCycleSame(kFALSE),
61         fESD(NULL), 
62         fESDTree(NULL), 
63         fFirst(qas.fFirst),  
64         fGAliceFileName(qas.fGAliceFileName), 
65         fRunNumber(qas.fRunNumber), 
66         fNumberOfEvents(qas.fNumberOfEvents), 
67         fRawReader(NULL), 
68         fRunLoader(NULL)  
69 {
70         // cpy ctor
71         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
72                 fLoader[iDet]      = qas.fLoader[iDet] ;
73                 fQADataMaker[iDet] = qas.fQADataMaker[iDet] ;
74                 fQACycles[iDet]    = qas.fQACycles[iDet] ;      
75         }
76 }
77
78 //_____________________________________________________________________________
79 AliQADataMakerSteer & AliQADataMakerSteer::operator = (const AliQADataMakerSteer & qas) 
80 {
81         // assignment operator
82   this->~AliQADataMakerSteer() ;
83   new(this) AliQADataMakerSteer(qas) ;
84   return *this ;
85 }
86
87 //_____________________________________________________________________________
88 AliQADataMakerSteer::~AliQADataMakerSteer() 
89 {
90         // dtor
91   for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
92     fLoader[iDet] = NULL;
93         if (fQADataMaker[iDet]) {
94                 (fQADataMaker[iDet])->Finish() ; 
95                 delete fQADataMaker[iDet] ;
96                 fQADataMaker[iDet] = NULL ;
97         }
98   }
99
100   fRunLoader = NULL ;
101   delete fRawReader ;
102   fRawReader = NULL ;
103 }
104
105 //_____________________________________________________________________________
106 AliLoader * AliQADataMakerSteer::GetLoader(Int_t iDet)
107 {
108         // get the loader for a detector
109         
110         TString detName = AliQA::GetDetName(iDet) ;
111     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
112         if (fLoader[iDet]) 
113                 return fLoader[iDet] ;
114         
115         // load the QA data maker object
116         TPluginManager* pluginManager = gROOT->GetPluginManager() ;
117         TString loaderName = "Ali" + detName + "Loader" ;
118
119         AliLoader * loader = NULL ;
120         // first check if a plugin is defined for the quality assurance data maker
121         TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
122         // if not, add a plugin for it
123         if (!pluginHandler) {
124                 AliDebug(1, Form("defining plugin for %s", loaderName.Data())) ;
125                 TString libs = gSystem->GetLibraries() ;
126                 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
127                         pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
128                 } else {
129                         pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
130                 }
131                 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
132         }
133         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
134                 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
135         }
136         if (loader) 
137                 fLoader[iDet] = loader ;
138         return loader ;
139 }
140
141 //_____________________________________________________________________________
142 AliQADataMaker * AliQADataMakerSteer::GetQADataMaker(Int_t iDet)
143 {
144         // get the quality assurance data maker for a detector
145         
146         if (fQADataMaker[iDet]) 
147                 return fQADataMaker[iDet] ;
148         
149         AliQADataMaker * qadm = NULL ;
150         
151         if (fFirst) {
152                 // load the QA data maker object
153                 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
154                 TString detName = AliQA::GetDetName(iDet) ;
155                 TString qadmName = "Ali" + detName + "QADataMaker" ;
156
157                 // first check if a plugin is defined for the quality assurance data maker
158                 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
159                 // if not, add a plugin for it
160                 if (!pluginHandler) {
161                         AliDebug(1, Form("defining plugin for %s", qadmName.Data())) ;
162                         TString libs = gSystem->GetLibraries() ;
163                         if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
164                                 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
165                         } else {
166                                 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
167                         }
168                         pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
169                 }
170                 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
171                         qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
172                 }
173                 if (qadm) 
174                         fQADataMaker[iDet] = qadm ;
175         }
176         return qadm ;
177 }
178
179 //_____________________________________________________________________________
180 Bool_t AliQADataMakerSteer::Init(const AliQA::TASKINDEX taskIndex, const  char * fileName )
181 {
182         // Initialize the event source and QA data makers
183         
184         if (taskIndex == AliQA::kRAWS) {
185                 fRawReader = new AliRawReaderRoot(fileName) ; 
186             if ( ! fRawReader ) 
187                         return kFALSE ; 
188                 fRawReader->NextEvent() ; 
189                 fRunNumber = fRawReader->GetRunNumber() ; 
190                 fRawReader->RewindEvents();
191                 fNumberOfEvents = 999999 ;
192         } else if (taskIndex == AliQA::kESDS) {
193                 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
194                         TFile * esdFile = TFile::Open("AliESDs.root") ;
195                         fESDTree = dynamic_cast<TTree *> (esdFile->Get("esdTree")) ; 
196                         fESD     = new AliESDEvent() ;
197                         fESD->ReadFromTree(fESDTree) ;
198                         fESDTree->GetEntry(0) ; 
199                         fRunNumber = fESD->GetRunNumber() ; 
200                         fNumberOfEvents = fESDTree->GetEntries() ;
201                 } else {
202                         AliError("AliESDs.root not found") ; 
203                         return kFALSE ; 
204                 }                       
205         } else {
206                 if ( !InitRunLoader() ) {
207                         AliError("Run Loader not found") ; 
208                         return kFALSE ; 
209                 } else {
210                         if (fRunLoader->GetAliRun()) 
211                                 fRunNumber      = fRunLoader->GetAliRun()->GetRunNumber() ; 
212                         fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
213                 }
214         }
215                 
216         // Initialize all QA data makers for all detectors
217         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
218                 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
219                 if (!qadm) {
220                         AliWarning(Form("AliQADataMaker not found for %s", AliQA::GetDetName(iDet))) ; 
221                 } else {
222                         AliInfo(Form("Data Maker found for %s", qadm->GetName())) ; 
223                         qadm->Init(taskIndex, fRunNumber, GetQACycles(iDet)) ;
224                         qadm->StartOfCycle(taskIndex, fCycleSame) ;
225                 }
226         } 
227         fFirst = kFALSE ;
228         return kTRUE ; 
229 }
230
231 //_____________________________________________________________________________
232 Bool_t AliQADataMakerSteer::InitRunLoader()
233 {
234         // get or create the run loader
235         if (fRunLoader) {
236                 fCycleSame = kTRUE ; 
237                 return kTRUE ;
238         } 
239                 
240         if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
241     // load all base libraries to get the loader classes
242                 TString libs = gSystem->GetLibraries() ;
243                 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
244                         TString detName = AliQA::GetDetName(iDet) ;
245                         if (detName == "HLT") 
246                                 continue;
247                         if (libs.Contains("lib" + detName + "base.so")) 
248                                 continue;
249                         gSystem->Load("lib" + detName + "base.so");
250                 }
251                 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
252                 if (!fRunLoader) {
253                         AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
254                         return kFALSE;
255                 }
256                 fRunLoader->CdGAFile();
257                 if (fRunLoader->LoadgAlice() == 0) {
258                         gAlice = fRunLoader->GetAliRun();
259                 }
260                 if (!gAlice) {
261                         AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
262                         return kFALSE;
263                 }
264
265         } else {               // galice.root does not exist
266                 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
267                 return kFALSE;
268     }
269
270   return kTRUE;
271 }
272
273 //_____________________________________________________________________________
274 Bool_t AliQADataMakerSteer::Finish(const AliQA::TASKINDEX taskIndex) 
275 {
276         // write output to file for all detectors
277         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
278                 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
279                 if (qadm) {
280                         qadm->EndOfCycle(taskIndex) ; 
281                 }
282         }
283         return kTRUE ; 
284 }
285
286 //_____________________________________________________________________________
287 Bool_t AliQADataMakerSteer::Merge() 
288 {
289         // Merge 
290
291 }
292
293 //_____________________________________________________________________________
294 void AliQADataMakerSteer::Reset()
295 {
296         // Reset the default data members
297         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
298                 fLoader[iDet] = NULL;
299                 if (fQADataMaker[iDet]) {
300                         (fQADataMaker[iDet])->Reset() ; 
301 //                      delete fQADataMaker[iDet] ;
302 //                      fQADataMaker[iDet] = NULL ;
303                 }
304         }
305
306         delete fRawReader ;
307         fRawReader      = NULL ;
308
309         fCycleSame      = kFALSE ; 
310         fESD            = NULL ; 
311         fESDTree        = NULL ; 
312         fFirst          = kTRUE ;   
313         fNumberOfEvents = 0 ;  
314 }
315
316 //_____________________________________________________________________________
317 Bool_t AliQADataMakerSteer::Run(const AliQA::TASKINDEX taskIndex, const  char * fileName )
318 {
319         // Runs all the QA data Maker for every detector
320         Bool_t rv = kFALSE ;
321         
322         if ( !Init(taskIndex, fileName) ) 
323                 return kFALSE ; 
324                 
325     // Fill QA data in event loop 
326         for (UInt_t iEvent = 0 ; iEvent < fNumberOfEvents ; iEvent++) {
327                 // Get the event
328                 AliDebug(1, Form("processing event %d", iEvent));
329                 if ( taskIndex == AliQA::kRAWS ) {
330                         if ( !fRawReader->NextEvent() )
331                                 break ;
332                 } else if ( taskIndex == AliQA::kESDS ) {
333                         fESDTree->GetEntry(iEvent) ;
334                 } else {
335                         fRunLoader->GetEvent(iEvent);
336                 }
337                 // loop over detectors
338                 TObjArray* detArray = NULL ; 
339                 if (fRunLoader) // check if RunLoader exists 
340                         if ( fRunLoader->GetAliRun() ) // check if AliRun exists in gAlice.root
341                                 detArray = fRunLoader->GetAliRun()->Detectors() ;
342                 for (UInt_t iDet = 0 ; iDet < fgkNDetectors ; iDet++) {
343                         if (detArray) {
344                                 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
345                                 if (!det || !det->IsActive()) 
346                                         continue;
347                         }
348                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
349                         if (!qadm) {
350                                 rv = kFALSE ;
351                         } else {
352                                 if ( qadm->IsCycleDone() ) {
353                                         qadm->EndOfCycle(AliQA::kRAWS) ;
354                                         qadm->StartOfCycle(AliQA::kRAWS) ;
355                                 }
356                                 TTree * data ; 
357                                 switch (taskIndex) {
358                                         case AliQA::kRAWS :
359                                                 qadm->Exec(taskIndex, fRawReader) ; 
360                                         break ; 
361                                         case AliQA::kHITS :
362                                                 GetLoader(iDet)->LoadHits() ; 
363                                                 data = GetLoader(iDet)->TreeH() ; 
364                                                 if ( ! data ) {
365                                                         AliWarning(Form(" Hit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
366                                                 } else {
367                                                         qadm->Exec(taskIndex, data) ;
368                                                 } 
369                                         break ;
370                                         case AliQA::kSDIGITS :
371                                                 GetLoader(iDet)->LoadSDigits() ; 
372                                                 data = GetLoader(iDet)->TreeS() ; 
373                                                 if ( ! data ) {
374                                                         AliWarning(Form(" SDigit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
375                                                 } else {
376                                                         qadm->Exec(taskIndex, data) ; 
377                                                 }
378                                         break; 
379                                         case AliQA::kDIGITS :
380                                                 GetLoader(iDet)->LoadDigits() ; 
381                                                 data = GetLoader(iDet)->TreeD() ; 
382                                                 if ( ! data ) {
383                                                         AliWarning(Form(" Digit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
384                                                 } else {
385                                                         qadm->Exec(taskIndex, data) ;
386                                                 }
387                                         break; 
388                                         case AliQA::kRECPOINTS :
389                                                 GetLoader(iDet)->LoadRecPoints() ; 
390                                                 data = GetLoader(iDet)->TreeR() ; 
391                                                 if (!data) {
392                                                         AliWarning(Form("RecPoints not found for %s", AliQA::GetDetName(iDet))) ; 
393                                                 } else {
394                                                         qadm->Exec(taskIndex, data) ; 
395                                                 }
396                                         break; 
397                                         case AliQA::kTRACKSEGMENTS :
398                                         break; 
399                                         case AliQA::kRECPARTICLES :
400                                         break; 
401                                         case AliQA::kESDS :
402                                                 qadm->Exec(taskIndex, fESD) ;
403                                         break; 
404                                 } //task switch
405                                 qadm->Increment() ; 
406                         } //data maker exist
407                 } // detector loop
408         } // event loop 
409         // Save QA data for all detectors
410         rv = Finish(taskIndex) ;
411         return rv ; 
412 }