1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliQAManager.cxx 30894 2009-02-05 13:46:48Z schutz $ */
17 ///////////////////////////////////////////////////////////////////////////////
19 // class for running the QA makers //
21 // AliQAManager qas; //
22 // qas.Run(AliQAv1::kRAWS, rawROOTFileName); //
23 // qas.Run(AliQAv1::kHITS); //
24 // qas.Run(AliQAv1::kSDIGITS); //
25 // qas.Run(AliQAv1::kDIGITS); //
26 // qas.Run(AliQAv1::kRECPOINTS); //
27 // qas.Run(AliQAv1::kESDS); //
29 ///////////////////////////////////////////////////////////////////////////////
33 #include <TFileMerger.h>
35 #include <TGridCollection.h>
36 #include <TGridResult.h>
37 #include <TPluginManager.h>
42 #include "AliCDBManager.h"
43 #include "AliCDBEntry.h"
45 #include "AliCDBMetaData.h"
46 #include "AliCodeTimer.h"
47 #include "AliCorrQADataMakerRec.h"
48 #include "AliDetectorRecoParam.h"
49 #include "AliESDEvent.h"
50 #include "AliGeomManager.h"
51 #include "AliGlobalQADataMaker.h"
52 #include "AliHeader.h"
54 #include "AliModule.h"
56 #include "AliQADataMakerRec.h"
57 #include "AliQADataMakerSim.h"
58 #include "AliQAManager.h"
59 #include "AliRawReaderDate.h"
60 #include "AliRawReaderFile.h"
61 #include "AliRawReaderRoot.h"
63 #include "AliRunLoader.h"
64 #include "AliRunTag.h"
66 ClassImp(AliQAManager)
67 AliQAManager* AliQAManager::fgQAInstance = 0x0;
69 //_____________________________________________________________________________
70 AliQAManager::AliQAManager() :
82 fNumberOfEvents(999999),
86 fRawReaderDelete(kTRUE),
91 fMaxEvents = fNumberOfEvents ;
92 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
93 if (IsSelected(AliQAv1::GetDetName(iDet))) {
94 fLoader[iDet] = NULL ;
95 fQADataMaker[iDet] = NULL ;
96 fQACycles[iDet] = 999999 ;
102 //_____________________________________________________________________________
103 AliQAManager::AliQAManager(const Char_t * mode, const Char_t* gAliceFilename) :
111 fGAliceFileName(gAliceFilename),
115 fNumberOfEvents(999999),
119 fRawReaderDelete(kTRUE),
124 fMaxEvents = fNumberOfEvents ;
125 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
126 if (IsSelected(AliQAv1::GetDetName(iDet))) {
127 fLoader[iDet] = NULL ;
128 fQADataMaker[iDet] = NULL ;
129 fQACycles[iDet] = 999999 ;
135 //_____________________________________________________________________________
136 AliQAManager::AliQAManager(const AliQAManager & qas) :
138 fCurrentEvent(qas.fCurrentEvent),
140 fDetectors(qas.fDetectors),
141 fDetectorsW(qas.fDetectorsW),
144 fGAliceFileName(qas.fGAliceFileName),
145 fFirstEvent(qas.fFirstEvent),
146 fMaxEvents(qas.fMaxEvents),
148 fNumberOfEvents(qas.fNumberOfEvents),
150 fRunNumber(qas.fRunNumber),
152 fRawReaderDelete(kTRUE),
157 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
158 fLoader[iDet] = qas.fLoader[iDet] ;
159 fQADataMaker[iDet] = qas.fQADataMaker[iDet] ;
160 fQACycles[iDet] = qas.fQACycles[iDet] ;
161 fQAWriteExpert[iDet] = qas.fQAWriteExpert[iDet] ;
165 //_____________________________________________________________________________
166 AliQAManager & AliQAManager::operator = (const AliQAManager & qas)
168 // assignment operator
169 this->~AliQAManager() ;
170 new(this) AliQAManager(qas) ;
174 //_____________________________________________________________________________
175 AliQAManager::~AliQAManager()
178 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
179 if (IsSelected(AliQAv1::GetDetName(iDet))) {
180 fLoader[iDet] = NULL;
181 if (fQADataMaker[iDet]) {
182 (fQADataMaker[iDet])->Finish() ;
183 delete fQADataMaker[iDet] ;
187 if (fRawReaderDelete) {
194 //_____________________________________________________________________________
195 Bool_t AliQAManager::DoIt(const AliQAv1::TASKINDEX_t taskIndex)
197 // Runs all the QA data Maker for every detector
200 // Fill QA data in event loop
201 for (UInt_t iEvent = fFirstEvent ; iEvent < (UInt_t)fMaxEvents ; iEvent++) {
204 if ( iEvent%10 == 0 )
205 AliDebug(AliQAv1::GetQADebugLevel(), Form("processing event %d", iEvent));
206 if ( taskIndex == AliQAv1::kRAWS ) {
207 if ( !fRawReader->NextEvent() )
209 } else if ( taskIndex == AliQAv1::kESDS ) {
210 if ( fESDTree->GetEntry(iEvent) == 0 )
213 if ( fRunLoader->GetEvent(iEvent) != 0 )
216 // loop over active loaders
217 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
218 if (IsSelected(AliQAv1::GetDetName(iDet))) {
219 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
220 if (!qadm) continue; // This detector doesn't have any QA (for example, HLT)
221 if ( qadm->IsCycleDone() ) {
222 qadm->EndOfCycle(taskIndex) ;
224 TTree * data = NULL ;
225 AliLoader* loader = GetLoader(qadm->GetUniqueID());
227 case AliQAv1::kNULLTASKINDEX :
229 case AliQAv1::kRAWS :
230 qadm->Exec(taskIndex, fRawReader) ;
232 case AliQAv1::kHITS :
235 data = loader->TreeH() ;
237 AliWarning(Form(" Hit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
241 qadm->Exec(taskIndex, data) ;
243 case AliQAv1::kSDIGITS :
245 loader->LoadSDigits() ;
246 data = loader->TreeS() ;
248 AliWarning(Form(" SDigit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
252 qadm->Exec(taskIndex, data) ;
254 case AliQAv1::kDIGITS :
256 loader->LoadDigits() ;
257 data = loader->TreeD() ;
259 AliWarning(Form(" Digit Tree not found for %s", AliQAv1::GetDetName(iDet))) ;
263 qadm->Exec(taskIndex, data) ;
265 case AliQAv1::kRECPOINTS :
267 loader->LoadRecPoints() ;
268 data = loader->TreeR() ;
270 AliWarning(Form("RecPoints not found for %s", AliQAv1::GetDetName(iDet))) ;
274 qadm->Exec(taskIndex, data) ;
276 case AliQAv1::kTRACKSEGMENTS :
278 case AliQAv1::kRECPARTICLES :
280 case AliQAv1::kESDS :
281 qadm->Exec(taskIndex, fESD) ;
283 case AliQAv1::kNTASKINDEX :
290 // Save QA data for all detectors
291 rv = Finish(taskIndex) ;
293 if ( taskIndex == AliQAv1::kRAWS )
294 fRawReader->RewindEvents() ;
299 //_____________________________________________________________________________
300 Bool_t AliQAManager::Finish(const AliQAv1::TASKINDEX_t taskIndex)
302 // write output to file for all detectors
303 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
304 if (IsSelected(AliQAv1::GetDetName(iDet))) {
305 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
307 qadm->EndOfCycle(taskIndex) ;
313 //_____________________________________________________________________________
314 TObjArray * AliQAManager::GetFromOCDB(AliQAv1::DETECTORINDEX_t det, AliQAv1::TASKINDEX_t task, const char * year) const
316 // Retrieve the list of QA data for a given detector and a given task
317 TObjArray * rv = NULL ;
318 if ( !strlen(AliQAv1::GetQARefStorage()) ) {
319 AliError("No storage defined, use AliQAv1::SetQARefStorage") ;
322 if ( ! IsDefaultStorageSet() ) {
323 TString tmp(AliQAv1::GetQARefDefaultStorage()) ;
326 Instance()->SetDefaultStorage(tmp.Data()) ;
327 Instance()->SetSpecificStorage(Form("%s/*", AliQAv1::GetQAName()), AliQAv1::GetQARefStorage()) ;
329 TString detOCDBDir(Form("%s/%s/%s", AliQAv1::GetQAName(), AliQAv1::GetDetName((Int_t)det), AliQAv1::GetRefOCDBDirName())) ;
330 AliDebug(AliQAv1::GetQADebugLevel(), Form("Retrieving reference data from %s/%s for %s", AliQAv1::GetQARefStorage(), detOCDBDir.Data(), AliQAv1::GetTaskName(task).Data())) ;
331 AliCDBEntry* entry = QAManager()->Get(detOCDBDir.Data(), 0) ; //FIXME 0 --> Run Number
332 TList * listDetQAD = dynamic_cast<TList *>(entry->GetObject()) ;
334 rv = dynamic_cast<TObjArray *>(listDetQAD->FindObject(AliQAv1::GetTaskName(task))) ;
338 //_____________________________________________________________________________
339 AliLoader * AliQAManager::GetLoader(Int_t iDet)
341 // get the loader for a detector
343 if ( !fRunLoader || iDet == AliQAv1::kCORR)
346 TString detName = AliQAv1::GetDetName(iDet) ;
347 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
349 return fLoader[iDet] ;
351 // load the QA data maker object
352 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
353 TString loaderName = "Ali" + detName + "Loader" ;
355 AliLoader * loader = NULL ;
356 // first check if a plugin is defined for the quality assurance data maker
357 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
358 // if not, add a plugin for it
359 if (!pluginHandler) {
360 AliDebug(AliQAv1::GetQADebugLevel(), Form("defining plugin for %s", loaderName.Data())) ;
361 TString libs = gSystem->GetLibraries() ;
362 if (libs.Contains("lib" + detName + "base.so") || (gSystem->Load("lib" + detName + "base.so") >= 0)) {
363 pluginManager->AddHandler("AliQADataMaker", detName, loaderName, detName + "loader", loaderName + "()") ;
365 pluginManager->AddHandler("AliLoader", detName, loaderName, detName, loaderName + "()") ;
367 pluginHandler = pluginManager->FindHandler("AliLoader", detName) ;
369 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
370 loader = (AliLoader *) pluginHandler->ExecPlugin(0) ;
373 fLoader[iDet] = loader ;
377 //_____________________________________________________________________________
378 AliQAv1 * AliQAManager::GetQA(UInt_t run, UInt_t evt)
380 // retrieves the QA object stored in a file named "Run{run}.Event{evt}_1.ESD.tag.root"
381 char * fileName = Form("Run%d.Event%d_1.ESD.tag.root", run, evt) ;
382 TFile * tagFile = TFile::Open(fileName) ;
384 AliError(Form("File %s not found", fileName)) ;
387 TTree * tagTree = dynamic_cast<TTree *>(tagFile->Get("T")) ;
389 AliError(Form("Tree T not found in %s", fileName)) ;
393 AliRunTag * tag = new AliRunTag ;
394 tagTree->SetBranchAddress("AliTAG", &tag) ;
395 tagTree->GetEntry(evt) ;
396 AliQAv1 * qa = AliQAv1::Instance(tag->GetQALength(), tag->GetQAArray(), tag->GetESLength(), tag->GetEventSpecies()) ;
401 //_____________________________________________________________________________
402 AliQADataMaker * AliQAManager::GetQADataMaker(const Int_t iDet)
404 // get the quality assurance data maker for a detector
406 if (fQADataMaker[iDet]) {
407 fQADataMaker[iDet]->SetEventSpecie(fQADataMaker[iDet]->GetRecoParam()->GetEventSpecie()) ;
408 return fQADataMaker[iDet] ;
411 AliQADataMaker * qadm = NULL ;
413 if (iDet == AliQAv1::kGLOBAL) { //Global QA
414 qadm = new AliGlobalQADataMaker();
415 qadm->SetName(AliQAv1::GetDetName(iDet));
416 qadm->SetUniqueID(iDet);
417 fQADataMaker[iDet] = qadm;
418 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
422 if (iDet == AliQAv1::kCORR) { //the data maker for correlations among detectors
423 qadm = new AliCorrQADataMakerRec(fQADataMaker) ;
424 qadm->SetName(AliQAv1::GetDetName(iDet));
425 qadm->SetUniqueID(iDet);
426 fQADataMaker[iDet] = qadm;
427 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
431 // load the QA data maker object
432 TPluginManager* pluginManager = gROOT->GetPluginManager() ;
433 TString detName = AliQAv1::GetDetName(iDet) ;
435 if (tmp.Contains("sim"))
436 tmp.ReplaceAll("s", "S") ;
437 else if (tmp.Contains("rec"))
438 tmp.ReplaceAll("r", "R") ;
440 TString qadmName = "Ali" + detName + "QADataMaker" + tmp ;
442 // first check if a plugin is defined for the quality assurance data maker
443 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
444 // if not, add a plugin for it
445 if (!pluginHandler) {
446 AliDebug(AliQAv1::GetQADebugLevel(), Form("defining plugin for %s", qadmName.Data())) ;
447 TString libs = gSystem->GetLibraries() ;
448 if (libs.Contains("lib" + detName + fMode + ".so") || (gSystem->Load("lib" + detName + fMode + ".so") >= 0)) {
449 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName + "qadm", qadmName + "()") ;
451 pluginManager->AddHandler("AliQADataMaker", detName, qadmName, detName, qadmName + "()") ;
453 pluginHandler = pluginManager->FindHandler("AliQADataMaker", detName) ;
455 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
456 qadm = (AliQADataMaker *) pluginHandler->ExecPlugin(0) ;
459 qadm->SetName(AliQAv1::GetDetName(iDet));
460 qadm->SetUniqueID(iDet);
461 fQADataMaker[iDet] = qadm ;
462 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
468 //_____________________________________________________________________________
469 void AliQAManager::EndOfCycle(TObjArray * detArray)
471 // End of cycle QADataMakers
473 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
474 if (IsSelected(AliQAv1::GetDetName(iDet))) {
475 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
478 // skip non active detectors
480 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQAv1::GetDetName(iDet))) ;
481 if (!det || !det->IsActive())
484 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
485 if ( fTasks.Contains(Form("%d", taskIndex)) )
486 qadm->EndOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex))) ;
493 //_____________________________________________________________________________
494 void AliQAManager::EndOfCycle(TString detectors)
496 // End of cycle QADataMakers
498 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
499 if (IsSelected(AliQAv1::GetDetName(iDet))) {
500 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
503 // skip non active detectors
504 if (!detectors.Contains(AliQAv1::GetDetName(iDet)))
506 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
507 if ( fTasks.Contains(Form("%d", taskIndex)) )
508 qadm->EndOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex))) ;
515 //_____________________________________________________________________________
516 void AliQAManager::Increment()
518 // Increments the cycle counter for all QA Data Makers
519 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
520 if (IsSelected(AliQAv1::GetDetName(iDet))) {
521 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
528 //_____________________________________________________________________________
529 Bool_t AliQAManager::InitQA(const AliQAv1::TASKINDEX_t taskIndex, const char * input )
531 // Initialize the event source and QA data makers
533 fTasks += Form("%d", taskIndex) ;
535 if (taskIndex == AliQAv1::kRAWS) {
537 fRawReader = AliRawReader::Create(input);
541 fRawReaderDelete = kTRUE ;
542 fRawReader->NextEvent() ;
543 fRunNumber = fRawReader->GetRunNumber() ;
545 fRawReader->RewindEvents();
546 fNumberOfEvents = 999999 ;
547 if ( fMaxEvents < 0 )
548 fMaxEvents = fNumberOfEvents ;
549 } else if (taskIndex == AliQAv1::kESDS) {
550 fTasks = AliQAv1::GetTaskName(AliQAv1::kESDS) ;
551 if (!gSystem->AccessPathName("AliESDs.root")) { // AliESDs.root exists
552 TFile * esdFile = TFile::Open("AliESDs.root") ;
553 fESDTree = dynamic_cast<TTree *> (esdFile->Get("esdTree")) ;
555 AliError("esdTree not found") ;
558 fESD = new AliESDEvent() ;
559 fESD->ReadFromTree(fESDTree) ;
560 fESDTree->GetEntry(0) ;
561 fRunNumber = fESD->GetRunNumber() ;
562 fNumberOfEvents = fESDTree->GetEntries() ;
563 if ( fMaxEvents < 0 )
564 fMaxEvents = fNumberOfEvents ;
567 AliError("AliESDs.root not found") ;
571 if ( !InitRunLoader() ) {
572 AliWarning("No Run Loader not found") ;
574 fNumberOfEvents = fRunLoader->GetNumberOfEvents() ;
575 if ( fMaxEvents < 0 )
576 fMaxEvents = fNumberOfEvents ;
581 TObjArray* detArray = NULL ;
582 if (fRunLoader) // check if RunLoader exists
583 if ( fRunLoader->GetAliRun() ) { // check if AliRun exists in gAlice.root
584 detArray = fRunLoader->GetAliRun()->Detectors() ;
585 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
588 // Initialize all QA data makers for all detectors
589 fRunNumber = AliCDBManager::Instance()->GetRun() ;
590 if ( ! AliGeomManager::GetGeometry() )
591 AliGeomManager::LoadGeometry() ;
593 InitQADataMaker(fRunNumber, detArray) ; //, fCycleSame, kTRUE, detArray) ;
597 //_____________________________________________________________________________
598 void AliQAManager::InitQADataMaker(UInt_t run, TObjArray * detArray)
600 // Initializes The QADataMaker for all active detectors and for all active tasks
601 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
602 if (IsSelected(AliQAv1::GetDetName(iDet))) {
603 AliQADataMaker * qadm = GetQADataMaker(iDet) ;
605 AliError(Form("AliQADataMaker not found for %s", AliQAv1::GetDetName(iDet))) ;
606 fDetectorsW.ReplaceAll(AliQAv1::GetDetName(iDet), "") ;
608 if (fQAWriteExpert[iDet])
609 qadm->SetWriteExpert() ;
610 AliDebug(AliQAv1::GetQADebugLevel(), Form("Data Maker found for %s %d", qadm->GetName(), qadm->WriteExpert())) ;
611 // skip non active detectors
613 AliModule* det = static_cast<AliModule*>(detArray->FindObject(AliQAv1::GetDetName(iDet))) ;
614 if (!det || !det->IsActive())
617 // Set default reco params
618 Bool_t sameCycle = kFALSE ;
619 for (UInt_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
620 if ( fTasks.Contains(Form("%d", taskIndex)) ) {
621 qadm->Init(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex)), GetQACycles(qadm->GetUniqueID())) ;
622 qadm->StartOfCycle(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(taskIndex)), run, sameCycle) ;
632 //_____________________________________________________________________________
633 Bool_t AliQAManager::InitRunLoader()
635 // get or create the run loader
639 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
640 // load all base libraries to get the loader classes
641 TString libs = gSystem->GetLibraries() ;
642 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
643 if (!IsSelected(AliQAv1::GetDetName(iDet)))
645 TString detName = AliQAv1::GetDetName(iDet) ;
646 if (detName == "HLT")
648 if (libs.Contains("lib" + detName + "base.so"))
650 gSystem->Load("lib" + detName + "base.so");
652 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
654 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
657 fRunLoader->CdGAFile();
658 if (fRunLoader->LoadgAlice() == 0) {
659 gAlice = fRunLoader->GetAliRun();
663 AliError(Form("no gAlice object found in file %s", fGAliceFileName.Data()));
667 } else { // galice.root does not exist
668 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
674 fRunLoader->LoadHeader();
675 fRunNumber = fRunLoader->GetHeader()->GetRun() ;
680 //_____________________________________________________________________________
681 Bool_t AliQAManager::IsSelected(const char * det)
683 // check whether detName is contained in detectors
684 // if yes, it is removed from detectors
687 const TString detName(det) ;
688 // always activates Correlation
689 if ( detName.Contains(AliQAv1::GetDetName(AliQAv1::kCORR))) {
692 // check if all detectors are selected
693 if (fDetectors.Contains("ALL")) {
696 } else if ((fDetectors.CompareTo(detName) == 0) ||
697 fDetectors.BeginsWith(detName+" ") ||
698 fDetectors.EndsWith(" "+detName) ||
699 fDetectors.Contains(" "+detName+" ")) {
706 //_____________________________________________________________________________
707 Bool_t AliQAManager::Merge(Int_t runNumber) const
709 // Merge data from all detectors from a given run in one single file
710 // Merge the QA results from all the data chunks in one run
711 if ( runNumber == -1)
712 runNumber = fRunNumber ;
713 Bool_t rv = MergeData(runNumber) ;
714 //rv *= MergeResults(runNumber) ; // not needed for the time being
718 //______________________________________________________________________
719 Bool_t AliQAManager::MergeXML(const char * collectionFile, const char * subFile, const char * outFile)
721 // merges files listed in a xml collection
722 // usage Merge(collection, outputFile))
723 // collection: is a xml collection
727 if ( strstr(collectionFile, ".xml") == 0 ) {
728 AliError("Input collection file must be an \".xml\" file\n") ;
733 TGrid::Connect("alien://");
737 // Open the file collection
738 printf("*** Create Collection ***\n");
739 printf("*** Wk-Dir = |%s| \n",gSystem->WorkingDirectory());
740 printf("*** Coll = |%s| \n",collectionFile);
742 TGridCollection * collection = (TGridCollection*)gROOT->ProcessLine(Form("TAlienCollection::Open(\"%s\")",collectionFile));
743 TGridResult* result = collection->GetGridResult("", 0, 0);
747 TFileMerger merger(kFALSE) ;
749 TString tempo(collectionFile) ;
751 tempo.ReplaceAll(".xml", subFile) ;
753 tempo.ReplaceAll(".xml", "_Merged.root") ;
754 outFile = tempo.Data() ;
756 merger.OutputFile(outFile) ;
758 while ( (turl = result->GetKey(index, "turl")) ) {
761 file = Form("%s#%s", turl, subFile) ;
763 file = Form("%s", turl) ;
765 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s\n", file)) ;
766 merger.AddFile(file) ;
773 AliDebug(AliQAv1::GetQADebugLevel(), Form("Files merged into %s\n", outFile)) ;
779 //_____________________________________________________________________________
780 void AliQAManager::MergeCustom() const
782 // Custom Merge of QA data from all detectors for all runs in one single file
783 // search all the run numbers
784 // search all the run numbers
785 gROOT->ProcessLine(".! ls *QA*.root > QAtempo.txt") ;
787 FILE * QAfiles = fopen("QAtempo.txt", "r") ;
790 TIter nextRun(&srunList) ;
791 TObjString * srun = NULL ;
792 Int_t loRun = 999999999 ;
794 while ( QAfile.Gets(QAfiles) ) {
795 Bool_t runExist = kFALSE ;
796 TString srunNew(QAfile(QAfile.Index("QA.")+3, QAfile.Index(".root")-(QAfile.Index("QA.")+3))) ;
797 Int_t cuRun = srunNew.Atoi() ;
802 while ( (srun = dynamic_cast<TObjString *> (nextRun())) ) {
803 if ( cuRun == (srun->String()).Atoi() ) {
810 srunList.Add(new TObjString(srunNew.Data()));
813 Int_t runNumber = 0 ;
814 TFile mergedFile(Form("Merged.%s.Data.root", AliQAv1::GetQADataFileName()), "RECREATE") ;
815 TH1I * hisRun = new TH1I("hLMR", "List of merged runs", hiRun-loRun+10, loRun, hiRun+10) ;
816 // create the structure into the merged file
817 for (Int_t iDet = 0; iDet < AliQAv1::kNDET ; iDet++) {
818 TDirectory * detDir = mergedFile.mkdir(AliQAv1::GetDetName(iDet)) ;
819 for (Int_t taskIndex = 0; taskIndex < AliQAv1::kNTASKINDEX; taskIndex++) {
821 TDirectory * taskDir = gDirectory->mkdir(AliQAv1::GetTaskName(taskIndex)) ;
822 for (Int_t es = 0 ; es < AliRecoParam::kNSpecies ; es++) {
824 TDirectory * esDir = gDirectory->mkdir(AliRecoParam::GetEventSpecieName(es)) ;
826 gDirectory->mkdir(AliQAv1::GetExpert()) ;
830 while ( (srun = dynamic_cast<TObjString *> (nextRun())) ) {
831 runNumber = (srun->String()).Atoi() ;
832 hisRun->Fill(runNumber) ;
833 AliDebug(AliQAv1::GetQADebugLevel(), Form("Merging run number %d", runNumber)) ;
834 // search all QA files for runNumber in the current directory
835 char * fileList[AliQAv1::kNDET] ;
837 for (Int_t iDet = 0; iDet < AliQAv1::kNDET ; iDet++) {
838 char * file = gSystem->Which(gSystem->WorkingDirectory(), Form("%s.%s.%d.root", AliQAv1::GetDetName(iDet), AliQAv1::GetQADataFileName(), runNumber));
840 fileList[index++] = file ;
843 AliError("No QA data file found\n") ;
846 for ( Int_t i = 0 ; i < index ; i++) {
847 TFile * inFile = TFile::Open(fileList[i]) ;
848 TList * listOfKeys =inFile->GetListOfKeys() ;
849 TIter nextkey(listOfKeys) ;
851 TString dirName("") ;
852 while ( (obj1 = nextkey()) ) {
853 TDirectory * directoryDet = inFile->GetDirectory(obj1->GetName()) ;
854 if ( directoryDet ) {
855 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s dir = %s", inFile->GetName(), directoryDet->GetName())) ;
856 dirName += Form("%s/", directoryDet->GetName() ) ;
858 TList * listOfTasks = directoryDet->GetListOfKeys() ;
859 TIter nextTask(listOfTasks) ;
861 while ( (obj2 = nextTask()) ) {
862 TDirectory * directoryTask = directoryDet->GetDirectory(obj2->GetName()) ;
863 if ( directoryTask ) {
864 dirName += Form("%s", obj2->GetName()) ;
865 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s", dirName.Data())) ;
866 directoryTask->cd() ;
867 TList * listOfEventSpecie = directoryTask->GetListOfKeys() ;
868 TIter nextEventSpecie(listOfEventSpecie) ;
870 while ( (obj3 = nextEventSpecie()) ) {
871 TDirectory * directoryEventSpecie = directoryTask->GetDirectory(obj3->GetName()) ;
872 if ( directoryEventSpecie ) {
873 dirName += Form("/%s/", obj3->GetName()) ;
874 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s\n", dirName.Data())) ;
875 directoryEventSpecie->cd() ;
876 // histograms are here
877 TDirectory * mergedDirectory = mergedFile.GetDirectory(dirName.Data()) ;
878 TList * listOfData = directoryEventSpecie->GetListOfKeys() ;
879 TIter nextData(listOfData) ;
881 while ( (key = dynamic_cast<TKey *>(nextData())) ) {
882 TString className(key->GetClassName()) ;
883 if ( className.Contains("TH") || className.Contains("TProfile") ) {
884 TH1 * histIn = dynamic_cast<TH1*> (key->ReadObj()) ;
885 TH1 * histOu = dynamic_cast<TH1*> (mergedDirectory->FindObjectAny(histIn->GetName())) ;
886 AliDebug(AliQAv1::GetQADebugLevel(), Form("%s %x %x\n", key->GetName(), histIn, histOu)) ;
887 mergedDirectory->cd() ;
891 histOu->Add(histIn) ;
892 histOu->Write(histOu->GetName(), kOverwrite) ;
895 else if ( className.Contains("TDirectoryFile") ) {
896 TDirectory * dirExpert = directoryEventSpecie->GetDirectory(key->GetName()) ;
898 TDirectory * mergedDirectoryExpert = mergedDirectory->GetDirectory(dirExpert->GetName()) ;
899 TList * listOfExpertData = dirExpert->GetListOfKeys() ;
900 TIter nextExpertData(listOfExpertData) ;
902 while ( (keykey = dynamic_cast<TKey *>(nextExpertData())) ) {
903 TString classNameExpert(keykey->GetClassName()) ;
904 if (classNameExpert.Contains("TH")) {
905 TH1 * histInExpert = dynamic_cast<TH1*> (keykey->ReadObj()) ;
906 TH1 * histOuExpert = dynamic_cast<TH1*> (mergedDirectory->FindObjectAny(histInExpert->GetName())) ;
907 mergedDirectoryExpert->cd() ;
908 if ( ! histOuExpert ) {
909 histInExpert->Write() ;
911 histOuExpert->Add(histInExpert) ;
912 histOuExpert->Write(histOuExpert->GetName(), kOverwrite) ;
917 AliError(Form("No merge done for this object %s in %s", key->GetName(), dirName.Data())) ;
920 dirName.ReplaceAll(Form("/%s/",obj3->GetName()), "") ;
923 dirName.ReplaceAll(obj2->GetName(), "") ;
937 //_____________________________________________________________________________
938 Bool_t AliQAManager::MergeData(const Int_t runNumber) const
940 // Merge QA data from all detectors for a given run in one single file
943 TString outFileName = Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()) ;
944 merger.OutputFile(outFileName.Data()) ;
945 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
946 char * file = gSystem->Which(gSystem->WorkingDirectory(), Form("%s.%s.%d.root", AliQAv1::GetDetName(iDet), AliQAv1::GetQADataFileName(), runNumber));
948 merger.AddFile(file) ;
954 //_____________________________________________________________________________
955 Bool_t AliQAManager::MergeResults(const Int_t runNumber) const
957 // Merge the QA result from all the data chunks in a run
958 // to be revised whwn it will be used (see MergeData)
960 cmd = Form(".! ls %s*.root > tempo.txt", AliQAv1::GetQADataFileName()) ;
961 gROOT->ProcessLine(cmd.Data()) ;
962 ifstream in("tempo.txt") ;
963 const Int_t chunkMax = 100 ;
964 TString fileList[chunkMax] ;
969 in >> fileList[index] ;
972 AliDebug(AliQAv1::GetQADebugLevel(), Form("index = %d file = %s", index, (fileList[index].Data()))) ;
977 AliError("No QA Result File found") ;
982 TString outFileName ;
984 outFileName = Form("Merged.%s.Result.%d.root",AliQAv1::GetQADataFileName(),runNumber);
986 outFileName = Form("Merged.%s.Result.root",AliQAv1::GetQADataFileName());
987 merger.OutputFile(outFileName.Data()) ;
988 for (Int_t ifile = 0 ; ifile < index ; ifile++) {
989 TString file = fileList[ifile] ;
990 merger.AddFile(file) ;
997 //_____________________________________________________________________________
998 void AliQAManager::Reset(const Bool_t sameCycle)
1000 // Reset the default data members
1002 for (UInt_t iDet = 0; iDet < fgkNDetectors ; iDet++) {
1003 if (IsSelected(AliQAv1::GetDetName(iDet))) {
1004 AliQADataMaker * qadm = GetQADataMaker(iDet);
1008 if (fRawReaderDelete) {
1013 fCycleSame = sameCycle ;
1017 fNumberOfEvents = 999999 ;
1020 //_____________________________________________________________________________
1021 AliQAManager * AliQAManager::QAManager(const Char_t * mode, TMap *entryCache, Int_t run)
1023 // returns AliQAManager instance (singleton)
1025 if (!fgQAInstance) {
1026 fgQAInstance = new AliQAManager(mode) ;
1028 fgQAInstance->Init();
1030 fgQAInstance->InitFromCache(entryCache,run);
1032 return fgQAInstance;
1035 //_____________________________________________________________________________
1036 TString AliQAManager::Run(const char * detectors, AliRawReader * rawReader, const Bool_t sameCycle)
1038 //Runs all the QA data Maker for Raws only
1040 fCycleSame = sameCycle ;
1041 fRawReader = rawReader ;
1042 fDetectors = detectors ;
1043 fDetectorsW = detectors ;
1045 AliCDBManager* man = AliCDBManager::Instance() ;
1047 if ( man->GetRun() == -1 ) {// check if run number not set previously and set it from raw data
1048 rawReader->NextEvent() ;
1049 man->SetRun(fRawReader->GetRunNumber()) ;
1050 rawReader->RewindEvents() ;
1054 if ( !InitQA(AliQAv1::kRAWS) )
1056 fRawReaderDelete = kFALSE ;
1058 DoIt(AliQAv1::kRAWS) ;
1059 return fDetectorsW ;
1062 //_____________________________________________________________________________
1063 TString AliQAManager::Run(const char * detectors, const char * fileName, const Bool_t sameCycle)
1065 //Runs all the QA data Maker for Raws only
1067 fCycleSame = sameCycle ;
1068 fDetectors = detectors ;
1069 fDetectorsW = detectors ;
1071 AliCDBManager* man = AliCDBManager::Instance() ;
1072 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
1073 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
1075 AliFatal("galice.root file not found in current directory") ;
1079 if ( ! rl->GetAliRun() ) {
1080 AliFatal("AliRun not found in galice.root") ;
1083 man->SetRun(rl->GetHeader()->GetRun());
1089 if ( !InitQA(AliQAv1::kRAWS, fileName) )
1092 DoIt(AliQAv1::kRAWS) ;
1093 return fDetectorsW ;
1096 //_____________________________________________________________________________
1097 TString AliQAManager::Run(const char * detectors, const AliQAv1::TASKINDEX_t taskIndex, Bool_t const sameCycle, const char * fileName )
1099 // Runs all the QA data Maker for every detector
1101 fCycleSame = sameCycle ;
1102 fDetectors = detectors ;
1103 fDetectorsW = detectors ;
1105 AliCDBManager* man = AliCDBManager::Instance() ;
1106 if ( man->GetRun() == -1 ) { // check if run number not set previously and set it from AliRun
1107 AliRunLoader * rl = AliRunLoader::Open("galice.root") ;
1109 AliFatal("galice.root file not found in current directory") ;
1113 if ( ! rl->GetAliRun() ) {
1114 AliDebug(AliQAv1::GetQADebugLevel(), "AliRun not found in galice.root") ;
1117 man->SetRun(rl->GetHeader()->GetRun()) ;
1123 if ( taskIndex == AliQAv1::kNULLTASKINDEX) {
1124 for (UInt_t task = 0; task < AliQAv1::kNTASKINDEX; task++) {
1125 if ( fTasks.Contains(Form("%d", task)) ) {
1127 if ( !InitQA(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(task)), fileName) )
1129 DoIt(AliQAv1::GetTaskIndex(AliQAv1::GetTaskName(task))) ;
1134 if ( !InitQA(taskIndex, fileName) )
1139 return fDetectorsW ;
1143 //_____________________________________________________________________________
1144 void AliQAManager::RunOneEvent(AliRawReader * rawReader)
1146 //Runs all the QA data Maker for Raws only and on one event only (event loop done by calling method)
1149 AliCodeTimerAuto("") ;
1150 if (fTasks.Contains(Form("%d", AliQAv1::kRAWS))){
1151 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1152 if (!IsSelected(AliQAv1::GetDetName(iDet)))
1154 AliQADataMaker *qadm = GetQADataMaker(iDet);
1157 if ( qadm->IsCycleDone() ) {
1158 qadm->EndOfCycle() ;
1160 AliCodeTimerStart(Form("running RAW quality assurance data maker for %s", AliQAv1::GetDetName(iDet)));
1161 qadm->SetEventSpecie(qadm->GetRecoParam()->GetEventSpecie()) ;
1162 qadm->Exec(AliQAv1::kRAWS, rawReader) ;
1163 AliCodeTimerStop(Form("running RAW quality assurance data maker for %s", AliQAv1::GetDetName(iDet)));
1168 //_____________________________________________________________________________
1169 void AliQAManager::RunOneEvent(AliESDEvent *& esd)
1171 //Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1173 AliCodeTimerAuto("") ;
1174 if (fTasks.Contains(Form("%d", AliQAv1::kESDS))) {
1175 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1176 if (!IsSelected(AliQAv1::GetDetName(iDet)))
1178 AliQADataMaker *qadm = GetQADataMaker(iDet);
1181 if ( qadm->IsCycleDone() ) {
1182 qadm->EndOfCycle() ;
1184 AliCodeTimerStart(Form("running ESD quality assurance data maker for %s", AliQAv1::GetDetName(iDet)));
1185 qadm->Exec(AliQAv1::kESDS, esd) ;
1186 AliCodeTimerStop(Form("running ESD quality assurance data maker for %s", AliQAv1::GetDetName(iDet)));
1191 //_____________________________________________________________________________
1192 void AliQAManager::RunOneEventInOneDetector(Int_t det, TTree * tree)
1194 // Runs all the QA data Maker for ESDs only and on one event only (event loop done by calling method)
1195 AliCodeTimerAuto("") ;
1196 if (fTasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
1197 if (IsSelected(AliQAv1::GetDetName(det))) {
1198 AliQADataMaker *qadm = GetQADataMaker(det);
1200 if ( qadm->IsCycleDone() ) {
1201 qadm->EndOfCycle() ;
1203 AliCodeTimerStart(Form("running RecPoints quality assurance data maker for %s", AliQAv1::GetDetName(det)));
1204 qadm->Exec(AliQAv1::kRECPOINTS, tree) ;
1205 AliCodeTimerStop(Form("running RecPoints quality assurance data maker for %s", AliQAv1::GetDetName(det)));
1211 //_____________________________________________________________________________
1212 Bool_t AliQAManager::Save2OCDB(const Int_t runNumber, AliRecoParam::EventSpecie_t es, const char * year, const char * detectors) const
1214 // take the locasl QA data merge into a single file and save in OCDB
1216 TString tmp(AliQAv1::GetQARefStorage()) ;
1217 if ( tmp.IsNull() ) {
1218 AliError("No storage defined, use AliQAv1::SetQARefStorage") ;
1221 if ( !(tmp.Contains(AliQAv1::GetLabLocalOCDB()) || tmp.Contains(AliQAv1::GetLabAliEnOCDB())) ) {
1222 AliError(Form("%s is a wrong storage, use %s or %s", AliQAv1::GetQARefStorage(), AliQAv1::GetLabLocalOCDB().Data(), AliQAv1::GetLabAliEnOCDB().Data())) ;
1225 TString sdet(detectors) ;
1228 if ( sdet.Contains("ALL") ) {
1229 rv = Merge(runNumber) ;
1232 TString inputFileName(Form("Merged.%s.Data.%d.root", AliQAv1::GetQADataFileName(), runNumber)) ;
1233 inputFile = TFile::Open(inputFileName.Data()) ;
1234 rv = SaveIt2OCDB(runNumber, inputFile, year, es) ;
1236 for (Int_t index = 0; index < AliQAv1::kNDET; index++) {
1237 if (sdet.Contains(AliQAv1::GetDetName(index))) {
1238 TString inputFileName(Form("%s.%s.%d.root", AliQAv1::GetDetName(index), AliQAv1::GetQADataFileName(), runNumber)) ;
1239 inputFile = TFile::Open(inputFileName.Data()) ;
1240 rv *= SaveIt2OCDB(runNumber, inputFile, year, es) ;
1247 //_____________________________________________________________________________
1248 Bool_t AliQAManager::SaveIt2OCDB(const Int_t runNumber, TFile * inputFile, const char * year, AliRecoParam::EventSpecie_t es) const
1250 // reads the TH1 from file and adds it to appropriate list before saving to OCDB
1252 AliDebug(AliQAv1::GetQADebugLevel(), Form("Saving TH1s in %s to %s", inputFile->GetName(), AliQAv1::GetQARefStorage())) ;
1253 if ( ! IsDefaultStorageSet() ) {
1254 TString tmp( AliQAv1::GetQARefStorage() ) ;
1255 if ( tmp.Contains(AliQAv1::GetLabLocalOCDB()) )
1256 Instance()->SetDefaultStorage(AliQAv1::GetQARefStorage()) ;
1258 TString tmp1(AliQAv1::GetQARefDefaultStorage()) ;
1260 tmp1.Append("?user=alidaq") ;
1261 Instance()->SetDefaultStorage(tmp1.Data()) ;
1264 Instance()->SetSpecificStorage("*", AliQAv1::GetQARefStorage()) ;
1266 Instance()->SetRun(runNumber);
1268 AliCDBMetaData mdr ;
1269 mdr.SetResponsible("yves schutz");
1271 for ( Int_t detIndex = 0 ; detIndex < AliQAv1::kNDET ; detIndex++) {
1272 TDirectory * detDir = inputFile->GetDirectory(AliQAv1::GetDetName(detIndex)) ;
1274 AliDebug(AliQAv1::GetQADebugLevel(), Form("Entering %s", detDir->GetName())) ;
1275 AliQAv1::SetQARefDataDirName(es) ;
1276 TString detOCDBDir(Form("%s/%s/%s", AliQAv1::GetDetName(detIndex), AliQAv1::GetRefOCDBDirName(), AliQAv1::GetRefDataDirName())) ;
1277 AliCDBId idr(detOCDBDir.Data(), runNumber, AliCDBRunRange::Infinity()) ;
1278 TList * listDetQAD = new TList() ;
1279 TString listName(Form("%s QA data Reference", AliQAv1::GetDetName(detIndex))) ;
1280 mdr.SetComment(Form("%s QA stuff", AliQAv1::GetDetName(detIndex)));
1281 listDetQAD->SetName(listName) ;
1282 TList * taskList = detDir->GetListOfKeys() ;
1283 TIter nextTask(taskList) ;
1285 while ( (taskKey = dynamic_cast<TKey*>(nextTask())) ) {
1286 TDirectory * taskDir = detDir->GetDirectory(taskKey->GetName()) ;
1287 TDirectory * esDir = taskDir->GetDirectory(AliRecoParam::GetEventSpecieName(es)) ;
1288 AliDebug(AliQAv1::GetQADebugLevel(), Form("Saving %s", esDir->GetName())) ;
1289 TObjArray * listTaskQAD = new TObjArray(100) ;
1290 listTaskQAD->SetName(Form("%s/%s", taskKey->GetName(), AliRecoParam::GetEventSpecieName(es))) ;
1291 listDetQAD->Add(listTaskQAD) ;
1292 TList * histList = esDir->GetListOfKeys() ;
1293 TIter nextHist(histList) ;
1295 while ( (histKey = dynamic_cast<TKey*>(nextHist())) ) {
1296 TObject * odata = esDir->Get(histKey->GetName()) ;
1298 AliError(Form("%s in %s/%s returns a NULL pointer !!", histKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1300 if ( AliQAv1::GetExpert() == histKey->GetName() ) {
1301 TDirectory * expertDir = esDir->GetDirectory(histKey->GetName()) ;
1302 TList * expertHistList = expertDir->GetListOfKeys() ;
1303 TIter nextExpertHist(expertHistList) ;
1304 TKey * expertHistKey ;
1305 while ( (expertHistKey = dynamic_cast<TKey*>(nextExpertHist())) ) {
1306 TObject * expertOdata = expertDir->Get(expertHistKey->GetName()) ;
1307 if ( !expertOdata ) {
1308 AliError(Form("%s in %s/%s/Expert returns a NULL pointer !!", expertHistKey->GetName(), detDir->GetName(), taskDir->GetName())) ;
1310 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", expertHistKey->GetName())) ;
1311 if ( expertOdata->IsA()->InheritsFrom("TH1") ) {
1312 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", expertHistKey->GetName())) ;
1313 TH1 * hExpertdata = static_cast<TH1*>(expertOdata) ;
1314 listTaskQAD->Add(hExpertdata) ;
1319 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", histKey->GetName())) ;
1320 if ( odata->IsA()->InheritsFrom("TH1") ) {
1321 AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding %s", histKey->GetName())) ;
1322 TH1 * hdata = static_cast<TH1*>(odata) ;
1323 listTaskQAD->Add(hdata) ;
1328 Instance()->Put(listDetQAD, idr, &mdr) ;
1334 //_____________________________________________________________________________
1335 void AliQAManager::SetEventSpecie(AliRecoParam::EventSpecie_t es)
1337 // set the current event specie and inform AliQAv1 that this event specie has been encountered
1338 AliQAv1::Instance()->SetEventSpecie(es) ;
1341 //_____________________________________________________________________________
1342 void AliQAManager::SetRecoParam(const Int_t det, const AliDetectorRecoParam *par)
1344 // Set custom reconstruction parameters for a given detector
1345 // Single set of parameters for all the events
1346 GetQADataMaker(det)->SetRecoParam(par) ;
1349 //_____________________________________________________________________________
1350 void AliQAManager::SetWriteExpert()
1352 // enable the writing of QA expert data
1353 for (UInt_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1354 if (IsSelected(AliQAv1::GetDetName(iDet)))
1355 fQAWriteExpert[iDet] = kTRUE ;