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 **************************************************************************/
17 // Author: Andrei Gheata, 31/05/2006
19 //==============================================================================
20 // AliAnalysisManager - Manager analysis class. Allows creation of several
21 // analysis tasks and data containers storing their input/output. Allows
22 // connecting/chaining tasks via shared data containers. Serializes the current
23 // event for all tasks depending only on initial input data.
24 //==============================================================================
26 //==============================================================================
28 #include "AliAnalysisManager.h"
31 #include <Riostream.h>
37 #include <TMethodCall.h>
42 #include <TStopwatch.h>
44 #include "AliAnalysisSelector.h"
45 #include "AliAnalysisGrid.h"
46 #include "AliAnalysisTask.h"
47 #include "AliAnalysisDataContainer.h"
48 #include "AliAnalysisDataSlot.h"
49 #include "AliVEventHandler.h"
50 #include "AliVEventPool.h"
51 #include "AliSysInfo.h"
52 #include "AliAnalysisStatistics.h"
54 ClassImp(AliAnalysisManager)
56 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
57 TString AliAnalysisManager::fgCommonFileName = "";
58 Int_t AliAnalysisManager::fPBUpdateFreq = 1;
60 //______________________________________________________________________________
61 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
64 fInputEventHandler(NULL),
65 fOutputEventHandler(NULL),
66 fMCtruthEventHandler(NULL),
70 fMode(kLocalAnalysis),
74 fSpecialOutputLocation(""),
87 fAutoBranchHandling(kTRUE),
96 // Default constructor.
97 fgAnalysisManager = this;
98 fgCommonFileName = "AnalysisResults.root";
99 fTasks = new TObjArray();
100 fTopTasks = new TObjArray();
101 fZombies = new TObjArray();
102 fContainers = new TObjArray();
103 fInputs = new TObjArray();
104 fOutputs = new TObjArray();
105 fParamCont = new TObjArray();
107 TObject::SetObjectStat(kFALSE);
110 //______________________________________________________________________________
111 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
114 fInputEventHandler(NULL),
115 fOutputEventHandler(NULL),
116 fMCtruthEventHandler(NULL),
121 fInitOK(other.fInitOK),
122 fIsRemote(other.fIsRemote),
123 fDebug(other.fDebug),
124 fSpecialOutputLocation(""),
137 fAutoBranchHandling(other.fAutoBranchHandling),
140 fNcalls(other.fNcalls),
141 fMaxEntries(other.fMaxEntries),
142 fStatisticsMsg(other.fStatisticsMsg),
143 fRequestedBranches(other.fRequestedBranches),
144 fStatistics(other.fStatistics)
147 fTasks = new TObjArray(*other.fTasks);
148 fTopTasks = new TObjArray(*other.fTopTasks);
149 fZombies = new TObjArray(*other.fZombies);
150 fContainers = new TObjArray(*other.fContainers);
151 fInputs = new TObjArray(*other.fInputs);
152 fOutputs = new TObjArray(*other.fOutputs);
153 fParamCont = new TObjArray(*other.fParamCont);
154 fgCommonFileName = "AnalysisResults.root";
155 fgAnalysisManager = this;
156 TObject::SetObjectStat(kFALSE);
159 //______________________________________________________________________________
160 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
163 if (&other != this) {
164 TNamed::operator=(other);
165 fInputEventHandler = other.fInputEventHandler;
166 fOutputEventHandler = other.fOutputEventHandler;
167 fMCtruthEventHandler = other.fMCtruthEventHandler;
168 fEventPool = other.fEventPool;
171 fNSysInfo = other.fNSysInfo;
173 fInitOK = other.fInitOK;
174 fIsRemote = other.fIsRemote;
175 fDebug = other.fDebug;
176 fTasks = new TObjArray(*other.fTasks);
177 fTopTasks = new TObjArray(*other.fTopTasks);
178 fZombies = new TObjArray(*other.fZombies);
179 fContainers = new TObjArray(*other.fContainers);
180 fInputs = new TObjArray(*other.fInputs);
181 fOutputs = new TObjArray(*other.fOutputs);
182 fParamCont = new TObjArray(*other.fParamCont);
184 fCommonOutput = NULL;
187 fExtraFiles = other.fExtraFiles;
188 fgCommonFileName = "AnalysisResults.root";
189 fgAnalysisManager = this;
190 fAutoBranchHandling = other.fAutoBranchHandling;
191 fTable.Clear("nodelete");
192 fRunFromPath = other.fRunFromPath;
193 fNcalls = other. fNcalls;
194 fMaxEntries = other.fMaxEntries;
195 fStatisticsMsg = other.fStatisticsMsg;
196 fRequestedBranches = other.fRequestedBranches;
197 fStatistics = other.fStatistics;
202 //______________________________________________________________________________
203 AliAnalysisManager::~AliAnalysisManager()
206 if (fTasks) {fTasks->Delete(); delete fTasks;}
207 if (fTopTasks) delete fTopTasks;
208 if (fZombies) delete fZombies;
209 if (fContainers) {fContainers->Delete(); delete fContainers;}
210 if (fInputs) delete fInputs;
211 if (fOutputs) delete fOutputs;
212 if (fParamCont) delete fParamCont;
213 if (fGridHandler) delete fGridHandler;
214 if (fInputEventHandler) delete fInputEventHandler;
215 if (fOutputEventHandler) delete fOutputEventHandler;
216 if (fMCtruthEventHandler) delete fMCtruthEventHandler;
217 if (fEventPool) delete fEventPool;
218 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
219 TObject::SetObjectStat(kTRUE);
222 //______________________________________________________________________________
223 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
225 // Read one entry of the tree or a whole branch.
226 fCurrentEntry = entry;
227 if (!fAutoBranchHandling)
229 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : -1;
232 //______________________________________________________________________________
233 Int_t AliAnalysisManager::GetRunFromAlienPath(const char *path)
235 // Attempt to extract run number from input data path. Works only for paths to
236 // alice data in alien.
237 // sim: /alice/sim/<production>/run_no/...
238 // data: /alice/data/year/period/000run_no/... (ESD or AOD)
242 Int_t index = s.Index("/alice/sim");
244 for (Int_t i=0; i<3; i++) {
245 index = s.Index("/", index+1);
246 if (index<0) return 0;
251 index = s.Index("/alice/data");
253 for (Int_t i=0; i<4; i++) {
254 index = s.Index("/", index+1);
255 if (index<0) return 0;
263 //______________________________________________________________________________
264 Bool_t AliAnalysisManager::Init(TTree *tree)
266 // The Init() function is called when the selector needs to initialize
267 // a new tree or chain. Typically here the branch addresses of the tree
268 // will be set. It is normaly not necessary to make changes to the
269 // generated code, but the routine can be extended by the user if needed.
270 // Init() will be called many times when running with PROOF.
271 Bool_t init = kFALSE;
272 if (!tree) return kFALSE; // Should not happen - protected in selector caller
274 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
276 // Call InitTree of EventHandler
277 if (fOutputEventHandler) {
278 if (fMode == kProofAnalysis) {
279 init = fOutputEventHandler->Init(0x0, "proof");
281 init = fOutputEventHandler->Init(0x0, "local");
284 Error("Init", "Output event handler failed to initialize");
289 if (fInputEventHandler) {
290 if (fMode == kProofAnalysis) {
291 init = fInputEventHandler->Init(tree, "proof");
293 init = fInputEventHandler->Init(tree, "local");
296 Error("Init", "Input event handler failed to initialize tree");
300 // If no input event handler we need to get the tree once
302 if(!tree->GetTree()) {
303 Long64_t readEntry = tree->LoadTree(0);
304 if (readEntry == -2) {
305 Error("Init", "Input tree has no entry. Exiting");
311 if (fMCtruthEventHandler) {
312 if (fMode == kProofAnalysis) {
313 init = fMCtruthEventHandler->Init(0x0, "proof");
315 init = fMCtruthEventHandler->Init(0x0, "local");
318 Error("Init", "MC event handler failed to initialize");
323 if (!fInitOK) InitAnalysis();
324 if (!fInitOK) return kFALSE;
327 AliAnalysisDataContainer *top = fCommonInput;
328 if (!top) top = (AliAnalysisDataContainer*)fInputs->At(0);
330 Error("Init","No top input container !");
334 CheckBranches(kFALSE);
336 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
341 //______________________________________________________________________________
342 void AliAnalysisManager::SlaveBegin(TTree *tree)
344 // The SlaveBegin() function is called after the Begin() function.
345 // When running with PROOF SlaveBegin() is called on each slave server.
346 // The tree argument is deprecated (on PROOF 0 is passed).
347 if (fDebug > 1) printf("->AliAnalysisManager::SlaveBegin()\n");
348 if (!CheckTasks()) Fatal("SlaveBegin", "Not all needed libraries were loaded");
349 static Bool_t isCalled = kFALSE;
350 Bool_t init = kFALSE;
351 Bool_t initOK = kTRUE;
353 TDirectory *curdir = gDirectory;
354 // Call SlaveBegin only once in case of mixing
355 if (isCalled && fMode==kMixingAnalysis) return;
357 // Call Init of EventHandler
358 if (fOutputEventHandler) {
359 if (fMode == kProofAnalysis) {
360 // Merging AOD's in PROOF via TProofOutputFile
361 if (fDebug > 1) printf(" Initializing AOD output file %s...\n", fOutputEventHandler->GetOutputFileName());
362 init = fOutputEventHandler->Init("proof");
363 if (!init) msg = "Failed to initialize output handler on worker";
365 init = fOutputEventHandler->Init("local");
366 if (!init) msg = "Failed to initialize output handler";
369 if (!fSelector) Error("SlaveBegin", "Selector not set");
370 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
373 if (fInputEventHandler) {
374 fInputEventHandler->SetInputTree(tree);
375 if (fMode == kProofAnalysis) {
376 init = fInputEventHandler->Init("proof");
377 if (!init) msg = "Failed to initialize input handler on worker";
379 init = fInputEventHandler->Init("local");
380 if (!init) msg = "Failed to initialize input handler";
383 if (!fSelector) Error("SlaveBegin", "Selector not set");
384 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
387 if (fMCtruthEventHandler) {
388 if (fMode == kProofAnalysis) {
389 init = fMCtruthEventHandler->Init("proof");
390 if (!init) msg = "Failed to initialize MC handler on worker";
392 init = fMCtruthEventHandler->Init("local");
393 if (!init) msg = "Failed to initialize MC handler";
396 if (!fSelector) Error("SlaveBegin", "Selector not set");
397 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
399 if (curdir) curdir->cd();
403 AliAnalysisTask *task;
404 // Call CreateOutputObjects for all tasks
405 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
406 Bool_t dirStatus = TH1::AddDirectoryStatus();
408 while ((task=(AliAnalysisTask*)next())) {
410 // Start with memory as current dir and make sure by default histograms do not get attached to files.
411 TH1::AddDirectory(kFALSE);
412 task->CreateOutputObjects();
413 if (!task->CheckPostData()) {
414 Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
415 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
416 ####### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ##########", task->GetName(), task->ClassName());
418 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
421 TH1::AddDirectory(dirStatus);
422 if (curdir) curdir->cd();
423 if (fDebug > 1) printf("<-AliAnalysisManager::SlaveBegin()\n");
426 //______________________________________________________________________________
427 Bool_t AliAnalysisManager::Notify()
429 // The Notify() function is called when a new file is opened. This
430 // can be either for a new TTree in a TChain or when when a new TTree
431 // is started when using PROOF. It is normaly not necessary to make changes
432 // to the generated code, but the routine can be extended by the
433 // user if needed. The return value is currently not used.
434 if (!fTree) return kFALSE;
436 fTable.Clear("nodelete"); // clearing the hash table may not be needed -> C.L.
437 if (fMode == kProofAnalysis) fIsRemote = kTRUE;
439 TFile *curfile = fTree->GetCurrentFile();
441 Error("Notify","No current file");
445 if (fDebug > 1) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
446 Int_t run = AliAnalysisManager::GetRunFromAlienPath(curfile->GetName());
447 if (run) SetRunFromPath(run);
448 if (fDebug > 1) printf(" ### run found from path: %d\n", run);
450 AliAnalysisTask *task;
452 // Call Notify of the event handlers
453 if (fInputEventHandler) {
454 fInputEventHandler->Notify(curfile->GetName());
457 if (fOutputEventHandler) {
458 fOutputEventHandler->Notify(curfile->GetName());
461 if (fMCtruthEventHandler) {
462 fMCtruthEventHandler->Notify(curfile->GetName());
465 // Call Notify for all tasks
466 while ((task=(AliAnalysisTask*)next()))
469 if (fDebug > 1) printf("<-AliAnalysisManager::Notify()\n");
473 //______________________________________________________________________________
474 Bool_t AliAnalysisManager::Process(Long64_t entry)
476 // The Process() function is called for each entry in the tree (or possibly
477 // keyed object in the case of PROOF) to be processed. The entry argument
478 // specifies which entry in the currently loaded tree is to be processed.
479 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
480 // to read either all or the required parts of the data. When processing
481 // keyed objects with PROOF, the object is already loaded and is available
482 // via the fObject pointer.
484 // This function should contain the "body" of the analysis. It can contain
485 // simple or elaborate selection criteria, run algorithms on the data
486 // of the event and typically fill histograms.
488 // WARNING when a selector is used with a TChain, you must use
489 // the pointer to the current TTree to call GetEntry(entry).
490 // The entry is always the local entry number in the current tree.
491 // Assuming that fChain is the pointer to the TChain being processed,
492 // use fChain->GetTree()->GetEntry(entry).
493 if (fDebug > 1) printf("->AliAnalysisManager::Process(%lld)\n", entry);
495 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
496 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
497 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
501 if (fInputEventHandler) fInputEventHandler ->GetEntry();
504 if (fDebug > 1) printf("<-AliAnalysisManager::Process()\n");
508 //______________________________________________________________________________
509 void AliAnalysisManager::PackOutput(TList *target)
511 // Pack all output data containers in the output list. Called at SlaveTerminate
512 // stage in PROOF case for each slave.
513 if (fDebug > 1) printf("->AliAnalysisManager::PackOutput()\n");
515 Error("PackOutput", "No target. Exiting.");
518 TDirectory *cdir = gDirectory;
520 if (fInputEventHandler) fInputEventHandler ->Terminate();
521 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
522 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
525 // Call FinishTaskOutput() for each event loop task (not called for
526 // post-event loop tasks - use Terminate() fo those)
527 TIter nexttask(fTasks);
528 AliAnalysisTask *task;
529 while ((task=(AliAnalysisTask*)nexttask())) {
530 if (!task->IsPostEventLoop()) {
531 if (fDebug > 1) printf("->FinishTaskOutput: task %s\n", task->GetName());
532 task->FinishTaskOutput();
534 if (fDebug > 1) printf("<-FinishTaskOutput: task %s\n", task->GetName());
537 // Write statistics message on the workers.
538 if (fStatistics) WriteStatisticsMsg(fNcalls);
540 if (fMode == kProofAnalysis) {
541 TIter next(fOutputs);
542 AliAnalysisDataContainer *output;
543 Bool_t isManagedByHandler = kFALSE;
546 while ((output=(AliAnalysisDataContainer*)next())) {
547 // Do not consider outputs of post event loop tasks
548 isManagedByHandler = kFALSE;
549 if (output->GetProducer() && output->GetProducer()->IsPostEventLoop()) continue;
550 const char *filename = output->GetFileName();
551 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
552 isManagedByHandler = kTRUE;
553 printf("#### Handler output. Extra: %s\n", fExtraFiles.Data());
554 filename = fOutputEventHandler->GetOutputFileName();
556 // Check if data was posted to this container. If not, issue an error.
557 if (!output->GetData() && !isManagedByHandler) {
558 Error("PackOutput", "No data for output container %s. Forgot to PostData ?", output->GetName());
561 if (!output->IsSpecialOutput()) {
563 if (strlen(filename) && !isManagedByHandler) {
564 // Backup current folder
565 TDirectory *opwd = gDirectory;
566 // File resident outputs.
567 // Check first if the file exists.
568 TString openoption = "RECREATE";
569 Bool_t firsttime = kTRUE;
570 if (filestmp.FindObject(output->GetFileName())) {
573 filestmp.Add(new TNamed(output->GetFileName(),""));
575 if (!gSystem->AccessPathName(output->GetFileName()) && !firsttime) openoption = "UPDATE";
576 // TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
577 // Save data to file, then close.
578 if (output->GetData()->InheritsFrom(TCollection::Class())) {
579 // If data is a collection, we set the name of the collection
580 // as the one of the container and we save as a single key.
581 TCollection *coll = (TCollection*)output->GetData();
582 coll->SetName(output->GetName());
583 // coll->Write(output->GetName(), TObject::kSingleKey);
585 if (output->GetData()->InheritsFrom(TTree::Class())) {
586 TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
587 // Save data to file, then close.
588 TTree *tree = (TTree*)output->GetData();
589 // Check if tree is in memory
590 if (tree->GetDirectory()==gROOT) tree->SetDirectory(gDirectory);
594 // output->GetData()->Write();
597 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
599 // printf(" file %s listing content:\n", filename);
602 // Clear file list to release object ownership to user.
605 output->SetFile(NULL);
606 // Restore current directory
607 if (opwd) opwd->cd();
609 // Memory-resident outputs
610 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
612 AliAnalysisDataWrapper *wrap = 0;
613 if (isManagedByHandler) {
614 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
615 wrap->SetName(output->GetName());
617 else wrap =output->ExportData();
618 // Output wrappers must NOT delete data after merging - the user owns them
619 wrap->SetDeleteData(kFALSE);
622 // Special outputs. The file must be opened and connected to the container.
623 TDirectory *opwd = gDirectory;
624 TFile *file = output->GetFile();
626 AliAnalysisTask *producer = output->GetProducer();
628 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
629 output->GetFileName(), output->GetName(), producer->ClassName());
632 TString outFilename = file->GetName();
633 if (fDebug > 1) printf("PackOutput %s: special output\n", output->GetName());
634 if (isManagedByHandler) {
635 // Terminate IO for files managed by the output handler
636 // file->Write() moved to AOD handler (A.G. 11.01.10)
637 // if (file) file->Write();
638 if (file && fDebug > 2) {
639 printf(" handled file %s listing content:\n", file->GetName());
642 fOutputEventHandler->TerminateIO();
645 // Release object ownership to users after writing data to file
646 if (output->GetData()->InheritsFrom(TCollection::Class())) {
647 // If data is a collection, we set the name of the collection
648 // as the one of the container and we save as a single key.
649 TCollection *coll = (TCollection*)output->GetData();
650 coll->SetName(output->GetName());
651 coll->Write(output->GetName(), TObject::kSingleKey);
653 if (output->GetData()->InheritsFrom(TTree::Class())) {
654 TTree *tree = (TTree*)output->GetData();
655 tree->SetDirectory(file);
658 output->GetData()->Write();
662 printf(" file %s listing content:\n", output->GetFileName());
665 // Clear file list to release object ownership to user.
668 output->SetFile(NULL);
670 // Restore current directory
671 if (opwd) opwd->cd();
672 // Check if a special output location was provided or the output files have to be merged
673 if (strlen(fSpecialOutputLocation.Data())) {
674 TString remote = fSpecialOutputLocation;
676 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
677 if (remote.BeginsWith("alien:")) {
678 gROOT->ProcessLine("TGrid::Connect(\"alien:\", gProofServ->GetUser());");
679 remote += outFilename;
680 remote.ReplaceAll(".root", Form("_%d.root", gid));
682 remote += Form("%s_%d_", gSystem->HostName(), gid);
683 remote += outFilename;
686 Info("PackOutput", "Output file for container %s to be copied \n at: %s. No merging.",
687 output->GetName(), remote.Data());
688 TFile::Cp ( outFilename.Data(), remote.Data() );
689 // Copy extra outputs
690 if (fExtraFiles.Length() && isManagedByHandler) {
691 TObjArray *arr = fExtraFiles.Tokenize(" ");
693 TIter nextfilename(arr);
694 while ((os=(TObjString*)nextfilename())) {
695 outFilename = os->GetString();
696 remote = fSpecialOutputLocation;
698 if (remote.BeginsWith("alien://")) {
699 remote += outFilename;
700 remote.ReplaceAll(".root", Form("_%d.root", gid));
702 remote += Form("%s_%d_", gSystem->HostName(), gid);
703 remote += outFilename;
706 Info("PackOutput", "Extra AOD file %s to be copied \n at: %s. No merging.",
707 outFilename.Data(), remote.Data());
708 TFile::Cp ( outFilename.Data(), remote.Data() );
713 // No special location specified-> use TProofOutputFile as merging utility
714 // The file at this output slot must be opened in CreateOutputObjects
715 if (fDebug > 1) printf(" File for container %s to be merged via file merger...\n", output->GetName());
721 if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
724 //______________________________________________________________________________
725 void AliAnalysisManager::ImportWrappers(TList *source)
727 // Import data in output containers from wrappers coming in source.
728 if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
729 TIter next(fOutputs);
730 AliAnalysisDataContainer *cont;
731 AliAnalysisDataWrapper *wrap;
733 Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
734 TDirectory *cdir = gDirectory;
735 while ((cont=(AliAnalysisDataContainer*)next())) {
737 if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
738 if (cont->IsRegisterDataset()) continue;
739 const char *filename = cont->GetFileName();
740 Bool_t isManagedByHandler = kFALSE;
741 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
742 isManagedByHandler = kTRUE;
743 filename = fOutputEventHandler->GetOutputFileName();
745 if (cont->IsSpecialOutput() || inGrid) {
746 if (strlen(fSpecialOutputLocation.Data())) continue;
747 // Copy merged file from PROOF scratch space.
748 // In case of grid the files are already in the current directory.
750 if (isManagedByHandler && fExtraFiles.Length()) {
751 // Copy extra registered dAOD files.
752 TObjArray *arr = fExtraFiles.Tokenize(" ");
754 TIter nextfilename(arr);
755 while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
758 if (!GetFileFromWrapper(filename, source)) continue;
760 // Normally we should connect data from the copied file to the
761 // corresponding output container, but it is not obvious how to do this
762 // automatically if several objects in file...
763 TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
764 if (!f) f = TFile::Open(filename, "READ");
766 Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
770 // Cd to the directory pointed by the container
771 TString folder = cont->GetFolderName();
772 if (!folder.IsNull()) f->cd(folder);
773 // Try to fetch first an object having the container name.
774 obj = gDirectory->Get(cont->GetName());
776 Warning("ImportWrappers", "Could not import object of type:%s for container %s in file %s:%s.\n Object will not be available in Terminate(). Try if possible to name the output object as the container (%s) or to embed it in a TList",
777 cont->GetType()->GetName(), cont->GetName(), filename, cont->GetFolderName(), cont->GetName());
780 wrap = new AliAnalysisDataWrapper(obj);
781 wrap->SetDeleteData(kFALSE);
783 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
785 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
790 printf(" Importing data for container %s\n", cont->GetName());
791 if (strlen(filename)) printf(" -> file %s\n", filename);
794 cont->ImportData(wrap);
796 if (cdir) cdir->cd();
797 if (fDebug > 1) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
800 //______________________________________________________________________________
801 void AliAnalysisManager::UnpackOutput(TList *source)
803 // Called by AliAnalysisSelector::Terminate only on the client.
804 if (fDebug > 1) printf("->AliAnalysisManager::UnpackOutput()\n");
806 Error("UnpackOutput", "No target. Exiting.");
809 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
811 if (fMode == kProofAnalysis) ImportWrappers(source);
813 TIter next(fOutputs);
814 AliAnalysisDataContainer *output;
815 while ((output=(AliAnalysisDataContainer*)next())) {
816 if (!output->GetData()) continue;
817 // Check if there are client tasks that run post event loop
818 if (output->HasConsumers()) {
819 // Disable event loop semaphore
820 output->SetPostEventLoop(kTRUE);
821 TObjArray *list = output->GetConsumers();
822 Int_t ncons = list->GetEntriesFast();
823 for (Int_t i=0; i<ncons; i++) {
824 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
825 task->CheckNotify(kTRUE);
826 // If task is active, execute it
827 if (task->IsPostEventLoop() && task->IsActive()) {
828 if (fDebug > 1) printf("== Executing post event loop task %s\n", task->GetName());
834 if (fDebug > 1) printf("<-AliAnalysisManager::UnpackOutput()\n");
837 //______________________________________________________________________________
838 void AliAnalysisManager::Terminate()
840 // The Terminate() function is the last function to be called during
841 // a query. It always runs on the client, it can be used to present
842 // the results graphically.
843 if (fDebug > 1) printf("->AliAnalysisManager::Terminate()\n");
844 TDirectory *cdir = gDirectory;
846 AliAnalysisTask *task;
847 AliAnalysisDataContainer *output;
850 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
851 // Call Terminate() for tasks
853 while (!IsSkipTerminate() && (task=(AliAnalysisTask*)next())) {
854 // Save all the canvases produced by the Terminate
855 TString pictname = Form("%s_%s", task->GetName(), task->ClassName());
859 AliSysInfo::AddStamp(Form("%s_TERMINATE",task->ClassName()),0, itask, 2);
861 if (TObject::TestBit(kSaveCanvases)) {
862 if (!gROOT->IsBatch()) {
863 if (fDebug>1) printf("Waiting 5 sec for %s::Terminate() to finish drawing ...\n", task->ClassName());
865 while (timer.CpuTime()<5) {
867 gSystem->ProcessEvents();
870 Int_t iend = gROOT->GetListOfCanvases()->GetEntries();
871 if (iend==0) continue;
873 for (Int_t ipict=0; ipict<iend; ipict++) {
874 canvas = (TCanvas*)gROOT->GetListOfCanvases()->At(ipict);
875 if (!canvas) continue;
876 canvas->SaveAs(Form("%s_%02d.gif", pictname.Data(),ipict));
878 gROOT->GetListOfCanvases()->Delete();
882 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
883 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
884 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
886 TObjArray *allOutputs = new TObjArray();
888 for (icont=0; icont<fOutputs->GetEntriesFast(); icont++) allOutputs->Add(fOutputs->At(icont));
889 if (!IsSkipTerminate())
890 for (icont=0; icont<fParamCont->GetEntriesFast(); icont++) allOutputs->Add(fParamCont->At(icont));
891 TIter next1(allOutputs);
892 TString handlerFile = "";
893 TString extraOutputs = "";
894 if (fOutputEventHandler) {
895 handlerFile = fOutputEventHandler->GetOutputFileName();
896 extraOutputs = fOutputEventHandler->GetExtraOutputs();
900 while ((output=(AliAnalysisDataContainer*)next1())) {
901 // Special outputs or grid files have the files already closed and written.
903 if (fMode == kGridAnalysis && icont<=fOutputs->GetEntriesFast()) continue;
904 if (fMode == kProofAnalysis) {
905 if (output->IsSpecialOutput() || output->IsRegisterDataset()) continue;
907 const char *filename = output->GetFileName();
908 TString openoption = "RECREATE";
909 if (!(strcmp(filename, "default"))) continue;
910 if (!strlen(filename)) continue;
911 if (!output->GetData()) continue;
912 TDirectory *opwd = gDirectory;
913 TFile *file = output->GetFile();
914 if (!file) file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
916 //if (handlerFile == filename && !gSystem->AccessPathName(filename)) openoption = "UPDATE";
917 Bool_t firsttime = kTRUE;
918 if (filestmp.FindObject(filename) || extraOutputs.Contains(filename)) {
921 filestmp.Add(new TNamed(filename,""));
923 if (!gSystem->AccessPathName(filename) && !firsttime) openoption = "UPDATE";
924 if (fDebug>1) printf("Opening file: %s option=%s\n",filename, openoption.Data());
925 file = new TFile(filename, openoption);
927 if (fDebug>1) printf("File <%s> already opened with option: <%s> \n", filename, file->GetOption());
928 openoption = file->GetOption();
929 if (openoption == "READ") {
930 if (fDebug>1) printf("...reopening in UPDATE mode\n");
931 file->ReOpen("UPDATE");
934 if (file->IsZombie()) {
935 Error("Terminate", "Cannot open output file %s", filename);
938 output->SetFile(file);
940 // Check for a folder request
941 TString dir = output->GetFolderName();
943 if (!file->GetDirectory(dir)) file->mkdir(dir);
946 if (fDebug > 1) printf("...writing container %s to file %s:%s\n", output->GetName(), file->GetName(), output->GetFolderName());
947 if (output->GetData()->InheritsFrom(TCollection::Class())) {
948 // If data is a collection, we set the name of the collection
949 // as the one of the container and we save as a single key.
950 TCollection *coll = (TCollection*)output->GetData();
951 coll->SetName(output->GetName());
952 coll->Write(output->GetName(), TObject::kSingleKey);
954 if (output->GetData()->InheritsFrom(TTree::Class())) {
955 TTree *tree = (TTree*)output->GetData();
956 tree->SetDirectory(gDirectory);
959 output->GetData()->Write();
962 if (opwd) opwd->cd();
966 while ((output=(AliAnalysisDataContainer*)next1())) {
967 // Close all files at output
968 TDirectory *opwd = gDirectory;
969 if (output->GetFile()) {
970 // Clear file list to release object ownership to user.
971 // output->GetFile()->Clear();
972 output->GetFile()->Close();
973 output->SetFile(NULL);
974 // Copy merged outputs in alien if requested
975 if (fSpecialOutputLocation.Length() &&
976 fSpecialOutputLocation.BeginsWith("alien://")) {
977 Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data());
978 TFile::Cp(output->GetFile()->GetName(),
979 Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
982 if (opwd) opwd->cd();
985 //Write statistics information on the client
986 if (fStatistics) WriteStatisticsMsg(fNcalls);
988 TDirectory *crtdir = gDirectory;
989 TFile f("syswatch.root", "RECREATE");
993 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
994 tree->SetName("syswatch");
995 tree->SetMarkerStyle(kCircle);
996 tree->SetMarkerColor(kBlue);
997 tree->SetMarkerSize(0.5);
998 if (!gROOT->IsBatch()) {
999 tree->SetAlias("event", "id0");
1000 tree->SetAlias("task", "id1");
1001 tree->SetAlias("stage", "id2");
1002 // Already defined aliases
1003 // tree->SetAlias("deltaT","stampSec-stampOldSec");
1004 // tree->SetAlias("T","stampSec-first");
1005 // tree->SetAlias("deltaVM","(pI.fMemVirtual-pIOld.fMemVirtual)");
1006 // tree->SetAlias("VM","pI.fMemVirtual");
1007 TCanvas *canvas = new TCanvas("SysInfo","SysInfo",10,10,1200,800);
1008 Int_t npads = 1 /*COO plot for all tasks*/ +
1009 fTopTasks->GetEntries() /*Exec plot per task*/ +
1010 1 /*Terminate plot for all tasks*/ +
1013 Int_t iopt = (Int_t)TMath::Sqrt((Double_t)npads);
1014 if (npads<iopt*(iopt+1))
1015 canvas->Divide(iopt, iopt+1, 0.01, 0.01);
1017 canvas->Divide(iopt+1, iopt+1, 0.01, 0.01);
1019 // draw the plot of deltaVM for Exec for each task
1020 for (itask=0; itask<fTopTasks->GetEntriesFast(); itask++) {
1021 task = (AliAnalysisTask*)fTopTasks->At(itask);
1023 cut = Form("task==%d && stage==1", itask);
1024 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1025 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1027 hist->SetTitle(Form("%s: Exec dVM[kB]/event", task->GetName()));
1028 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1031 // Draw the plot of deltaVM for CreateOutputObjects for all tasks
1033 tree->SetMarkerStyle(kFullTriangleUp);
1034 tree->SetMarkerColor(kRed);
1035 tree->SetMarkerSize(0.8);
1036 cut = "task>=0 && task<1000 && stage==0";
1037 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1038 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1040 hist->SetTitle("Memory in CreateOutputObjects()");
1041 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1042 hist->GetXaxis()->SetTitle("task");
1044 // draw the plot of deltaVM for Terminate for all tasks
1046 tree->SetMarkerStyle(kOpenSquare);
1047 tree->SetMarkerColor(kMagenta);
1048 cut = "task>=0 && task<1000 && stage==2";
1049 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1050 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1052 hist->SetTitle("Memory in Terminate()");
1053 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1054 hist->GetXaxis()->SetTitle("task");
1058 tree->SetMarkerStyle(kFullCircle);
1059 tree->SetMarkerColor(kGreen);
1060 cut = Form("task==%d && stage==1",fTopTasks->GetEntriesFast()-1);
1061 tree->Draw("VM:event",cut,"", 1234567890, 0);
1062 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1064 hist->SetTitle("Virtual memory");
1065 hist->GetYaxis()->SetTitle("VM [MB]");
1069 tree->SetMarkerStyle(kCircle);
1070 tree->SetMarkerColor(kBlue);
1071 tree->SetMarkerSize(0.5);
1076 if (crtdir) crtdir->cd();
1078 // Validate the output files
1079 if (ValidateOutputFiles() && fIsRemote && fMode!=kProofAnalysis) {
1081 out.open("outputs_valid", ios::out);
1085 if (fDebug > 1) printf("<-AliAnalysisManager::Terminate()\n");
1087 //______________________________________________________________________________
1088 void AliAnalysisManager::ProfileTask(Int_t itop, const char *option) const
1090 // Profiles the task having the itop index in the list of top (first level) tasks.
1091 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->At(itop);
1093 Error("ProfileTask", "There are only %d top tasks in the manager", fTopTasks->GetEntries());
1096 ProfileTask(task->GetName(), option);
1099 //______________________________________________________________________________
1100 void AliAnalysisManager::ProfileTask(const char *name, const char */*option*/) const
1102 // Profile a managed task after the execution of the analysis in case NSysInfo
1104 if (gSystem->AccessPathName("syswatch.root")) {
1105 Error("ProfileTask", "No file syswatch.root found in the current directory");
1108 if (gROOT->IsBatch()) return;
1109 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->FindObject(name);
1111 Error("ProfileTask", "No top task named %s known by the manager.", name);
1114 Int_t itop = fTopTasks->IndexOf(task);
1115 Int_t itask = fTasks->IndexOf(task);
1116 // Create canvas with 2 pads: first draw COO + Terminate, second Exec
1117 TDirectory *cdir = gDirectory;
1118 TFile f("syswatch.root");
1119 TTree *tree = (TTree*)f.Get("syswatch");
1121 Error("ProfileTask", "No tree named <syswatch> found in file syswatch.root");
1124 if (fDebug > 1) printf("=== Profiling task %s (class %s)\n", name, task->ClassName());
1125 TCanvas *canvas = new TCanvas(Form("profile_%d",itop),Form("Profile of task %s (class %s)",name,task->ClassName()),10,10,800,600);
1126 canvas->Divide(2, 2, 0.01, 0.01);
1130 // VM profile for COO and Terminate methods
1132 cut = Form("task==%d && (stage==0 || stage==2)",itask);
1133 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1134 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1136 hist->SetTitle("Alocated VM[kB] for COO and Terminate");
1137 hist->GetYaxis()->SetTitle("deltaVM [kB]");
1138 hist->GetXaxis()->SetTitle("method");
1140 // CPU profile per event
1142 cut = Form("task==%d && stage==1",itop);
1143 tree->Draw("deltaT:event",cut,"", 1234567890, 0);
1144 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1146 hist->SetTitle("Execution time per event");
1147 hist->GetYaxis()->SetTitle("CPU/event [s]");
1149 // VM profile for Exec
1151 cut = Form("task==%d && stage==1",itop);
1152 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1153 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1155 hist->SetTitle("Alocated VM[kB] per event");
1156 hist->GetYaxis()->SetTitle("deltaVM [kB]");
1161 if (cdir) cdir->cd();
1164 //______________________________________________________________________________
1165 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1167 // Adds a user task to the global list of tasks.
1168 if (fTasks->FindObject(task)) {
1169 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1172 task->SetActive(kFALSE);
1176 //______________________________________________________________________________
1177 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1179 // Retreive task by name.
1180 if (!fTasks) return NULL;
1181 return (AliAnalysisTask*)fTasks->FindObject(name);
1184 //______________________________________________________________________________
1185 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
1186 TClass *datatype, EAliAnalysisContType type, const char *filename)
1188 // Create a data container of a certain type. Types can be:
1189 // kExchangeContainer = 0, used to exchange data between tasks
1190 // kInputContainer = 1, used to store input data
1191 // kOutputContainer = 2, used for writing result to a file
1192 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1193 // the output object to a folder inside the output file
1194 if (fContainers->FindObject(name)) {
1195 Error("CreateContainer","A container named %s already defined !",name);
1198 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1199 fContainers->Add(cont);
1201 case kInputContainer:
1204 case kOutputContainer:
1205 fOutputs->Add(cont);
1206 if (filename && strlen(filename)) {
1207 cont->SetFileName(filename);
1208 cont->SetDataOwned(kFALSE); // data owned by the file
1211 case kParamContainer:
1212 fParamCont->Add(cont);
1213 if (filename && strlen(filename)) {
1214 cont->SetFileName(filename);
1215 cont->SetDataOwned(kFALSE); // data owned by the file
1218 case kExchangeContainer:
1224 //______________________________________________________________________________
1225 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1226 AliAnalysisDataContainer *cont)
1228 // Connect input of an existing task to a data container.
1230 Error("ConnectInput", "Task pointer is NULL");
1233 if (!fTasks->FindObject(task)) {
1235 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1237 Bool_t connected = task->ConnectInput(islot, cont);
1241 //______________________________________________________________________________
1242 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1243 AliAnalysisDataContainer *cont)
1245 // Connect output of an existing task to a data container.
1247 Error("ConnectOutput", "Task pointer is NULL");
1250 if (!fTasks->FindObject(task)) {
1252 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1254 Bool_t connected = task->ConnectOutput(islot, cont);
1258 //______________________________________________________________________________
1259 void AliAnalysisManager::CleanContainers()
1261 // Clean data from all containers that have already finished all client tasks.
1262 TIter next(fContainers);
1263 AliAnalysisDataContainer *cont;
1264 while ((cont=(AliAnalysisDataContainer *)next())) {
1265 if (cont->IsOwnedData() &&
1266 cont->IsDataReady() &&
1267 cont->ClientsExecuted()) cont->DeleteData();
1271 //______________________________________________________________________________
1272 Bool_t AliAnalysisManager::InitAnalysis()
1274 // Initialization of analysis chain of tasks. Should be called after all tasks
1275 // and data containers are properly connected
1276 // Reset flag and remove valid_outputs file if exists
1278 if (!gSystem->AccessPathName("outputs_valid"))
1279 gSystem->Unlink("outputs_valid");
1280 // Check for top tasks (depending only on input data containers)
1281 if (!fTasks->First()) {
1282 Error("InitAnalysis", "Analysis has no tasks !");
1286 AliAnalysisTask *task;
1287 AliAnalysisDataContainer *cont;
1290 Bool_t iszombie = kFALSE;
1291 Bool_t istop = kTRUE;
1293 while ((task=(AliAnalysisTask*)next())) {
1296 Int_t ninputs = task->GetNinputs();
1297 for (i=0; i<ninputs; i++) {
1298 cont = task->GetInputSlot(i)->GetContainer();
1302 fZombies->Add(task);
1306 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
1307 i, task->GetName());
1309 if (iszombie) continue;
1310 // Check if cont is an input container
1311 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1312 // Connect to parent task
1316 fTopTasks->Add(task);
1320 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1323 // Check now if there are orphan tasks
1324 for (i=0; i<ntop; i++) {
1325 task = (AliAnalysisTask*)fTopTasks->At(i);
1330 while ((task=(AliAnalysisTask*)next())) {
1331 if (!task->IsUsed()) {
1333 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1336 // Check the task hierarchy (no parent task should depend on data provided
1337 // by a daughter task)
1338 for (i=0; i<ntop; i++) {
1339 task = (AliAnalysisTask*)fTopTasks->At(i);
1340 if (task->CheckCircularDeps()) {
1341 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1346 // Check that all containers feeding post-event loop tasks are in the outputs list
1347 TIter nextcont(fContainers); // loop over all containers
1348 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1349 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1350 if (cont->HasConsumers()) {
1351 // Check if one of the consumers is post event loop
1352 TIter nextconsumer(cont->GetConsumers());
1353 while ((task=(AliAnalysisTask*)nextconsumer())) {
1354 if (task->IsPostEventLoop()) {
1355 fOutputs->Add(cont);
1362 // Check if all special output containers have a file name provided
1363 TIter nextout(fOutputs);
1364 while ((cont=(AliAnalysisDataContainer*)nextout())) {
1365 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1366 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1370 // Initialize requested branch list if needed
1371 if (!fAutoBranchHandling) {
1373 while ((task=(AliAnalysisTask*)next())) {
1374 if (!task->HasBranches()) {
1375 Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches.\nUse: fBranchNames = \"ESD:br1,br2,...,brN AOD:bra1,bra2,...,braM\"",
1376 task->GetName(), task->ClassName());
1379 if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1380 Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1383 TString taskbranches;
1384 task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1385 if (taskbranches.IsNull()) {
1386 Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1387 task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1390 AddBranches(taskbranches);
1397 //______________________________________________________________________________
1398 void AliAnalysisManager::AddBranches(const char *branches)
1400 // Add branches to the existing fRequestedBranches.
1401 TString br(branches);
1402 TObjArray *arr = br.Tokenize(",");
1405 while ((obj=next())) {
1406 if (!fRequestedBranches.Contains(obj->GetName())) {
1407 if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1408 fRequestedBranches += obj->GetName();
1414 //______________________________________________________________________________
1415 void AliAnalysisManager::CheckBranches(Bool_t load)
1417 // The method checks the input branches to be loaded during the analysis.
1418 if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;
1419 TObjArray *arr = fRequestedBranches.Tokenize(",");
1422 while ((obj=next())) {
1423 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1425 br = fTree->GetBranch(obj->GetName());
1427 Error("CheckBranches", "Could not find branch %s",obj->GetName());
1432 if (load && br->GetReadEntry()!=GetCurrentEntry()) br->GetEntry(GetCurrentEntry());
1437 //______________________________________________________________________________
1438 Bool_t AliAnalysisManager::CheckTasks() const
1440 // Check consistency of tasks.
1441 // Get the pointer to AliAnalysisTaskSE::Class()
1442 TClass *badptr = (TClass*)gROOT->ProcessLine("AliAnalysisTaskSE::Class()");
1443 // Loop all tasks to check if their corresponding library was loaded
1446 while ((obj=next())) {
1447 if (obj->IsA() == badptr) {
1448 Error("CheckTasks", "##################\n \
1449 Class for task %s NOT loaded. You probably forgot to load the library for this task (or compile it dynamically).\n###########################\n",obj->GetName());
1456 //______________________________________________________________________________
1457 void AliAnalysisManager::PrintStatus(Option_t *option) const
1459 // Print task hierarchy.
1461 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1464 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1466 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1467 TIter next(fTopTasks);
1468 AliAnalysisTask *task;
1469 while ((task=(AliAnalysisTask*)next()))
1470 task->PrintTask(option);
1472 if (!fAutoBranchHandling && !fRequestedBranches.IsNull())
1473 printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1475 TString sopt(option);
1478 if (sopt.Contains("ALL"))
1480 if ( fOutputEventHandler )
1482 cout << TString('_',78) << endl;
1483 cout << "OutputEventHandler:" << endl;
1484 fOutputEventHandler->Print(" ");
1489 //______________________________________________________________________________
1490 void AliAnalysisManager::ResetAnalysis()
1492 // Reset all execution flags and clean containers.
1496 //______________________________________________________________________________
1497 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1499 // Start analysis having a grid handler.
1500 if (!fGridHandler) {
1501 Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1502 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1506 return StartAnalysis(type, tree, nentries, firstentry);
1509 //______________________________________________________________________________
1510 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1512 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1513 // MIX. Process nentries starting from firstentry
1515 // Backup current directory and make sure gDirectory points to gROOT
1516 TDirectory *cdir = gDirectory;
1519 Error("StartAnalysis","Analysis manager was not initialized !");
1523 if (!CheckTasks()) Fatal("StartAnalysis", "Not all needed libraries were loaded");
1524 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1525 fMaxEntries = nentries;
1527 TString anaType = type;
1529 fMode = kLocalAnalysis;
1530 Bool_t runlocalinit = kTRUE;
1531 if (anaType.Contains("file")) {
1532 runlocalinit = kFALSE;
1535 if (anaType.Contains("proof")) fMode = kProofAnalysis;
1536 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1537 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
1539 if (fMode == kGridAnalysis) {
1541 if (!anaType.Contains("terminate")) {
1542 if (!fGridHandler) {
1543 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1544 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1548 // Write analysis manager in the analysis file
1549 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1550 // run local task configuration
1551 TIter nextTask(fTasks);
1552 AliAnalysisTask *task;
1553 while ((task=(AliAnalysisTask*)nextTask())) {
1557 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1558 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1563 // Terminate grid analysis
1564 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1565 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1566 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1567 if (!fGridHandler->MergeOutputs()) {
1568 // Return if outputs could not be merged or if it alien handler
1569 // was configured for offline mode or local testing.
1574 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1575 ImportWrappers(NULL);
1581 SetEventLoop(kFALSE);
1582 // Enable event loop mode if a tree was provided
1583 if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1586 TString ttype = "TTree";
1587 if (tree && tree->IsA() == TChain::Class()) {
1588 chain = (TChain*)tree;
1589 if (!chain || !chain->GetListOfFiles()->First()) {
1590 Error("StartAnalysis", "Cannot process null or empty chain...");
1597 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1598 if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1599 // Initialize locally all tasks (happens for all modes)
1601 AliAnalysisTask *task;
1603 while ((task=(AliAnalysisTask*)next())) {
1607 if (getsysInfo) AliSysInfo::AddStamp("LocalInit_all", 0);
1611 case kLocalAnalysis:
1612 if (!tree && !fGridHandler) {
1613 TIter nextT(fTasks);
1614 // Call CreateOutputObjects for all tasks
1616 Bool_t dirStatus = TH1::AddDirectoryStatus();
1617 while ((task=(AliAnalysisTask*)nextT())) {
1618 TH1::AddDirectory(kFALSE);
1619 task->CreateOutputObjects();
1620 if (!task->CheckPostData()) {
1621 Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
1622 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
1623 ########### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ###########", task->GetName(), task->ClassName());
1625 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1629 TH1::AddDirectory(dirStatus);
1630 if (IsExternalLoop()) {
1631 Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1632 \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1639 fSelector = new AliAnalysisSelector(this);
1640 // Check if a plugin handler is used
1642 // Get the chain from the plugin
1643 TString dataType = "esdTree";
1644 if (fInputEventHandler) {
1645 dataType = fInputEventHandler->GetDataType();
1649 chain = fGridHandler->GetChainForTestMode(dataType);
1651 Error("StartAnalysis", "No chain for test mode. Aborting.");
1654 cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1655 retv = chain->Process(fSelector, "", nentries, firstentry);
1658 // Run tree-based analysis via AliAnalysisSelector
1659 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1660 retv = tree->Process(fSelector, "", nentries, firstentry);
1662 case kProofAnalysis:
1664 // Check if the plugin is used
1666 return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1668 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1669 Error("StartAnalysis", "No PROOF!!! Exiting.");
1673 line = Form("gProof->AddInput((TObject*)%p);", this);
1674 gROOT->ProcessLine(line);
1677 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1678 retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1680 Error("StartAnalysis", "No chain!!! Exiting.");
1687 if (!anaType.Contains("terminate")) {
1688 if (!fGridHandler) {
1689 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1690 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1694 // Write analysis manager in the analysis file
1695 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1696 // Start the analysis via the handler
1697 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1698 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1703 // Terminate grid analysis
1704 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1705 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1706 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1707 if (!fGridHandler->MergeOutputs()) {
1708 // Return if outputs could not be merged or if it alien handler
1709 // was configured for offline mode or local testing.
1714 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1715 ImportWrappers(NULL);
1719 case kMixingAnalysis:
1720 // Run event mixing analysis
1722 Error("StartAnalysis", "Cannot run event mixing without event pool");
1726 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1727 fSelector = new AliAnalysisSelector(this);
1728 while ((chain=fEventPool->GetNextChain())) {
1730 // Call NotifyBinChange for all tasks
1731 while ((task=(AliAnalysisTask*)next()))
1732 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1733 retv = chain->Process(fSelector);
1735 Error("StartAnalysis", "Mixing analysis failed");
1740 PackOutput(fSelector->GetOutputList());
1747 //______________________________________________________________________________
1748 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1750 // Start analysis for this manager on a given dataset. Analysis task can be:
1751 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1753 Error("StartAnalysis","Analysis manager was not initialized !");
1757 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1758 TString anaType = type;
1760 if (!anaType.Contains("proof")) {
1761 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1764 fMode = kProofAnalysis;
1766 SetEventLoop(kTRUE);
1767 // Set the dataset flag
1768 TObject::SetBit(kUseDataSet);
1771 // Start proof analysis using the grid handler
1772 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1773 Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
1776 // Check if the plugin is in test mode
1777 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
1778 dataset = "test_collection";
1780 dataset = fGridHandler->GetProofDataSet();
1784 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1785 Error("StartAnalysis", "No PROOF!!! Exiting.");
1789 // Initialize locally all tasks
1791 AliAnalysisTask *task;
1792 while ((task=(AliAnalysisTask*)next())) {
1796 line = Form("gProof->AddInput((TObject*)%p);", this);
1797 gROOT->ProcessLine(line);
1799 line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1800 dataset, nentries, firstentry);
1801 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1802 retv = (Long_t)gROOT->ProcessLine(line);
1806 //______________________________________________________________________________
1807 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
1809 // Opens according the option the file specified by cont->GetFileName() and changes
1810 // current directory to cont->GetFolderName(). If the file was already opened, it
1811 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
1812 // be optionally ignored.
1813 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1814 TString filename = cont->GetFileName();
1816 if (filename.IsNull()) {
1817 ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
1820 if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
1822 f = mgr->OpenProofFile(cont,option);
1824 // Check first if the file is already opened
1825 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1827 // Check if option "UPDATE" was preserved
1828 TString opt(option);
1830 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1831 ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1833 f = TFile::Open(filename, option);
1836 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1840 // Check for a folder request
1841 TString dir = cont->GetFolderName();
1842 if (!dir.IsNull()) {
1843 if (!f->GetDirectory(dir)) f->mkdir(dir);
1848 ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
1849 cont->SetFile(NULL);
1853 //______________________________________________________________________________
1854 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
1856 // Opens a special output file used in PROOF.
1858 TString filename = cont->GetFileName();
1859 if (cont == fCommonOutput) {
1860 if (fOutputEventHandler) {
1861 if (strlen(extaod)) filename = extaod;
1862 filename = fOutputEventHandler->GetOutputFileName();
1864 else Fatal("OpenProofFile","No output container. Exiting.");
1867 if (fMode!=kProofAnalysis || !fSelector) {
1868 Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
1871 if (fSpecialOutputLocation.Length()) {
1872 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1874 // Check if option "UPDATE" was preserved
1875 TString opt(option);
1877 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1878 ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1880 f = new TFile(filename, option);
1882 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1886 // Check for a folder request
1887 TString dir = cont->GetFolderName();
1889 if (!f->GetDirectory(dir)) f->mkdir(dir);
1894 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1895 cont->SetFile(NULL);
1898 // Check if there is already a proof output file in the output list
1899 TObject *pof = fSelector->GetOutputList()->FindObject(filename);
1901 // Get the actual file
1902 line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
1903 filename = (const char*)gROOT->ProcessLine(line);
1905 printf("File: %s already booked via TProofOutputFile\n", filename.Data());
1907 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1909 Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
1912 // Check if option "UPDATE" was preserved
1913 TString opt(option);
1915 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1916 Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
1918 if (cont->IsRegisterDataset()) {
1919 TString dsetName = filename;
1920 dsetName.ReplaceAll(".root", cont->GetTitle());
1921 dsetName.ReplaceAll(":","_");
1922 if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
1923 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
1925 if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
1926 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
1928 if (fDebug > 1) printf("=== %s\n", line.Data());
1929 gROOT->ProcessLine(line);
1930 line = Form("pf->OpenFile(\"%s\");", option);
1931 gROOT->ProcessLine(line);
1934 gROOT->ProcessLine("pf->Print()");
1935 printf(" == proof file name: %s", f->GetName());
1937 // Add to proof output list
1938 line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
1939 if (fDebug > 1) printf("=== %s\n", line.Data());
1940 gROOT->ProcessLine(line);
1942 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1946 // Check for a folder request
1947 TString dir = cont->GetFolderName();
1948 if (!dir.IsNull()) {
1949 if (!f->GetDirectory(dir)) f->mkdir(dir);
1954 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1955 cont->SetFile(NULL);
1959 //______________________________________________________________________________
1960 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1962 // Execute analysis.
1963 static Long64_t nentries = 0;
1964 static TTree *lastTree = 0;
1965 static TStopwatch *timer = new TStopwatch();
1966 if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
1968 if (fTree && (fTree != lastTree)) {
1969 nentries += fTree->GetEntries();
1972 if (!fNcalls) timer->Start();
1973 if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, TMath::Min(fMaxEntries,nentries), timer, kFALSE);
1976 TDirectory *cdir = gDirectory;
1977 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1978 if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
1980 Error("ExecAnalysis", "Analysis manager was not initialized !");
1985 AliAnalysisTask *task;
1986 // Check if the top tree is active.
1988 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1989 AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
1991 // De-activate all tasks
1992 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1993 AliAnalysisDataContainer *cont = fCommonInput;
1994 if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
1996 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
2000 cont->SetData(fTree); // This will notify all consumers
2001 Long64_t entry = fTree->GetTree()->GetReadEntry();
2003 // Call BeginEvent() for optional input/output and MC services
2004 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
2005 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
2006 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
2008 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2009 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2011 // Execute the tasks
2012 // TIter next1(cont->GetConsumers());
2013 TIter next1(fTopTasks);
2015 while ((task=(AliAnalysisTask*)next1())) {
2017 cout << " Executing task " << task->GetName() << endl;
2019 task->ExecuteTask(option);
2021 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2022 AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
2026 // Call FinishEvent() for optional output and MC services
2027 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2028 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2029 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2030 // Gather system information if requested
2031 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2032 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
2036 // The event loop is not controlled by TSelector
2038 // Call BeginEvent() for optional input/output and MC services
2039 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
2040 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
2041 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
2043 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2044 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2045 TIter next2(fTopTasks);
2046 while ((task=(AliAnalysisTask*)next2())) {
2047 task->SetActive(kTRUE);
2049 cout << " Executing task " << task->GetName() << endl;
2051 task->ExecuteTask(option);
2055 // Call FinishEvent() for optional output and MC services
2056 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2057 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2058 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2059 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2060 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2064 //______________________________________________________________________________
2065 Bool_t AliAnalysisManager::IsPipe(std::ostream &out)
2067 // Check if the stdout is connected to a pipe (C.Holm)
2068 Bool_t ispipe = kFALSE;
2069 out.seekp(0, std::ios_base::cur);
2072 if (errno == ESPIPE) ispipe = kTRUE;
2077 //______________________________________________________________________________
2078 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2080 // Set the input event handler and create a container for it.
2081 fInputEventHandler = handler;
2082 fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2085 //______________________________________________________________________________
2086 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2088 // Set the input event handler and create a container for it.
2089 fOutputEventHandler = handler;
2090 fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2091 fCommonOutput->SetSpecialOutput();
2094 //______________________________________________________________________________
2095 void AliAnalysisManager::SetDebugLevel(UInt_t level)
2097 // Set verbosity of the analysis manager. If the progress bar is used, the call is ignored
2098 if (TObject::TestBit(kUseProgressBar)) {
2099 Info("SetDebugLevel","Ignored. Disable the progress bar first.");
2105 //______________________________________________________________________________
2106 void AliAnalysisManager::SetUseProgressBar(Bool_t flag, Int_t freq)
2108 // Enable a text mode progress bar. Resets debug level to 0.
2109 Info("SetUseProgressBar", "Progress bar enabled, updated every %d events.\n ### NOTE: Debug level reset to 0 ###", freq);
2110 TObject::SetBit(kUseProgressBar,flag);
2111 fPBUpdateFreq = freq;
2115 //______________________________________________________________________________
2116 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2118 // This method is used externally to register output files which are not
2119 // connected to any output container, so that the manager can properly register,
2120 // retrieve or merge them when running in distributed mode. The file names are
2121 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2122 if (fExtraFiles.Contains(fname)) return;
2123 if (fExtraFiles.Length()) fExtraFiles += " ";
2124 fExtraFiles += fname;
2127 //______________________________________________________________________________
2128 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2130 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2134 TObject *pof = source->FindObject(filename);
2135 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2136 Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2139 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2140 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2141 TString clientUrl(chUrl);
2142 TString fullPath_str(fullPath);
2143 if (clientUrl.Contains("localhost")){
2144 TObjArray* array = fullPath_str.Tokenize ( "//" );
2145 TObjString *strobj = ( TObjString *)array->At(1);
2146 TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2147 TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2148 fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2149 fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2150 if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2154 else if (clientUrl.Contains("__lite__")) {
2155 // Special case for ProofLite environement - get file info and copy.
2156 gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2157 fullPath_str = Form("%s/%s", tmp, fullPath);
2160 Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2161 Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename);
2163 Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2167 //______________________________________________________________________________
2168 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2170 // Fill analysis type in the provided string.
2172 case kLocalAnalysis:
2175 case kProofAnalysis:
2181 case kMixingAnalysis:
2186 //______________________________________________________________________________
2187 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2189 // Validate all output files.
2190 TIter next(fOutputs);
2191 AliAnalysisDataContainer *output;
2192 TDirectory *cdir = gDirectory;
2193 TString openedFiles;
2194 while ((output=(AliAnalysisDataContainer*)next())) {
2195 if (output->IsRegisterDataset()) continue;
2196 TString filename = output->GetFileName();
2197 if (filename == "default") {
2198 if (!fOutputEventHandler) continue;
2199 filename = fOutputEventHandler->GetOutputFileName();
2200 // Main AOD may not be there
2201 if (gSystem->AccessPathName(filename)) continue;
2203 // Check if the file is closed
2204 if (openedFiles.Contains(filename)) continue;;
2205 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2207 Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2208 // Clear file list to release object ownership to user.
2212 file = TFile::Open(filename);
2213 if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2214 Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2219 openedFiles += filename;
2226 //______________________________________________________________________________
2227 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2229 // Implements a nice text mode progress bar.
2230 static Long64_t icount = 0;
2231 static TString oname;
2232 static TString nname;
2233 static Long64_t ocurrent = 0;
2234 static Long64_t osize = 0;
2235 static Int_t oseconds = 0;
2236 static TStopwatch *owatch = 0;
2237 static Bool_t oneoftwo = kFALSE;
2238 static Int_t nrefresh = 0;
2239 static Int_t nchecks = 0;
2240 static char lastChar = 0;
2241 const char symbol[4] = {'-','\\','|','/'};
2243 if (!lastChar) lastChar = (IsPipe(std::cerr))?'\r':'\n';
2249 ocurrent = TMath::Abs(current);
2250 osize = TMath::Abs(size);
2251 if (ocurrent > osize) ocurrent=osize;
2256 if ((current % fPBUpdateFreq) != 0) return;
2258 char progress[11] = " ";
2259 Int_t ichar = icount%4;
2264 if (owatch && !last) {
2266 time = owatch->RealTime();
2267 seconds = int(time) % 60;
2268 minutes = (int(time) / 60) % 60;
2269 hours = (int(time) / 60 / 60);
2271 if (oseconds==seconds) {
2275 oneoftwo = !oneoftwo;
2279 if (refresh && oneoftwo) {
2281 if (nchecks <= 0) nchecks = nrefresh+1;
2282 Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2283 oname = Form(" == %d%% ==", pctdone);
2285 Double_t percent = 100.0*ocurrent/osize;
2286 Int_t nchar = Int_t(percent/10);
2287 if (nchar>10) nchar=10;
2289 for (i=0; i<nchar; i++) progress[i] = '=';
2290 progress[nchar] = symbol[ichar];
2291 for (i=nchar+1; i<10; i++) progress[i] = ' ';
2292 progress[10] = '\0';
2295 if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2296 else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2297 else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2299 Int_t full = Int_t(ocurrent > 0 ?
2300 time * (float(osize)/ocurrent) + .5 :
2302 Int_t remain = Int_t(full - time);
2303 Int_t rsec = remain % 60;
2304 Int_t rmin = (remain / 60) % 60;
2305 Int_t rhour = (remain / 60 / 60);
2306 fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d ETA %.2d:%.2d:%.2d%c",
2307 percent, hours, minutes, seconds, rhour, rmin, rsec, lastChar);
2309 else fprintf(stderr, "[%6.2f %%]%c", percent, lastChar);
2310 if (refresh && oneoftwo) oname = nname;
2311 if (owatch) owatch->Continue();
2320 fprintf(stderr, "\n");
2324 //______________________________________________________________________________
2325 void AliAnalysisManager::DoLoadBranch(const char *name)
2327 // Get tree and load branch if needed.
2328 static Long64_t crtEntry = -100;
2330 if (fAutoBranchHandling || !fTree)
2333 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2335 br = fTree->GetBranch(name);
2337 Error("DoLoadBranch", "Could not find branch %s",name);
2342 if (br->GetReadEntry()==fCurrentEntry) return;
2343 Int_t ret = br->GetEntry(GetCurrentEntry());
2345 Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2346 if (crtEntry != fCurrentEntry) {
2347 CountEvent(1,0,1,0);
2348 crtEntry = fCurrentEntry;
2351 if (crtEntry != fCurrentEntry) {
2352 CountEvent(1,1,0,0);
2353 crtEntry = fCurrentEntry;
2358 //______________________________________________________________________________
2359 void AliAnalysisManager::AddStatisticsTask(UInt_t offlineMask)
2361 // Add the statistics task to the manager.
2363 Info("AddStatisticsTask", "Already added");
2366 TString line = Form("AliAnalysisTaskStat::AddToManager(%u);", offlineMask);
2367 gROOT->ProcessLine(line);
2370 //______________________________________________________________________________
2371 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2373 // Bookkeep current event;
2374 if (!fStatistics) return;
2375 fStatistics->AddInput(ninput);
2376 fStatistics->AddProcessed(nprocessed);
2377 fStatistics->AddFailed(nfailed);
2378 fStatistics->AddAccepted(naccepted);
2381 //______________________________________________________________________________
2382 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2384 // Add a line in the statistics message. If available, the statistics message is written
2385 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2387 if (!strlen(line)) return;
2388 if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2389 fStatisticsMsg += line;
2392 //______________________________________________________________________________
2393 void AliAnalysisManager::WriteStatisticsMsg(Int_t)
2395 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2396 static Bool_t done = kFALSE;
2399 if (!fStatistics) return;
2401 AddStatisticsMsg(Form("Number of input events: %lld",fStatistics->GetNinput()));
2402 AddStatisticsMsg(Form("Number of processed events: %lld",fStatistics->GetNprocessed()));
2403 AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2404 AddStatisticsMsg(Form("Number of accepted events for mask %s: %lld", AliAnalysisStatistics::GetMaskAsString(fStatistics->GetOfflineMask()), fStatistics->GetNaccepted()));
2405 out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2406 fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2407 fStatistics->GetNaccepted()), ios::out);
2408 out << fStatisticsMsg << endl;
2412 //______________________________________________________________________________
2413 const char* AliAnalysisManager::GetOADBPath()
2415 // returns the path of the OADB
2416 // this static function just depends on environment variables
2418 static TString oadbPath;
2420 if (gSystem->Getenv("OADB_PATH"))
2421 oadbPath = gSystem->Getenv("OADB_PATH");
2422 else if (gSystem->Getenv("ALICE_ROOT"))
2423 oadbPath.Form("%s/OADB", gSystem->Getenv("ALICE_ROOT"));
2425 ::Fatal("AliAnalysisManager::GetOADBPath", "Cannot figure out AODB path. Define ALICE_ROOT or OADB_PATH!");