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