]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMakerSteer.cxx
Adding VZERO track references (Brigitte)
[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 "AliCodeTimer.h"
44 #include "AliCorrQADataMakerRec.h"
45 #include "AliDetectorRecoParam.h"
46 #include "AliESDEvent.h"
47 #include "AliGeomManager.h"
48 #include "AliGlobalQADataMaker.h"
49 #include "AliHeader.h"
50 #include "AliLog.h"
51 #include "AliModule.h"
52 #include "AliQA.h"
53 #include "AliQADataMakerRec.h"
54 #include "AliQADataMakerSim.h"
55 #include "AliQADataMakerSteer.h" 
56 #include "AliRawReaderDate.h"
57 #include "AliRawReaderFile.h"
58 #include "AliRawReaderRoot.h"
59 #include "AliRun.h"
60 #include "AliRunLoader.h"
61 #include "AliRunTag.h"
62
63 ClassImp(AliQADataMakerSteer) 
64
65 //_____________________________________________________________________________
66 AliQADataMakerSteer::AliQADataMakerSteer(const Char_t * mode, const Char_t* gAliceFilename, const Char_t * name, const Char_t * title) :
67         TNamed(name, title), 
68         fCurrentEvent(0),  
69         fCycleSame(kFALSE),
70         fDetectors("ALL"), 
71         fDetectorsW("ALL"), 
72         fESD(NULL), 
73         fESDTree(NULL),
74         fGAliceFileName(gAliceFilename), 
75         fFirstEvent(0),        
76         fMaxEvents(0),   
77   fMode(mode), 
78         fNumberOfEvents(999999), 
79   fRecoParam(),
80         fRunNumber(0), 
81         fRawReader(NULL), 
82         fRawReaderDelete(kTRUE), 
83         fRunLoader(NULL), 
84   fTasks(""), 
85   fEventSpecie(AliRecoParam::kDefault)
86 {
87         // default ctor
88         fMaxEvents = fNumberOfEvents ; 
89         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
90                 if (IsSelected(AliQA::GetDetName(iDet))) {
91                         fLoader[iDet]      = NULL ;
92                         fQADataMaker[iDet] = NULL ;
93                         fQACycles[iDet]    = 999999 ;
94       fQAWriteExpert[iDet] = kTRUE ;
95                 }
96         }       
97 }
98
99 //_____________________________________________________________________________
100 AliQADataMakerSteer::AliQADataMakerSteer(const AliQADataMakerSteer & qas) : 
101         TNamed(qas), 
102         fCurrentEvent(qas.fCurrentEvent),  
103         fCycleSame(kFALSE),
104         fDetectors(qas.fDetectors), 
105         fDetectorsW(qas.fDetectorsW), 
106         fESD(NULL), 
107         fESDTree(NULL), 
108         fGAliceFileName(qas.fGAliceFileName), 
109         fFirstEvent(qas.fFirstEvent),        
110         fMaxEvents(qas.fMaxEvents),    
111         fMode(qas.fMode), 
112         fNumberOfEvents(qas.fNumberOfEvents), 
113         fRecoParam(),           
114         fRunNumber(qas.fRunNumber), 
115         fRawReader(NULL), 
116         fRawReaderDelete(kTRUE), 
117         fRunLoader(NULL), 
118   fTasks(qas.fTasks),
119   fEventSpecie(qas.fEventSpecie)
120 {
121         // cpy ctor
122         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
123                 fLoader[iDet]         = qas.fLoader[iDet] ;
124                 fQADataMaker[iDet]    = qas.fQADataMaker[iDet] ;
125                 fQACycles[iDet]       = qas.fQACycles[iDet] ;   
126     fQAWriteExpert[iDet] = qas.fQAWriteExpert[iDet] ;
127         }
128 }
129
130 //_____________________________________________________________________________
131 AliQADataMakerSteer & AliQADataMakerSteer::operator = (const AliQADataMakerSteer & qas) 
132 {
133         // assignment operator
134   this->~AliQADataMakerSteer() ;
135   new(this) AliQADataMakerSteer(qas) ;
136   return *this ;
137 }
138
139 //_____________________________________________________________________________
140 AliQADataMakerSteer::~AliQADataMakerSteer() 
141 {
142         // dtor
143         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
144                 if (IsSelected(AliQA::GetDetName(iDet))) {
145                   fLoader[iDet] = NULL;
146                   if (fQADataMaker[iDet]) {
147                           (fQADataMaker[iDet])->Finish() ; 
148                                 delete fQADataMaker[iDet] ;
149                   }
150                 }
151         }
152         if (fRawReaderDelete) { 
153                 fRunLoader = NULL ;
154                 delete fRawReader ;
155                 fRawReader = NULL ;
156         }
157 }
158
159 //_____________________________________________________________________________
160 Bool_t AliQADataMakerSteer::DoIt(const AliQA::TASKINDEX_t taskIndex)
161 {
162         // Runs all the QA data Maker for every detector
163                 
164         Bool_t rv = kFALSE ;
165     // Fill QA data in event loop 
166         for (UInt_t iEvent = fFirstEvent ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
167                 fCurrentEvent++ ; 
168                 // Get the event
169                 if ( iEvent%10 == 0  ) 
170                         AliInfo(Form("processing event %d", iEvent));
171                 if ( taskIndex == AliQA::kRAWS ) {
172                         if ( !fRawReader->NextEvent() )
173                                 break ;
174                 } else if ( taskIndex == AliQA::kESDS ) {
175                         if ( fESDTree->GetEntry(iEvent) == 0 )
176                                 break ;
177                 } else {
178                         if ( fRunLoader->GetEvent(iEvent) != 0 )
179                                 break ;
180                 }
181                 // loop  over active loaders
182                 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
183                         if (IsSelected(AliQA::GetDetName(iDet))) {
184                                 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
185                                 if (!qadm) continue; // This detector doesn't have any QA (for example, HLT)
186                                 if ( qadm->IsCycleDone() ) {
187           qadm->EndOfCycle(taskIndex) ;
188                                 }
189                                 TTree * data = NULL ; 
190                                 AliLoader* loader = GetLoader(qadm->GetUniqueID());
191                                 switch (taskIndex) {
192                                         case AliQA::kNULLTASKINDEX : 
193                                                 break ; 
194                                         case AliQA::kRAWS :
195                                                 qadm->Exec(taskIndex, fRawReader) ; 
196                                                 break ; 
197                                         case AliQA::kHITS :
198             if( loader ) {
199               loader->LoadHits() ; 
200               data = loader->TreeH() ; 
201               if ( ! data ) {
202                 AliWarning(Form(" Hit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
203                 break ; 
204               } 
205             } 
206             qadm->Exec(taskIndex, data) ;
207                                                 break ;
208                                                 case AliQA::kSDIGITS :
209             if( loader ) {      
210               loader->LoadSDigits() ; 
211               data = loader->TreeS() ; 
212               if ( ! data ) {
213                 AliWarning(Form(" SDigit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
214                 break ; 
215               } 
216             }
217             qadm->Exec(taskIndex, data) ; 
218                                                 break; 
219                                                 case AliQA::kDIGITS :
220             if( loader ) {      
221               loader->LoadDigits() ; 
222               data = loader->TreeD() ; 
223               if ( ! data ) {
224                 AliWarning(Form(" Digit Tree not found for  %s", AliQA::GetDetName(iDet))) ; 
225                 break ; 
226               } 
227             }
228             qadm->Exec(taskIndex, data) ;
229                                                 break; 
230                                                 case AliQA::kRECPOINTS :
231             if( loader ) {      
232               loader->LoadRecPoints() ; 
233               data = loader->TreeR() ; 
234               if (!data) {
235                 AliWarning(Form("RecPoints not found for %s", AliQA::GetDetName(iDet))) ; 
236                 break ; 
237               } 
238             }
239             qadm->Exec(taskIndex, data) ; 
240             break; 
241                                                 case AliQA::kTRACKSEGMENTS :
242                                                 break; 
243                                                 case AliQA::kRECPARTICLES :
244                                                 break; 
245                                                 case AliQA::kESDS :
246                                                 qadm->Exec(taskIndex, fESD) ;
247                                                 break; 
248                                                 case AliQA::kNTASKINDEX :
249                                                 break; 
250                                 } //task switch
251                         }
252                 } // detector loop
253     Increment() ; 
254         } // event loop 
255         // Save QA data for all detectors
256         rv = Finish(taskIndex) ;
257         
258         if ( taskIndex == AliQA::kRAWS ) 
259                 fRawReader->RewindEvents() ;
260
261         return rv ; 
262 }
263
264 //_____________________________________________________________________________
265 Bool_t AliQADataMakerSteer::Finish(const AliQA::TASKINDEX_t taskIndex) 
266 {
267         // write output to file for all detectors
268         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
269                 if (IsSelected(AliQA::GetDetName(iDet))) {
270                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
271                         if (qadm) 
272         qadm->EndOfCycle(taskIndex) ;
273                 }
274         }
275         return kTRUE ; 
276 }
277
278 //_____________________________________________________________________________
279 TObjArray * AliQADataMakerSteer::GetFromOCDB(AliQA::DETECTORINDEX_t det, AliQA::TASKINDEX_t task, const char * year) const 
280 {
281         // Retrieve the list of QA data for a given detector and a given task 
282         TObjArray * rv = NULL ;
283         if ( !strlen(AliQA::GetQARefStorage()) ) { 
284                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
285                 return NULL ; 
286         }       
287         AliCDBManager* man = AliCDBManager::Instance() ; 
288         if ( ! man->IsDefaultStorageSet() ) {
289                 TString tmp(AliQA::GetQARefDefaultStorage()) ; 
290                 tmp.Append(year) ; 
291                 tmp.Append("/") ; 
292                 man->SetDefaultStorage(tmp.Data()) ;            
293                 man->SetSpecificStorage(Form("%s/*", AliQA::GetQAName()), AliQA::GetQARefStorage()) ;
294         }
295         TString detOCDBDir(Form("%s/%s/%s", AliQA::GetQAName(), AliQA::GetDetName((Int_t)det), AliQA::GetRefOCDBDirName())) ; 
296         AliInfo(Form("Retrieving reference data from %s/%s for %s", AliQA::GetQARefStorage(), detOCDBDir.Data(), AliQA::GetTaskName(task).Data())) ; 
297         AliCDBEntry* entry = man->Get(detOCDBDir.Data(), 0) ; //FIXME 0 --> Run Number
298         TList * listDetQAD = dynamic_cast<TList *>(entry->GetObject()) ;
299         if ( listDetQAD ) 
300                 rv = dynamic_cast<TObjArray *>(listDetQAD->FindObject(AliQA::GetTaskName(task))) ; 
301         return rv ; 
302 }
303
304 //_____________________________________________________________________________
305 AliLoader * AliQADataMakerSteer::GetLoader(Int_t iDet)
306 {
307         // get the loader for a detector
308
309         if ( !fRunLoader || iDet == AliQA::kCORR) 
310                 return NULL ; 
311         
312         TString detName = AliQA::GetDetName(iDet) ;
313     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
314         if (fLoader[iDet]) 
315                 return fLoader[iDet] ;
316         
317         // load the QA data maker object
318         TPluginManager* pluginManager = gROOT->GetPluginManager() ;
319         TString loaderName = "Ali" + detName + "Loader" ;
320
321         AliLoader * loader = NULL ;
322         // first check if a plugin is defined for the quality assurance data maker
323         TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
324         // if not, add a plugin for it
325         if (!pluginHandler) {
326                 AliDebug(1, Form("defining plugin for %s", loaderName.Data())) ;
327                 TString libs = gSystem->GetLibraries() ;
328                 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
329                         pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
330                 } else {
331                         pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
332                 }
333                 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
334         }
335         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
336                 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
337         }
338         if (loader) 
339                 fLoader[iDet] = loader ;
340         return loader ;
341 }
342
343 //_____________________________________________________________________________
344 AliQA * AliQADataMakerSteer::GetQA(UInt_t run, UInt_t evt) 
345 {
346 // retrieves the QA object stored in a file named "Run{run}.Event{evt}_1.ESD.tag.root"  
347   char * fileName = Form("Run%d.Event%d_1.ESD.tag.root", run, evt) ; 
348   TFile * tagFile = TFile::Open(fileName) ;
349   if ( !tagFile ) {
350     AliError(Form("File %s not found", fileName)) ;
351     return NULL ; 
352   }
353   TTree * tagTree = dynamic_cast<TTree *>(tagFile->Get("T")) ; 
354   if ( !tagTree ) {
355     AliError(Form("Tree T not found in %s", fileName)) ; 
356     tagFile->Close() ; 
357     return NULL ; 
358   }
359   AliRunTag * tag = new AliRunTag ; 
360   tagTree->SetBranchAddress("AliTAG", &tag) ; 
361   tagTree->GetEntry(evt) ; 
362   AliQA * qa = AliQA::Instance(tag->GetQALength(), tag->GetQA(), tag->GetESLength(), tag->GetEventSpecies()) ; 
363   tagFile->Close() ; 
364   return qa ; 
365 }
366
367 //_____________________________________________________________________________
368 AliQADataMaker * AliQADataMakerSteer::GetQADataMaker(const Int_t iDet) 
369 {
370         // get the quality assurance data maker for a detector
371         
372         if (fQADataMaker[iDet]) {
373     fQADataMaker[iDet]->SetEventSpecie(fEventSpecie) ; 
374                 return fQADataMaker[iDet] ;
375   }
376         
377         AliQADataMaker * qadm = NULL ;
378         
379         if (iDet == AliQA::kGLOBAL) { //Global QA
380                 qadm = new AliGlobalQADataMaker();
381                 qadm->SetName(AliQA::GetDetName(iDet));
382                 qadm->SetUniqueID(iDet);
383                 fQADataMaker[iDet] = qadm;
384     qadm->SetEventSpecie(fEventSpecie) ; 
385                 return qadm;
386         }
387
388         if (iDet == AliQA::kCORR) { //the data maker for correlations among detectors
389     qadm = new AliCorrQADataMakerRec(fQADataMaker) ; 
390                 qadm->SetName(AliQA::GetDetName(iDet));
391                 qadm->SetUniqueID(iDet);
392                 fQADataMaker[iDet] = qadm;
393     qadm->SetEventSpecie(fEventSpecie) ; 
394                 return qadm;
395   }
396
397         // load the QA data maker object
398         TPluginManager* pluginManager = gROOT->GetPluginManager() ;
399         TString detName = AliQA::GetDetName(iDet) ;
400         TString tmp(fMode) ; 
401         if (tmp.Contains("sim")) 
402                 tmp.ReplaceAll("s", "S") ; 
403         else if (tmp.Contains("rec")) 
404                 tmp.ReplaceAll("r", "R") ; 
405
406         TString qadmName = "Ali" + detName + "QADataMaker" + tmp ;
407
408         // first check if a plugin is defined for the quality assurance data maker
409         TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
410         // if not, add a plugin for it
411         if (!pluginHandler) {
412                 AliDebug(1, Form("defining plugin for %s", qadmName.Data())) ;
413                 TString libs = gSystem->GetLibraries() ;
414                 if (libs.Contains("lib" + detName + fMode + ".so") || (gSystem->Load("lib" + detName + fMode + ".so") >= 0)) {
415                         pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
416                 } else {
417                         pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
418                 }
419                 pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
420         }
421         if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
422                 qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
423         }
424         if (qadm) {
425                 qadm->SetName(AliQA::GetDetName(iDet));
426                 qadm->SetUniqueID(iDet);
427                 fQADataMaker[iDet] = qadm ;
428     qadm->SetEventSpecie(fEventSpecie) ; 
429         }
430
431   return qadm ;
432 }
433
434 //_____________________________________________________________________________
435 void  AliQADataMakerSteer::EndOfCycle(TObjArray * detArray) 
436 {
437         // End of cycle QADataMakers 
438         
439         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
440                 if (IsSelected(AliQA::GetDetName(iDet))) {
441                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
442                         if (!qadm) 
443                                 continue ;      
444                         // skip non active detectors
445                         if (detArray) {
446                                 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
447                                 if (!det || !det->IsActive())  
448                                         continue ;
449                         }
450                         for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
451                                 if ( fTasks.Contains(Form("%d", taskIndex)) ) 
452                                         qadm->EndOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex))) ;
453                         }
454                         qadm->Finish();
455                 }
456         }
457 }
458
459 //_____________________________________________________________________________
460 void  AliQADataMakerSteer::EndOfCycle(TString detectors) 
461 {
462         // End of cycle QADataMakers 
463         
464         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
465                 if (IsSelected(AliQA::GetDetName(iDet))) {
466                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
467                         if (!qadm) 
468                                 continue ;      
469                         // skip non active detectors
470       if (!detectors.Contains(AliQA::GetDetName(iDet))) 
471         continue ;
472                 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
473                                 if ( fTasks.Contains(Form("%d", taskIndex)) ) 
474                                         qadm->EndOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex))) ;
475                         }
476                         qadm->Finish();
477                 }
478         }
479 }
480
481 //_____________________________________________________________________________
482 void AliQADataMakerSteer::Increment()
483 {
484   // Increments the cycle counter for all QA Data Makers
485         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
486                 if (IsSelected(AliQA::GetDetName(iDet))) {
487                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
488                         if (qadm) 
489         qadm->Increment() ;
490     }
491   }
492 }
493   
494 //_____________________________________________________________________________
495 Bool_t AliQADataMakerSteer::Init(const AliQA::TASKINDEX_t taskIndex, const  char * input )
496 {
497         // Initialize the event source and QA data makers
498         
499         fTasks += Form("%d", taskIndex) ; 
500
501         if (taskIndex == AliQA::kRAWS) { 
502                 if (!fRawReader) {
503                         fRawReader = AliRawReader::Create(input);
504                 }
505                 if ( ! fRawReader ) 
506                         return kFALSE ; 
507                 fRawReaderDelete = kTRUE ; 
508                 fRawReader->NextEvent() ; 
509                 fRunNumber = fRawReader->GetRunNumber() ; 
510                 AliCDBManager::Instance()->SetRun(fRunNumber) ; 
511                 fRawReader->RewindEvents();
512                 fNumberOfEvents = 999999 ;
513                 if ( fMaxEvents < 0 ) 
514                         fMaxEvents = fNumberOfEvents ; 
515                 } else if (taskIndex == AliQA::kESDS) {
516                         fTasks = AliQA::GetTaskName(AliQA::kESDS) ; 
517       if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
518         TFile * esdFile = TFile::Open("AliESDs.root") ;
519         fESDTree = dynamic_cast<TTree *> (esdFile->Get("esdTree")) ; 
520         if ( !fESDTree ) {
521           AliError("esdTree not found") ; 
522           return kFALSE ; 
523         } else {
524           fESD     = new AliESDEvent() ;
525           fESD->ReadFromTree(fESDTree) ;
526           fESDTree->GetEntry(0) ; 
527           fRunNumber = fESD->GetRunNumber() ; 
528           fNumberOfEvents = fESDTree->GetEntries() ;
529           if ( fMaxEvents < 0 ) 
530             fMaxEvents = fNumberOfEvents ; 
531         }
532       } else {
533         AliError("AliESDs.root not found") ; 
534         return kFALSE ; 
535       }                 
536     } else {
537       if ( !InitRunLoader() ) { 
538         AliWarning("No Run Loader not found") ; 
539       } else {
540         fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
541         if ( fMaxEvents < 0 ) 
542           fMaxEvents = fNumberOfEvents ; 
543       }
544     }
545
546   // Get Detectors 
547   TObjArray* detArray = NULL ; 
548         if (fRunLoader) // check if RunLoader exists 
549                 if ( fRunLoader->GetAliRun() ) { // check if AliRun exists in gAlice.root
550                         detArray = fRunLoader->GetAliRun()->Detectors() ;
551                         fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
552                 }
553
554         // Initialize all QA data makers for all detectors
555         fRunNumber = AliCDBManager::Instance()->GetRun() ; 
556         if ( !  AliGeomManager::GetGeometry() ) 
557                 AliGeomManager::LoadGeometry() ; 
558         
559         InitQADataMaker(fRunNumber, detArray) ; //, fCycleSame, kTRUE, detArray) ; 
560         return kTRUE ; 
561 }
562
563 //_____________________________________________________________________________
564 void  AliQADataMakerSteer::InitQADataMaker(UInt_t run, TObjArray * detArray) 
565 {
566         // Initializes The QADataMaker for all active detectors and for all active tasks 
567         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
568                 if (IsSelected(AliQA::GetDetName(iDet))) {
569                         AliQADataMaker * qadm = GetQADataMaker(iDet) ;
570                         if (!qadm) {
571                                 AliError(Form("AliQADataMaker not found for %s", AliQA::GetDetName(iDet))) ; 
572                                 fDetectorsW.ReplaceAll(AliQA::GetDetName(iDet), "") ; 
573                         } else {
574         if (fQAWriteExpert[iDet])
575           qadm->SetWriteExpert() ; 
576                                 AliDebug(1, Form("Data Maker found for %s", qadm->GetName())) ; 
577                                 // skip non active detectors
578                                 if (detArray) {
579                                         AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQA::GetDetName(iDet))) ;
580                                         if (!det || !det->IsActive())  
581                                                 continue ;
582                                 }
583                                 if (fQAWriteExpert[iDet]) qadm->SetWriteExpert() ; 
584               // Set default reco params
585         Bool_t sameCycle = kFALSE ; 
586                                 for (UInt_t taskIndex = 0; taskIndex < AliQA::kNTASKINDEX; taskIndex++) {
587                                         if ( fTasks.Contains(Form("%d", taskIndex)) ) {
588                                                 qadm->Init(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), GetQACycles(qadm->GetUniqueID())) ;
589             qadm->StartOfCycle(AliQA::GetTaskIndex(AliQA::GetTaskName(taskIndex)), run,  sameCycle) ;
590             sameCycle = kTRUE ;
591                                         }
592                                 }
593                         }
594                 }
595         }
596 }
597
598
599 //_____________________________________________________________________________
600 Bool_t AliQADataMakerSteer::InitRunLoader()
601 {
602         // get or create the run loader
603         if (fRunLoader) {
604                 fCycleSame = kTRUE ; 
605         } else {
606                 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
607                         // load all base libraries to get the loader classes
608                         TString libs = gSystem->GetLibraries() ;
609                         for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
610                                 if (!IsSelected(AliQA::GetDetName(iDet))) 
611                                         continue ; 
612                                 TString detName = AliQA::GetDetName(iDet) ;
613                                 if (detName == "HLT") 
614                                         continue;
615                                 if (libs.Contains("lib" + detName + "base.so")) 
616                                         continue;
617                                 gSystem->Load("lib" + detName + "base.so");
618                         }
619                         fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
620                         if (!fRunLoader) {
621                                 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
622                                 return kFALSE;
623                         }
624                         fRunLoader->CdGAFile();
625                         if (fRunLoader->LoadgAlice() == 0) {
626                                 gAlice = fRunLoader->GetAliRun();
627                         }
628
629                         if (!gAlice) {
630                                 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
631                                 return kFALSE;
632                         }
633
634                 } else {               // galice.root does not exist
635                         AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
636                         return kFALSE;
637                 }
638         }
639
640         if (!fRunNumber) { 
641                 fRunLoader->LoadHeader();
642                 fRunNumber = fRunLoader->GetHeader()->GetRun() ; 
643         }
644         return kTRUE;
645 }
646
647 //_____________________________________________________________________________
648 Bool_t AliQADataMakerSteer::IsSelected(const char * det) 
649 {
650         // check whether detName is contained in detectors
651         // if yes, it is removed from detectors
652         
653         Bool_t rv = kFALSE;
654         const TString detName(det) ;
655   // always activates Correlation
656   if ( detName.Contains(AliQA::GetDetName(AliQA::kCORR))) {
657     rv = kTRUE ; 
658   } else {
659     // check if all detectors are selected
660     if (fDetectors.Contains("ALL")) {
661       fDetectors = "ALL";
662       rv = kTRUE;
663     } else if ((fDetectors.CompareTo(detName) == 0) ||
664                fDetectors.BeginsWith(detName+" ") ||
665                fDetectors.EndsWith(" "+detName) ||
666                fDetectors.Contains(" "+detName+" ")) {
667       rv = kTRUE;
668     }
669   }
670         return rv ;
671 }
672
673 //_____________________________________________________________________________
674 Bool_t AliQADataMakerSteer::Merge(const Int_t runNumber) const
675 {
676         // Merge data from all the cycles from all detectors in one single file per run
677         // Merge the QA results from all the data chunks in one run 
678  Bool_t rv = MergeData(runNumber) ; 
679  rv *= MergeResults(runNumber) ;
680  return rv ; 
681 }
682         
683         
684 //_____________________________________________________________________________
685 Bool_t AliQADataMakerSteer::MergeData(const Int_t runNumber) const
686 {
687         // Merge all the cycles from all detectors in one single file per run
688         TString cmd ;
689         if (runNumber == -1) 
690                 cmd = Form(".! ls *%s*.%d.root > tempo.txt", AliQA::GetQADataFileName(), runNumber) ; 
691         else 
692                 cmd = Form(".! ls *%s*.*.root > tempo.txt", AliQA::GetQADataFileName()) ; 
693         gROOT->ProcessLine(cmd.Data()) ;
694         ifstream in("tempo.txt") ; 
695         const Int_t runMax = 10 ;  
696         TString file[AliQA::kNDET*runMax] ;
697         
698         Int_t index = 0 ; 
699         while ( 1 ) {
700                 in >> file[index] ; 
701                 if ( !in.good() ) 
702                         break ; 
703                 AliInfo(Form("index = %d file = %s", index, (file[index]).Data())) ; 
704                 index++ ;
705         }
706         
707         if ( index == 0 ) { 
708                 AliError(Form("run number %d not found", runNumber)) ; 
709                 return kFALSE ; 
710         }
711
712   TFileMerger merger ; 
713   TString outFileName(Form("Merged.%s.Data.%d.root",AliQA::GetQADataFileName(),runNumber)); 
714   merger.OutputFile(outFileName.Data()) ; 
715   for (Int_t ifile = 0 ; ifile < index-1 ; ifile++) {
716     TString pattern(Form("%s.%d.", AliQA::GetQADataFileName(), runNumber)); 
717     TString tmp(file[ifile]) ; 
718     if (tmp.Contains(pattern)) {
719       merger.AddFile(tmp) ; 
720     }
721         }
722   merger.Merge() ; 
723         return kTRUE ; 
724 }
725
726 //_____________________________________________________________________________
727 Bool_t AliQADataMakerSteer::MergeResults(const Int_t runNumber) const
728 {
729         // Merge the QA result from all the data chunks in a run 
730         TString cmd ;
731         cmd = Form(".! ls %s*.root > tempo.txt", AliQA::GetQADataFileName()) ; 
732         gROOT->ProcessLine(cmd.Data()) ;
733         ifstream in("tempo.txt") ; 
734         const Int_t chunkMax = 100 ;  
735         TString fileList[chunkMax] ;
736         
737         Int_t index = 0 ; 
738         while ( 1 ) {
739                 TString file ; 
740                 in >> fileList[index] ; 
741                 if ( !in.good() ) 
742                         break ; 
743                 AliInfo(Form("index = %d file = %s", index, (fileList[index].Data()))) ; 
744                 index++ ;
745         }
746         
747         if ( index == 0 ) { 
748                 AliError("No QA Result File found") ; 
749                 return kFALSE ; 
750         }
751         
752         TFileMerger merger ; 
753         TString outFileName(Form("Merged.%s.Result.%d.root", AliQA::GetQADataFileName(), runNumber));           
754         merger.OutputFile(outFileName.Data()) ; 
755         for (Int_t ifile = 0 ; ifile < index ; ifile++) {
756                 TString file = fileList[ifile] ; 
757                 merger.AddFile(file) ; 
758         }
759         merger.Merge() ; 
760         
761         return kTRUE ; 
762 }
763
764 //_____________________________________________________________________________
765 void AliQADataMakerSteer::Reset(const Bool_t sameCycle)
766 {
767         // Reset the default data members
768
769         for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
770                 if (IsSelected(AliQA::GetDetName(iDet))) {
771                         AliQADataMaker * qadm = GetQADataMaker(iDet);
772                         qadm->Reset();
773                 }
774         } 
775         if (fRawReaderDelete) { 
776                 delete fRawReader ;
777                 fRawReader      = NULL ;
778         }
779
780         fCycleSame      = sameCycle ; 
781         fESD            = NULL ; 
782         fESDTree        = NULL ; 
783         //fFirst          = kTRUE ;   
784         fNumberOfEvents = 999999 ;  
785 }
786
787 //_____________________________________________________________________________
788 TString AliQADataMakerSteer::Run(const char * detectors, AliRawReader * rawReader, const Bool_t sameCycle) 
789 {
790         //Runs all the QA data Maker for Raws only
791         
792         fCycleSame       = sameCycle ;
793         fRawReader       = rawReader ;
794         fDetectors       = detectors ; 
795         fDetectorsW      = detectors ;  
796         
797         AliCDBManager* man = AliCDBManager::Instance() ; 
798
799         if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
800                 rawReader->NextEvent() ; 
801                 man->SetRun(fRawReader->GetRunNumber()) ;
802                 rawReader->RewindEvents() ;
803         }       
804         
805         if (!fCycleSame) 
806     if ( !Init(AliQA::kRAWS) ) 
807       return kFALSE ; 
808   fRawReaderDelete = kFALSE ; 
809
810         DoIt(AliQA::kRAWS) ; 
811         return  fDetectorsW ;
812 }
813
814 //_____________________________________________________________________________
815 TString AliQADataMakerSteer::Run(const char * detectors, const char * fileName, const Bool_t sameCycle) 
816 {
817         //Runs all the QA data Maker for Raws only
818
819         fCycleSame       = sameCycle ;
820         fDetectors       = detectors ; 
821         fDetectorsW      = detectors ;  
822         
823         AliCDBManager* man = AliCDBManager::Instance() ; 
824         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
825                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
826                 if ( ! rl ) {
827                         AliFatal("galice.root file not found in current directory") ; 
828                 } else {
829                         rl->CdGAFile() ; 
830                         rl->LoadgAlice() ;
831                         if ( ! rl->GetAliRun() ) {
832                                 AliFatal("AliRun not found in galice.root") ;
833                         } else {
834                                 rl->LoadHeader() ;
835                                 man->SetRun(rl->GetHeader()->GetRun());
836                         }
837                 }
838         }
839         
840         if (!fCycleSame) 
841     if ( !Init(AliQA::kRAWS, fileName) ) 
842       return kFALSE ; 
843         
844         DoIt(AliQA::kRAWS) ; 
845         return  fDetectorsW ;
846 }
847
848 //_____________________________________________________________________________
849 TString AliQADataMakerSteer::Run(const char * detectors, const AliQA::TASKINDEX_t taskIndex, Bool_t const sameCycle, const  char * fileName ) 
850 {
851         // Runs all the QA data Maker for every detector
852         
853         fCycleSame       = sameCycle ;
854         fDetectors       = detectors ; 
855         fDetectorsW      = detectors ;          
856         
857         AliCDBManager* man = AliCDBManager::Instance() ;        
858         if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
859                 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
860                 if ( ! rl ) {
861                         AliFatal("galice.root file not found in current directory") ; 
862                 } else {
863                         rl->CdGAFile() ; 
864                         rl->LoadgAlice() ;
865                         if ( ! rl->GetAliRun() ) {
866                                 AliInfo("AliRun not found in galice.root") ;
867                         } else {
868                                 rl->LoadHeader() ;
869                                 man->SetRun(rl->GetHeader()->GetRun()) ;
870                         }
871                 }
872         }
873         
874
875         if ( taskIndex == AliQA::kNULLTASKINDEX) { 
876                 for (UInt_t task = 0; task < AliQA::kNTASKINDEX; task++) {
877                         if ( fTasks.Contains(Form("%d", task)) ) {
878         if (!fCycleSame)
879           if ( !Init(AliQA::GetTaskIndex(AliQA::GetTaskName(task)), fileName) ) 
880             return kFALSE ;
881         DoIt(AliQA::GetTaskIndex(AliQA::GetTaskName(task))) ;
882                         }
883                 }
884         } else {
885     if (! fCycleSame )
886       if ( !Init(taskIndex, fileName) ) 
887         return kFALSE ; 
888       DoIt(taskIndex) ; 
889   }             
890         
891         return fDetectorsW ;
892
893 }
894
895 //_____________________________________________________________________________
896 void AliQADataMakerSteer::RunOneEvent(AliRawReader * rawReader) 
897 {
898         //Runs all the QA data Maker for Raws only and on one event only (event loop done by calling method)
899   if ( ! rawReader ) 
900     return ; 
901         AliCodeTimerAuto("") ;
902   if (fTasks.Contains(Form("%d", AliQA::kRAWS))){
903     for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
904       if (!IsSelected(AliQA::GetDetName(iDet))) 
905         continue;
906       AliQADataMaker *qadm = GetQADataMaker(iDet);  
907       if (!qadm) 
908         continue;
909       if ( qadm->IsCycleDone() ) {
910         qadm->EndOfCycle() ;
911       }
912       AliCodeTimerStart(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet))); 
913       qadm->SetEventSpecie(fEventSpecie) ; 
914                         qadm->Exec(AliQA::kRAWS, rawReader) ;
915       AliCodeTimerStop(Form("running RAW quality assurance data maker for %s", AliQA::GetDetName(iDet)));
916                 }
917   }
918 }
919
920 //_____________________________________________________________________________
921 void AliQADataMakerSteer::RunOneEvent(AliESDEvent *& esd) 
922 {
923         //Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
924         
925   AliCodeTimerAuto("") ;
926   if (fTasks.Contains(Form("%d", AliQA::kESDS))) {
927     for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
928       if (!IsSelected(AliQA::GetDetName(iDet))) 
929         continue;
930       AliQADataMaker *qadm = GetQADataMaker(iDet);  
931       if (!qadm) 
932         continue;
933       if ( qadm->IsCycleDone() ) {
934         qadm->EndOfCycle() ;
935       }
936       AliCodeTimerStart(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
937                         qadm->Exec(AliQA::kESDS, esd) ;
938       AliCodeTimerStop(Form("running ESD quality assurance data maker for %s", AliQA::GetDetName(iDet)));
939                 }
940         }
941 }
942
943 //_____________________________________________________________________________
944 void AliQADataMakerSteer::RunOneEventInOneDetector(Int_t det, TTree * tree) 
945 {
946         // Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
947         AliCodeTimerAuto("") ;
948   if (fTasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
949     if (IsSelected(AliQA::GetDetName(det))) {
950       AliQADataMaker *qadm = GetQADataMaker(det);  
951       if (qadm) { 
952         if ( qadm->IsCycleDone() ) {
953           qadm->EndOfCycle() ;
954         }
955         AliCodeTimerStart(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
956         qadm->Exec(AliQA::kRECPOINTS, tree) ;
957         AliCodeTimerStop(Form("running RecPoints quality assurance data maker for %s", AliQA::GetDetName(det)));
958       }
959     }
960   }
961 }
962
963 //_____________________________________________________________________________
964 Bool_t AliQADataMakerSteer::Save2OCDB(const Int_t runNumber, AliRecoParam::EventSpecie_t es, const char * year, const char * detectors) const
965 {
966         // take the locasl QA data merge into a single file and save in OCDB 
967         Bool_t rv = kTRUE ; 
968         TString tmp(AliQA::GetQARefStorage()) ; 
969         if ( tmp.IsNull() ) { 
970                 AliError("No storage defined, use AliQA::SetQARefStorage") ; 
971                 return kFALSE ; 
972         }
973         if ( !(tmp.Contains(AliQA::GetLabLocalOCDB()) || tmp.Contains(AliQA::GetLabAliEnOCDB())) ) {
974                 AliError(Form("%s is a wrong storage, use %s or %s", AliQA::GetQARefStorage(), AliQA::GetLabLocalOCDB().Data(), AliQA::GetLabAliEnOCDB().Data())) ; 
975                 return kFALSE ; 
976         }
977         TString sdet(detectors) ; 
978         sdet.ToUpper() ;
979         TFile * inputFile ; 
980         if ( sdet.Contains("ALL") ) {
981                 rv = Merge(runNumber) ; 
982                 if ( ! rv )
983                         return kFALSE ; 
984                 TString inputFileName(Form("Merged.%s.Data.%d.root", AliQA::GetQADataFileName(), runNumber)) ; 
985                 inputFile = TFile::Open(inputFileName.Data()) ; 
986                 rv = SaveIt2OCDB(runNumber, inputFile, year, es) ; 
987         } else {
988                 for (Int_t index = 0; index < AliQA::kNDET; index++) {
989                         if (sdet.Contains(AliQA::GetDetName(index))) {
990                                 TString inputFileName(Form("%s.%s.%d.root", AliQA::GetDetName(index), AliQA::GetQADataFileName(), runNumber)) ; 
991                                 inputFile = TFile::Open(inputFileName.Data()) ;                         
992                                 rv *= SaveIt2OCDB(runNumber, inputFile, year, es) ; 
993                         }
994                 }
995         }
996         return rv ; 
997 }
998
999 //_____________________________________________________________________________
1000 Bool_t AliQADataMakerSteer::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const char * year, AliRecoParam::EventSpecie_t es) const
1001 {
1002         // reads the TH1 from file and adds it to appropriate list before saving to OCDB
1003         Bool_t rv = kTRUE ;
1004         AliInfo(Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQA::GetQARefStorage())) ; 
1005         AliCDBManager* man = AliCDBManager::Instance() ; 
1006         if ( ! man->IsDefaultStorageSet() ) {
1007                 TString tmp( AliQA::GetQARefStorage() ) ; 
1008                 if ( tmp.Contains(AliQA::GetLabLocalOCDB()) ) 
1009                         man->SetDefaultStorage(AliQA::GetQARefStorage()) ;
1010                 else {
1011                         TString tmp1(AliQA::GetQARefDefaultStorage()) ; 
1012                         tmp1.Append(year) ; 
1013                         tmp1.Append("?user=alidaq") ; 
1014                         man->SetDefaultStorage(tmp1.Data()) ; 
1015                 }
1016         }
1017         man->SetSpecificStorage("*", AliQA::GetQARefStorage()) ; 
1018         if(man->GetRun() < 0) 
1019                 man->SetRun(runNumber);
1020
1021         AliCDBMetaData mdr ;
1022         mdr.SetResponsible("yves schutz");
1023
1024         for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++) {
1025                 TDirectory * detDir = inputFile->GetDirectory(AliQA::GetDetName(detIndex)) ; 
1026                 if ( detDir ) {
1027                         AliInfo(Form("Entering %s", detDir->GetName())) ;
1028       AliQA::SetQARefDataDirName(es) ;
1029                         TString detOCDBDir(Form("%s/%s/%s", AliQA::GetDetName(detIndex), AliQA::GetRefOCDBDirName(), AliQA::GetRefDataDirName())) ; 
1030                         AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity())  ;
1031                         TList * listDetQAD = new TList() ;
1032                         TString listName(Form("%s QA data Reference", AliQA::GetDetName(detIndex))) ; 
1033                         mdr.SetComment(Form("%s QA stuff", AliQA::GetDetName(detIndex)));
1034                         listDetQAD->SetName(listName) ; 
1035                         TList * taskList = detDir->GetListOfKeys() ; 
1036                         TIter nextTask(taskList) ; 
1037                         TKey * taskKey ; 
1038                         while ( (taskKey = dynamic_cast<TKey*>(nextTask())) ) {
1039                                 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ; 
1040         TDirectory * esDir   = taskDir->GetDirectory(AliRecoParam::GetEventSpecieName(es)) ; 
1041                                 AliInfo(Form("Saving %s", esDir->GetName())) ; 
1042                                 TObjArray * listTaskQAD = new TObjArray(100) ; 
1043                                 listTaskQAD->SetName(Form("%s/%s", taskKey->GetName(), AliRecoParam::GetEventSpecieName(es))) ;
1044                                 listDetQAD->Add(listTaskQAD) ; 
1045                                 TList * histList = esDir->GetListOfKeys() ; 
1046                                 TIter nextHist(histList) ; 
1047                                 TKey * histKey ; 
1048                                 while ( (histKey = dynamic_cast<TKey*>(nextHist())) ) {
1049                                         TObject * odata = esDir->Get(histKey->GetName()) ; 
1050                                         if ( !odata ) {
1051                                                 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1052                                         } else {
1053             if ( AliQA::GetExpert() == histKey->GetName() ) {
1054               TDirectory * expertDir   = esDir->GetDirectory(histKey->GetName()) ; 
1055               TList * expertHistList = expertDir->GetListOfKeys() ; 
1056               TIter nextExpertHist(expertHistList) ; 
1057               TKey * expertHistKey ; 
1058               while ( (expertHistKey = dynamic_cast<TKey*>(nextExpertHist())) ) {
1059                 TObject * expertOdata = expertDir->Get(expertHistKey->GetName()) ; 
1060                 if ( !expertOdata ) {
1061                   AliError(Form("%s in %s/%s/Expert returns a NULL pointer !!", expertHistKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1062                 } else {
1063                   AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1064                   if ( expertOdata->IsA()->InheritsFrom("TH1") ) {
1065                     AliInfo(Form("Adding %s", expertHistKey->GetName())) ;
1066                     TH1 * hExpertdata = static_cast<TH1*>(expertOdata) ; 
1067                     listTaskQAD->Add(hExpertdata) ; 
1068                   }                  
1069                 }                
1070               }
1071             }
1072                                                 AliInfo(Form("Adding %s", histKey->GetName())) ;
1073                                                 if ( odata->IsA()->InheritsFrom("TH1") ) {
1074                                                         AliInfo(Form("Adding %s", histKey->GetName())) ;
1075                                                         TH1 * hdata = static_cast<TH1*>(odata) ; 
1076                                                         listTaskQAD->Add(hdata) ; 
1077                                                 }
1078                                         }
1079                                 }
1080                         }
1081                         man->Put(listDetQAD, idr, &mdr) ;
1082                 }
1083         }
1084         return rv ; 
1085 }       
1086
1087 //_____________________________________________________________________________
1088 void AliQADataMakerSteer::SetEventSpecie(AliRecoParam::EventSpecie_t es) 
1089 {
1090   // set the current event specie and inform AliQA that this event specie has been encountered
1091   fEventSpecie = es ;
1092   AliQA::Instance()->SetEventSpecie(es) ; 
1093 }
1094
1095 //_____________________________________________________________________________
1096 void AliQADataMakerSteer::SetRecoParam(const Int_t det, const AliDetectorRecoParam *par) 
1097 {
1098   // Set custom reconstruction parameters for a given detector
1099   // Single set of parameters for all the events
1100   GetQADataMaker(det)->SetRecoParam(par) ; 
1101 }
1102
1103