]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMakerSteer.cxx
Include cstdlib (gcc 4.3.0)
[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 //                                                                           //
19 // class for running the QA makers                                           //
20 //                                                                           //
21 //   AliQADataMakerSteer qas;                                                //
22 //   qas.Run(AliQA::kRAWS, rawROOTFileName);                                 //
23 //   qas.Run(AliQA::kHITS);                                                  //
24 //   qas.Run(AliQA::kSDIGITS);                                               //
25 //   qas.Run(AliQA::kDIGITS);                                                //
26 //   qas.Run(AliQA::kRECPOINTS);                                             //
27 //   qas.Run(AliQA::kESDS);                                                  //
28 //                                                                           //
29 ///////////////////////////////////////////////////////////////////////////////
30
31 #include <TKey.h>
32 #include <TFile.h>
33 #include <TFileMerger.h>
34 #include <TPluginManager.h>
35 #include <TROOT.h>
36 #include <TString.h>
37 #include <TSystem.h>
38
39 #include "AliCDBManager.h"
40 #include "AliCDBEntry.h"
41 #include "AliCDBId.h"
42 #include "AliCDBMetaData.h"
43 #include "AliESDEvent.h"
44 #include "AliGeomManager.h"
45 #include "AliHeader.h"
46 #include "AliLog.h"
47 #include "AliModule.h"
48 #include "AliQA.h"
49 #include "AliQADataMakerRec.h"
50 #include "AliQADataMakerSim.h"
51 #include "AliQADataMakerSteer.h" 
52 #include "AliRawReaderDate.h"
53 #include "AliRawReaderFile.h"
54 #include "AliRawReaderRoot.h"
55 #include "AliRun.h"
56 #include "AliRunLoader.h"
57
58 ClassImp(AliQADataMakerSteer) 
59
60 //_____________________________________________________________________________
61 AliQADataMakerSteer::AliQADataMakerSteer(const char* gAliceFilename, const char * name, const char * title) :
62         TNamed(name, title), 
63         fCurrentEvent(0),  
64         fCycleSame(kFALSE),
65         fDetectors("ALL"), 
66         fDetectorsW("ALL"), 
67         fESD(NULL), 
68         fESDTree(NULL),
69         fFirst(kTRUE),  
70         fGAliceFileName(gAliceFilename), 
71         fFirstEvent(0),        
72         fMaxEvents(0),        
73         fNumberOfEvents(999999), 
74     fRunNumber(0), 
75         fRawReader(NULL), 
76         fRawReaderDelete(kTRUE), 
77         fRunLoader(NULL),
78         fQADataMakers()
79
80 {
81         // default ctor
82         fMaxEvents = fNumberOfEvents ; 
83         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
84                 if (IsSelected(AliQA::GetDetName(iDet))) {
85                         fLoader[iDet]      = NULL ;
86                         fQADataMaker[iDet] = NULL ;
87                         fQACycles[iDet]    = 999999 ;
88                 }
89         }       
90 }
91
92 //_____________________________________________________________________________
93 AliQADataMakerSteer::AliQADataMakerSteer(const AliQADataMakerSteer & qas) : 
94         TNamed(qas), 
95         fCurrentEvent(qas.fCurrentEvent),  
96         fCycleSame(kFALSE),
97         fDetectors(qas.fDetectors), 
98         fDetectorsW(qas.fDetectorsW), 
99         fESD(NULL), 
100         fESDTree(NULL), 
101         fFirst(qas.fFirst),  
102         fGAliceFileName(qas.fGAliceFileName), 
103         fFirstEvent(qas.fFirstEvent),        
104         fMaxEvents(qas.fMaxEvents),        
105         fNumberOfEvents(qas.fNumberOfEvents), 
106     fRunNumber(qas.fRunNumber), 
107         fRawReader(NULL), 
108         fRawReaderDelete(kTRUE), 
109         fRunLoader(NULL),
110         fQADataMakers()
111 {
112         // cpy ctor
113         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
114                 fLoader[iDet]      = qas.fLoader[iDet] ;
115                 fQADataMaker[iDet] = qas.fQADataMaker[iDet] ;
116                 fQACycles[iDet]    = qas.fQACycles[iDet] ;      
117         }
118 }
119
120 //_____________________________________________________________________________
121 AliQADataMakerSteer & AliQADataMakerSteer::operator = (const AliQADataMakerSteer & qas) 
122 {
123         // assignment operator
124   this->~AliQADataMakerSteer() ;
125   new(this) AliQADataMakerSteer(qas) ;
126   return *this ;
127 }
128
129 //_____________________________________________________________________________
130 AliQADataMakerSteer::~AliQADataMakerSteer() 
131 {
132         // dtor
133         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
134                 if (IsSelected(AliQA::GetDetName(iDet))) {
135                   fLoader[iDet] = NULL;
136                   if (fQADataMaker[iDet]) {
137                           (fQADataMaker[iDet])->Finish() ; 
138                           fQADataMaker[iDet] = NULL ;
139                   }
140                 }
141         }
142         if (fRawReaderDelete) { 
143                 fRunLoader = NULL ;
144                 delete fRawReader ;
145                 fRawReader = NULL ;
146         }
147 }
148
149 //_____________________________________________________________________________
150 Bool_t AliQADataMakerSteer::DoIt(const AliQA::TASKINDEX_t taskIndex, const char * mode)
151 {
152         // Runs all the QA data Maker for every detector
153                 
154         Bool_t rv = kFALSE ;
155     // Fill QA data in event loop 
156         for (UInt_t iEvent = fFirstEvent ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
157                 fCurrentEvent++ ; 
158                 // Get the event
159                 if ( iEvent%10 == 0  ) 
160                         AliInfo(Form("processing event %d", iEvent));
161                 if ( taskIndex == AliQA::kRAWS ) {
162                         if ( !fRawReader->NextEvent() )
163                                 break ;
164                 } else if ( taskIndex == AliQA::kESDS ) {
165                         if ( fESDTree->GetEntry(iEvent) == 0 )
166                                 break ;
167                 } else {
168                         if ( fRunLoader->GetEvent(iEvent) != 0 )
169                                 break ;
170                 }
171                 // loop  over active loaders
172                 for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
173                         AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i)) ;
174                         if ( qadm->IsCycleDone() ) {
175                                 qadm->EndOfCycle(AliQA::kRAWS) ;
176                                 qadm->StartOfCycle(AliQA::kRAWS) ;
177                         }
178                         TTree * data ; 
179                         Int_t iDet = qadm->GetUniqueID();
180                         AliLoader* loader = GetLoader(iDet);
181                         switch (taskIndex) {
182                                 case AliQA::kRAWS :
183                                         qadm->Exec(taskIndex, fRawReader) ; 
184                                         break ; 
185                                 case AliQA::kHITS :
186                                         loader->LoadHits() ; 
187                                         data = loader->TreeH() ; 
188                                         if ( ! data ) {
189                                                 AliWarning(Form(" Hit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
190                                         } else {
191                                                 qadm->Exec(taskIndex, data) ;
192                                         } 
193                                         break ;
194                                         case AliQA::kSDIGITS :
195                                         loader->LoadSDigits() ; 
196                                         data = loader->TreeS() ; 
197                                         if ( ! data ) {
198                                                 AliWarning(Form(" SDigit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
199                                         } else {
200                                                 qadm->Exec(taskIndex, data) ; 
201                                         }
202                                         break; 
203                                         case AliQA::kDIGITS :
204                                         loader->LoadDigits() ; 
205                                         data = loader->TreeD() ; 
206                                         if ( ! data ) {
207                                                 AliWarning(Form(" Digit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
208                                         } else {
209                                                 qadm->Exec(taskIndex, data) ;
210                                         }
211                                         break; 
212                                         case AliQA::kRECPOINTS :
213                                         loader->LoadRecPoints() ; 
214                                         data = loader->TreeR() ; 
215                                         if (!data) {
216                                                 AliWarning(Form("RecPoints not found for %s", AliQA::GetDetName(iDet))) ; 
217                                         } else {
218                                                 qadm->Exec(taskIndex, data) ; 
219                                         }
220                                         break; 
221                                         case AliQA::kTRACKSEGMENTS :
222                                         break; 
223                                         case AliQA::kRECPARTICLES :
224                                         break; 
225                                         case AliQA::kESDS :
226                                         qadm->Exec(taskIndex, fESD) ;
227                                         break; 
228                                         case AliQA::kNTASKINDEX :
229                                         break; 
230                         } //task switch
231                         //qadm->Increment() ; 
232                 } // detector loop
233         } // event loop 
234 //      // Save QA data for all detectors
235         rv = Finish(taskIndex, mode) ;
236         
237         if ( taskIndex == AliQA::kRAWS ) 
238                 fRawReader->RewindEvents() ;
239
240         return rv ; 
241 }
242
243 //_____________________________________________________________________________
244 Bool_t AliQADataMakerSteer::Finish(const AliQA::TASKINDEX_t taskIndex, const char * /*mode*/) 
245 {
246         // write output to file for all detectors
247         for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
248                 AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i));
249                 qadm->EndOfCycle(taskIndex) ; 
250         }    
251         return kTRUE ; 
252 }
253
254 //_____________________________________________________________________________
255 TObjArray * AliQADataMakerSteer::GetFromOCDB(AliQA::DETECTORINDEX_t det, AliQA::TASKINDEX_t task, const char * year) const 
256 {
257         // Retrieve the list of QA data for a given detector and a given task 
258         TObjArray * rv = NULL ;
259         if ( !strlen(AliQA::GetQARefStorage()) ) { 
260                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
261                 return NULL ; 
262         }       
263         AliCDBManager* man = AliCDBManager::Instance() ; 
264         if ( ! man->IsDefaultStorageSet() ) {
265                 TString tmp(AliQA::GetQARefDefaultStorage()) ; 
266                 tmp.Append(year) ; 
267                 tmp.Append("/") ; 
268                 man->SetDefaultStorage(tmp.Data()) ;            
269                 man->SetSpecificStorage(Form("%s/*", AliQA::GetQAName()), AliQA::GetQARefStorage()) ;
270         }
271         TString detOCDBDir(Form("%s/%s/%s", AliQA::GetQAName(), AliQA::GetDetName((Int_t)det), AliQA::GetRefOCDBDirName())) ; 
272         AliInfo(Form("Retrieving reference data from %s/%s for %s", AliQA::GetQARefStorage(), detOCDBDir.Data(), AliQA::GetTaskName(task).Data())) ; 
273         AliCDBEntry* entry = man->Get(detOCDBDir.Data(), 0) ; //FIXME 0 --> Run Number
274         TList * listDetQAD = dynamic_cast<TList *>(entry->GetObject()) ;
275         if ( listDetQAD ) 
276                 rv = dynamic_cast<TObjArray *>(listDetQAD->FindObject(AliQA::GetTaskName(task))) ; 
277         return rv ; 
278 }
279
280 //_____________________________________________________________________________
281 AliLoader * AliQADataMakerSteer::GetLoader(Int_t iDet)
282 {
283         // get the loader for a detector
284
285         if ( !fRunLoader ) 
286                 return NULL ; 
287         
288         TString detName = AliQA::GetDetName(iDet) ;
289     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
290         if (fLoader[iDet]) 
291                 return fLoader[iDet] ;
292         
293         // load the QA data maker object
294         TPluginManager* pluginManager = gROOT->GetPluginManager() ;
295         TString loaderName = "Ali" + detName + "Loader" ;
296
297         AliLoader * loader = NULL ;
298         // first check if a plugin is defined for the quality assurance data maker
299         TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
300         // if not, add a plugin for it
301         if (!pluginHandler) {
302                 AliDebug(1, Form("defining plugin for %s", loaderName.Data())) ;
303                 TString libs = gSystem->GetLibraries() ;
304                 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
305                         pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
306                 } else {
307                         pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
308                 }
309                 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
310         }
311         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
312                 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
313         }
314         if (loader) 
315                 fLoader[iDet] = loader ;
316         return loader ;
317 }
318
319 //_____________________________________________________________________________
320 AliQADataMaker * AliQADataMakerSteer::GetQADataMaker(const Int_t iDet, const char * mode )
321 {
322         // get the quality assurance data maker for a detector
323         
324         if (fQADataMaker[iDet]) 
325                 return fQADataMaker[iDet] ;
326         
327         AliQADataMaker * qadm = NULL ;
328         
329         if (fFirst) {
330                 // load the QA data maker object
331                 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
332                 TString detName = AliQA::GetDetName(iDet) ;
333                 TString tmp(mode) ; 
334                 if (tmp.Contains("sim")) 
335                         tmp.ReplaceAll("s", "S") ; 
336                 else if (tmp.Contains("rec")) 
337                         tmp.ReplaceAll("r", "R") ; 
338                 TString qadmName = "Ali" + detName + "QADataMaker" + tmp ;
339
340                 // first check if a plugin is defined for the quality assurance data maker
341                 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
342                 // if not, add a plugin for it
343                 if (!pluginHandler) {
344                         AliDebug(1, Form("defining plugin for %s", qadmName.Data())) ;
345                         TString libs = gSystem->GetLibraries() ;
346                         if (libs.Contains("lib" + detName + mode + ".so") || (gSystem->Load("lib" + detName + mode + ".so") >= 0)) {
347                                 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
348                         } else {
349                                 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
350                         }
351                         pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
352                 }
353                 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
354                         qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
355                 }
356                 if (qadm) 
357                         fQADataMaker[iDet] = qadm ;
358         }
359         return qadm ;
360 }
361
362 //_____________________________________________________________________________
363 Bool_t AliQADataMakerSteer::Init(const AliQA::TASKINDEX_t taskIndex, const char * mode, const  char * input )
364 {
365         // Initialize the event source and QA data makers
366         
367         if (taskIndex == AliQA::kRAWS) { 
368                 if (!fRawReader) {
369                         fRawReader = AliRawReader::Create(input);
370                 }
371             if ( ! fRawReader ) 
372                         return kFALSE ; 
373                 fRawReaderDelete = kTRUE ; 
374                 fRawReader->NextEvent() ; 
375                 fRunNumber = fRawReader->GetRunNumber() ; 
376                 AliCDBManager::Instance()->SetRun(fRunNumber) ; 
377                 fRawReader->RewindEvents();
378                 fNumberOfEvents = 999999 ;
379                 if ( fMaxEvents < 0 ) 
380                         fMaxEvents = fNumberOfEvents ; 
381                 } else if (taskIndex == AliQA::kESDS) {
382                 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
383                         TFile * esdFile = TFile::Open("AliESDs.root") ;
384                         fESDTree = dynamic_cast<TTree *> (esdFile->Get("esdTree")) ; 
385                         if ( !fESDTree ) {
386                                 AliError("esdTree not found") ; 
387                                 return kFALSE ; 
388                         } else {
389                                 fESD     = new AliESDEvent() ;
390                                 fESD->ReadFromTree(fESDTree) ;
391                                 fESDTree->GetEntry(0) ; 
392                                 fRunNumber = fESD->GetRunNumber() ; 
393                                 fNumberOfEvents = fESDTree->GetEntries() ;
394                                 if ( fMaxEvents < 0 ) 
395                                         fMaxEvents = fNumberOfEvents ; 
396                         }
397                 } else {
398                         AliError("AliESDs.root not found") ; 
399                         return kFALSE ; 
400                 }                       
401                 } else {
402                 if ( !InitRunLoader() ) { 
403                         AliWarning("No Run Loader not found") ; 
404                 } else {
405                         fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
406                         if ( fMaxEvents < 0 ) 
407                                 fMaxEvents = fNumberOfEvents ; 
408
409                 }
410         }
411
412         // Get Detectors 
413         TObjArray* detArray = NULL ; 
414         if (fRunLoader) // check if RunLoader exists 
415                 if ( fRunLoader->GetAliRun() ) { // check if AliRun exists in gAlice.root
416                         detArray = fRunLoader->GetAliRun()->Detectors() ;
417                         fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
418                 }
419         // Build array of QA data makers for all detectors
420         fQADataMakers.Clear();
421         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
422                 if (IsSelected(AliQA::GetDetName(iDet))) {
423                         AliQADataMaker * qadm = GetQADataMaker(iDet, mode) ;
424                         if (!qadm) {
425                                 AliError(Form("AliQADataMaker not found for %s", AliQA::GetDetName(iDet))) ; 
426                                 fDetectorsW.ReplaceAll(AliQA::GetDetName(iDet), "") ; 
427                         } else {
428                                 AliDebug(1, Form("Data Maker found for %s", qadm->GetName())) ; 
429                                 // skip non active detectors
430                                 if (detArray) {
431                                         AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
432                                         if (!det || !det->IsActive())  
433                                                 continue ;
434                                 }
435                                 qadm->SetName(AliQA::GetDetName(iDet));
436                                 qadm->SetUniqueID(iDet);
437                                 fQADataMakers.Add(qadm);
438                         }
439                 }
440         } 
441         // Initialize all QA data makers for all detectors
442         fRunNumber = AliCDBManager::Instance()->GetRun() ; 
443         if ( !  AliGeomManager::GetGeometry() ) 
444                 AliGeomManager::LoadGeometry() ; 
445         for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
446                 AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i));
447                 qadm->Init(taskIndex, fRunNumber, GetQACycles(qadm->GetUniqueID())) ;
448                 qadm->StartOfCycle(taskIndex, fCycleSame) ;
449         } 
450         fFirst = kFALSE ;
451         return kTRUE ; 
452 }
453
454 //_____________________________________________________________________________
455 Bool_t AliQADataMakerSteer::InitRunLoader()
456 {
457         // get or create the run loader
458         if (fRunLoader) {
459                 fCycleSame = kTRUE ; 
460         } else {
461                 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
462                         // load all base libraries to get the loader classes
463                         TString libs = gSystem->GetLibraries() ;
464                         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
465                                 if (!IsSelected(AliQA::GetDetName(iDet))) 
466                                         continue ; 
467                                 TString detName = AliQA::GetDetName(iDet) ;
468                                 if (detName == "HLT") 
469                                         continue;
470                                 if (libs.Contains("lib" + detName + "base.so")) 
471                                         continue;
472                                 gSystem->Load("lib" + detName + "base.so");
473                         }
474                         fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
475                         if (!fRunLoader) {
476                                 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
477                                 return kFALSE;
478                         }
479                         fRunLoader->CdGAFile();
480                         if (fRunLoader->LoadgAlice() == 0) {
481                                 gAlice = fRunLoader->GetAliRun();
482                         }
483
484                         if (!gAlice) {
485                                 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
486                                 return kFALSE;
487                         }
488
489                 } else {               // galice.root does not exist
490                         AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
491                         return kFALSE;
492                 }
493         }
494
495         if (!fRunNumber) { 
496                 fRunLoader->LoadHeader();
497                 fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
498         }
499         return kTRUE;
500 }
501
502 //_____________________________________________________________________________
503 Bool_t AliQADataMakerSteer::IsSelected(const char * det) 
504 {
505         // check whether detName is contained in detectors
506         // if yes, it is removed from detectors
507         
508         Bool_t rv = kFALSE;
509         const TString detName(det) ;
510         // check if all detectors are selected
511         if (fDetectors.Contains("ALL")) {
512                 fDetectors = "ALL";
513                 rv = kTRUE;
514         } else if ((fDetectors.CompareTo(detName) == 0) ||
515                            fDetectors.BeginsWith(detName+" ") ||
516                            fDetectors.EndsWith(" "+detName) ||
517                            fDetectors.Contains(" "+detName+" ")) {
518                 //              fDetectors.ReplaceAll(detName, "");
519                 rv = kTRUE;
520         }
521         
522         // clean up the detectors string
523         //      while (fDetectors.Contains("  ")) 
524         //              fDetectors.ReplaceAll("  ", " ");
525         //      while (fDetectors.BeginsWith(" ")) 
526         //              fDetectors.Remove(0, 1);
527         //      while (fDetectors.EndsWith(" ")) 
528         //              fDetectors.Remove(fDetectors.Length()-1, 1);
529         
530         return rv ;
531 }
532
533 //_____________________________________________________________________________
534 Bool_t AliQADataMakerSteer::Merge(const Int_t runNumber) const
535 {
536         // Merge all the cycles from all detectors in one single file per run
537         TString cmd ;
538         if ( runNumber == -1 )
539                 cmd = Form(".! ls *%s*.*.*.root > tempo.txt", AliQA::GetQADataFileName()) ; 
540         else 
541                 cmd = Form(".! ls *%s*.%d.*.root > tempo.txt", AliQA::GetQADataFileName(), runNumber) ; 
542         gROOT->ProcessLine(cmd.Data()) ;
543         ifstream in("tempo.txt") ; 
544         const Int_t runMax = 10 ;  
545         TString file[AliQA::kNDET*runMax] ;
546         Int_t run[AliQA::kNDET*runMax] = {-1} ;
547         
548         Int_t index = 0 ; 
549         while ( 1 ) {
550                 in >> file[index] ; 
551                 if ( !in.good() ) 
552                         break ; 
553                 AliInfo(Form("index = %d file = %s", index, (file[index]).Data())) ; 
554                 index++ ;
555         }
556         
557         if ( index == 0 ) { 
558                 AliError(Form("run number %d not found", runNumber)) ; 
559                 return kFALSE ; 
560         }
561         
562         Int_t runIndex    = 0 ;  
563         Int_t runIndexMax = 0 ; 
564         TString stmp(Form(".%s.", AliQA::GetQADataFileName())) ; 
565         for (Int_t ifile = 0 ; ifile < index ; ifile++) {
566                 TString tmp(file[ifile]) ; 
567                 tmp.ReplaceAll(".root", "") ; 
568                 TString det = tmp(0, tmp.Index(".")) ; 
569                 tmp.Remove(0, tmp.Index(stmp)+4) ; 
570                 TString ttmp = tmp(0, tmp.Index(".")) ; 
571                 Int_t newRun = ttmp.Atoi() ;
572                 for (Int_t irun = 0; irun <= runIndexMax; irun++) {
573                         if (newRun == run[irun]) 
574                                 break ; 
575                         run[runIndex] = newRun ; 
576                         runIndex++ ; 
577                 }
578                 runIndexMax = runIndex ; 
579                 ttmp = tmp(tmp.Index(".")+1, tmp.Length()) ; 
580                 Int_t cycle = ttmp.Atoi() ;  
581                 AliDebug(1, Form("%s : det = %s run = %d cycle = %d \n", file[ifile].Data(), det.Data(), newRun, cycle)) ; 
582         }
583         for (Int_t irun = 0 ; irun < runIndexMax ; irun++) {
584                 TFileMerger merger ; 
585                 TString outFileName(Form("Merged.%s.%d.root",AliQA::GetQADataFileName(),run[irun]));            
586                 merger.OutputFile(outFileName.Data()) ; 
587                 for (Int_t ifile = 0 ; ifile < index-1 ; ifile++) {
588                         TString pattern(Form("%s.%d.", AliQA::GetQADataFileName(), run[irun])) ; 
589                         TString tmp(file[ifile]) ; 
590                         if (tmp.Contains(pattern)) {
591                                 merger.AddFile(tmp) ; 
592                         }
593                 }
594                 merger.Merge() ; 
595         }
596         
597         return kTRUE ; 
598 }
599
600 //_____________________________________________________________________________
601 void AliQADataMakerSteer::Reset(const Bool_t sameCycle)
602 {
603         // Reset the default data members
604
605         for (Int_t i = 0; i < fQADataMakers.GetEntriesFast() ; i++) {
606                 AliQADataMaker * qadm = static_cast<AliQADataMaker *>(fQADataMakers.At(i));
607                 qadm->Reset(sameCycle);
608         } 
609         if (fRawReaderDelete) { 
610                 delete fRawReader ;
611                 fRawReader      = NULL ;
612         }
613
614         fCycleSame      = sameCycle ; 
615         fESD            = NULL ; 
616         fESDTree        = NULL ; 
617         fFirst          = kTRUE ;   
618         fNumberOfEvents = 999999 ;  
619 }
620
621 //_____________________________________________________________________________
622 TString AliQADataMakerSteer::Run(const char * detectors, AliRawReader * rawReader, const Bool_t sameCycle) 
623 {
624         //Runs all the QA data Maker for Raws only
625         
626         fCycleSame       = sameCycle ;
627         fRawReader       = rawReader ;
628         fDetectors       = detectors ; 
629         fDetectorsW      = detectors ;  
630         
631         AliCDBManager* man = AliCDBManager::Instance() ; 
632
633         if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
634                 rawReader->NextEvent() ; 
635                 man->SetRun(fRawReader->GetRunNumber()) ;
636                 rawReader->RewindEvents() ;
637         }       
638         
639         if ( !Init(AliQA::kRAWS, "rec") ) 
640                 return kFALSE ; 
641         fRawReaderDelete = kFALSE ; 
642
643         DoIt(AliQA::kRAWS, "rec") ; 
644         return  fDetectorsW ;
645 }
646
647 //_____________________________________________________________________________
648 TString AliQADataMakerSteer::Run(const char * detectors, const char * fileName, const Bool_t sameCycle) 
649 {
650         //Runs all the QA data Maker for Raws only
651
652         fCycleSame       = sameCycle ;
653         fDetectors       = detectors ; 
654         fDetectorsW      = detectors ;  
655         
656         AliCDBManager* man = AliCDBManager::Instance() ; 
657         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
658                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
659                 if ( ! rl ) {
660                         AliFatal("galice.root file not found in current directory") ; 
661                 } else {
662                         rl->CdGAFile() ; 
663                         rl->LoadgAlice() ;
664                         if ( ! rl->GetAliRun() ) {
665                                 AliFatal("AliRun not found in galice.root") ;
666                         } else {
667                                 rl->LoadHeader() ;
668                                 man->SetRun(rl->GetHeader()->GetRun());
669                         }
670                 }
671         }
672         
673         if ( !Init(AliQA::kRAWS, "rec", fileName) ) 
674                 return kFALSE ; 
675         
676         DoIt(AliQA::kRAWS, "rec") ; 
677         return  fDetectorsW ;
678 }
679
680 //_____________________________________________________________________________
681 TString AliQADataMakerSteer::Run(const char * detectors, const AliQA::TASKINDEX_t taskIndex, Bool_t const sameCycle, const  char * fileName ) 
682 {
683         // Runs all the QA data Maker for every detector
684         
685         fCycleSame       = sameCycle ;
686         fDetectors       = detectors ; 
687         fDetectorsW      = detectors ;          
688         
689         TString mode ; 
690         if ( (taskIndex == AliQA::kHITS) || (taskIndex == AliQA::kSDIGITS) || (taskIndex == AliQA::kDIGITS) ) 
691                 mode = "sim" ; 
692         else if ( (taskIndex == AliQA::kRAWS) || (taskIndex == AliQA::kRECPOINTS) || (taskIndex == AliQA::kESDS) )
693                 mode = "rec" ; 
694         else {
695                 AliError(Form("%s not implemented", AliQA::GetTaskName(taskIndex).Data())) ; 
696                 return "" ;
697         }
698
699         AliCDBManager* man = AliCDBManager::Instance() ;        
700         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
701                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
702                 if ( ! rl ) {
703                         AliFatal("galice.root file not found in current directory") ; 
704                 } else {
705                         rl->CdGAFile() ; 
706                         rl->LoadgAlice() ;
707                         if ( ! rl->GetAliRun() ) {
708                                 AliInfo("AliRun not found in galice.root") ;
709                         } else {
710                                 rl->LoadHeader() ;
711                                 man->SetRun(rl->GetHeader()->GetRun()) ;
712                         }
713                 }
714         }
715         
716         if ( !Init(taskIndex, mode.Data(), fileName) ) 
717                 return kFALSE ; 
718
719         DoIt(taskIndex, mode.Data()) ;
720         
721         return fDetectorsW ;
722
723 }
724
725 //_____________________________________________________________________________
726 Bool_t AliQADataMakerSteer::Save2OCDB(const Int_t runNumber, const char * year, const Int_t cycleNumber, const char * detectors) const
727 {
728         // take the locasl QA data merge into a single file and save in OCDB 
729         Bool_t rv = kTRUE ; 
730         TString tmp(AliQA::GetQARefStorage()) ; 
731         if ( tmp.IsNull() ) { 
732                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
733                 return kFALSE ; 
734         }
735         if ( !(tmp.Contains(AliQA::GetLabLocalOCDB()) || tmp.Contains(AliQA::GetLabAliEnOCDB())) ) {
736                 AliError(Form("%s is a wrong storage, use %s or %s", AliQA::GetQARefStorage(), AliQA::GetLabLocalOCDB().Data(), AliQA::GetLabAliEnOCDB().Data())) ; 
737                 return kFALSE ; 
738         }
739         TString sdet(detectors) ; 
740         sdet.ToUpper() ;
741         TFile * inputFile ; 
742         if ( sdet.Contains("ALL") ) {
743                 rv = Merge(runNumber) ; 
744                 if ( ! rv )
745                         return kFALSE ; 
746                 TString inputFileName(Form("Merged.%s.%d.root", AliQA::GetQADataFileName(), runNumber)) ; 
747                 inputFile = TFile::Open(inputFileName.Data()) ; 
748                 rv = SaveIt2OCDB(runNumber, inputFile, year) ; 
749         } else {
750                 for (Int_t index = 0; index < AliQA::kNDET; index++) {
751                         if (sdet.Contains(AliQA::GetDetName(index))) {
752                                 TString inputFileName(Form("%s.%s.%d.%d.root", AliQA::GetDetName(index), AliQA::GetQADataFileName(), runNumber, cycleNumber)) ; 
753                                 inputFile = TFile::Open(inputFileName.Data()) ;                         
754                                 rv *= SaveIt2OCDB(runNumber, inputFile, year) ; 
755                         }
756                 }
757         }
758         return rv ; 
759 }
760
761 //_____________________________________________________________________________
762 Bool_t AliQADataMakerSteer::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const char * year) const
763 {
764         // reads the TH1 from file and adds it to appropriate list before saving to OCDB
765         Bool_t rv = kTRUE ;
766         AliInfo(Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQA::GetQARefStorage())) ; 
767         AliCDBManager* man = AliCDBManager::Instance() ; 
768         if ( ! man->IsDefaultStorageSet() ) {
769                 TString tmp( AliQA::GetQARefStorage() ) ; 
770                 if ( tmp.Contains(AliQA::GetLabLocalOCDB()) ) 
771                         man->SetDefaultStorage(AliQA::GetQARefStorage()) ;
772                 else {
773                         TString tmp1(AliQA::GetQARefDefaultStorage()) ; 
774                         tmp1.Append(year) ; 
775                         tmp1.Append("?user=alidaq") ; 
776                         man->SetDefaultStorage(tmp1.Data()) ; 
777                 }
778         }
779         man->SetSpecificStorage("*", AliQA::GetQARefStorage()) ; 
780         if(man->GetRun() < 0) 
781                 man->SetRun(runNumber);
782
783         AliCDBMetaData mdr ;
784         mdr.SetResponsible("yves schutz");
785
786         for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++) {
787                 TDirectory * detDir = inputFile->GetDirectory(AliQA::GetDetName(detIndex)) ; 
788                 if ( detDir ) {
789                         AliInfo(Form("Entering %s", detDir->GetName())) ;
790                         TString detOCDBDir(Form("%s/%s/%s", AliQA::GetDetName(detIndex), AliQA::GetRefOCDBDirName(), AliQA::GetRefDataDirName())) ; 
791                         AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity())  ;
792                         TList * listDetQAD = new TList() ;
793                         TString listName(Form("%s QA data Reference", AliQA::GetDetName(detIndex))) ; 
794                         mdr.SetComment("HMPID QA stuff");
795                         listDetQAD->SetName(listName) ; 
796                         TList * taskList = detDir->GetListOfKeys() ; 
797                         TIter nextTask(taskList) ; 
798                         TKey * taskKey ; 
799                         while ( (taskKey = dynamic_cast<TKey*>(nextTask())) ) {
800                                 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ; 
801                                 AliInfo(Form("Saving %s", taskDir->GetName())) ; 
802                                 TObjArray * listTaskQAD = new TObjArray(100) ; 
803                                 listTaskQAD->SetName(taskKey->GetName()) ;
804                                 listDetQAD->Add(listTaskQAD) ; 
805                                 TList * histList = taskDir->GetListOfKeys() ; 
806                                 TIter nextHist(histList) ; 
807                                 TKey * histKey ; 
808                                 while ( (histKey = dynamic_cast<TKey*>(nextHist())) ) {
809                                         TObject * odata = taskDir->Get(histKey->GetName()) ; 
810                                         if ( !odata ) {
811                                                 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
812                                         } else {
813                                                 AliInfo(Form("Adding %s", histKey->GetName())) ;
814                                                 if ( odata->IsA()->InheritsFrom("TH1") ) {
815                                                         AliInfo(Form("Adding %s", histKey->GetName())) ;
816                                                         TH1 * hdata = static_cast<TH1*>(odata) ; 
817                                                         listTaskQAD->Add(hdata) ; 
818                                                 }
819                                         }
820                                 }
821                         }
822                         man->Put(listDetQAD, idr, &mdr) ;
823                 }
824         }
825         return rv ; 
826 }       
827