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 static Bool_t isCalled = kFALSE;
349 Bool_t init = kFALSE;
350 Bool_t initOK = kTRUE;
352 TDirectory *curdir = gDirectory;
353 // Call SlaveBegin only once in case of mixing
354 if (isCalled && fMode==kMixingAnalysis) return;
356 // Call Init of EventHandler
357 if (fOutputEventHandler) {
358 if (fMode == kProofAnalysis) {
359 // Merging AOD's in PROOF via TProofOutputFile
360 if (fDebug > 1) printf(" Initializing AOD output file %s...\n", fOutputEventHandler->GetOutputFileName());
361 init = fOutputEventHandler->Init("proof");
362 if (!init) msg = "Failed to initialize output handler on worker";
364 init = fOutputEventHandler->Init("local");
365 if (!init) msg = "Failed to initialize output handler";
368 if (!fSelector) Error("SlaveBegin", "Selector not set");
369 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
372 if (fInputEventHandler) {
373 fInputEventHandler->SetInputTree(tree);
374 if (fMode == kProofAnalysis) {
375 init = fInputEventHandler->Init("proof");
376 if (!init) msg = "Failed to initialize input handler on worker";
378 init = fInputEventHandler->Init("local");
379 if (!init) msg = "Failed to initialize input handler";
382 if (!fSelector) Error("SlaveBegin", "Selector not set");
383 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
386 if (fMCtruthEventHandler) {
387 if (fMode == kProofAnalysis) {
388 init = fMCtruthEventHandler->Init("proof");
389 if (!init) msg = "Failed to initialize MC handler on worker";
391 init = fMCtruthEventHandler->Init("local");
392 if (!init) msg = "Failed to initialize MC handler";
395 if (!fSelector) Error("SlaveBegin", "Selector not set");
396 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
398 if (curdir) curdir->cd();
402 AliAnalysisTask *task;
403 // Call CreateOutputObjects for all tasks
404 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
405 Bool_t dirStatus = TH1::AddDirectoryStatus();
407 while ((task=(AliAnalysisTask*)next())) {
409 // Start with memory as current dir and make sure by default histograms do not get attached to files.
410 TH1::AddDirectory(kFALSE);
411 task->CreateOutputObjects();
412 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
415 TH1::AddDirectory(dirStatus);
416 if (curdir) curdir->cd();
417 if (fDebug > 1) printf("<-AliAnalysisManager::SlaveBegin()\n");
420 //______________________________________________________________________________
421 Bool_t AliAnalysisManager::Notify()
423 // The Notify() function is called when a new file is opened. This
424 // can be either for a new TTree in a TChain or when when a new TTree
425 // is started when using PROOF. It is normaly not necessary to make changes
426 // to the generated code, but the routine can be extended by the
427 // user if needed. The return value is currently not used.
428 if (!fTree) return kFALSE;
430 fTable.Clear("nodelete"); // clearing the hash table may not be needed -> C.L.
431 if (fMode == kProofAnalysis) fIsRemote = kTRUE;
433 TFile *curfile = fTree->GetCurrentFile();
435 Error("Notify","No current file");
439 if (fDebug > 1) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
440 Int_t run = AliAnalysisManager::GetRunFromAlienPath(curfile->GetName());
441 if (run) SetRunFromPath(run);
442 if (fDebug > 1) printf(" ### run found from path: %d\n", run);
444 AliAnalysisTask *task;
446 // Call Notify of the event handlers
447 if (fInputEventHandler) {
448 fInputEventHandler->Notify(curfile->GetName());
451 if (fOutputEventHandler) {
452 fOutputEventHandler->Notify(curfile->GetName());
455 if (fMCtruthEventHandler) {
456 fMCtruthEventHandler->Notify(curfile->GetName());
459 // Call Notify for all tasks
460 while ((task=(AliAnalysisTask*)next()))
463 if (fDebug > 1) printf("<-AliAnalysisManager::Notify()\n");
467 //______________________________________________________________________________
468 Bool_t AliAnalysisManager::Process(Long64_t entry)
470 // The Process() function is called for each entry in the tree (or possibly
471 // keyed object in the case of PROOF) to be processed. The entry argument
472 // specifies which entry in the currently loaded tree is to be processed.
473 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
474 // to read either all or the required parts of the data. When processing
475 // keyed objects with PROOF, the object is already loaded and is available
476 // via the fObject pointer.
478 // This function should contain the "body" of the analysis. It can contain
479 // simple or elaborate selection criteria, run algorithms on the data
480 // of the event and typically fill histograms.
482 // WARNING when a selector is used with a TChain, you must use
483 // the pointer to the current TTree to call GetEntry(entry).
484 // The entry is always the local entry number in the current tree.
485 // Assuming that fChain is the pointer to the TChain being processed,
486 // use fChain->GetTree()->GetEntry(entry).
487 if (fDebug > 1) printf("->AliAnalysisManager::Process(%lld)\n", entry);
489 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
490 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
491 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
495 if (fInputEventHandler) fInputEventHandler ->GetEntry();
498 if (fDebug > 1) printf("<-AliAnalysisManager::Process()\n");
502 //______________________________________________________________________________
503 void AliAnalysisManager::PackOutput(TList *target)
505 // Pack all output data containers in the output list. Called at SlaveTerminate
506 // stage in PROOF case for each slave.
507 if (fDebug > 1) printf("->AliAnalysisManager::PackOutput()\n");
509 Error("PackOutput", "No target. Exiting.");
512 TDirectory *cdir = gDirectory;
514 if (fInputEventHandler) fInputEventHandler ->Terminate();
515 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
516 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
519 // Call FinishTaskOutput() for each event loop task (not called for
520 // post-event loop tasks - use Terminate() fo those)
521 TIter nexttask(fTasks);
522 AliAnalysisTask *task;
523 while ((task=(AliAnalysisTask*)nexttask())) {
524 if (!task->IsPostEventLoop()) {
525 if (fDebug > 1) printf("->FinishTaskOutput: task %s\n", task->GetName());
526 task->FinishTaskOutput();
528 if (fDebug > 1) printf("<-FinishTaskOutput: task %s\n", task->GetName());
531 // Write statistics message on the workers.
532 if (fStatistics) WriteStatisticsMsg(fNcalls);
534 if (fMode == kProofAnalysis) {
535 TIter next(fOutputs);
536 AliAnalysisDataContainer *output;
537 Bool_t isManagedByHandler = kFALSE;
540 while ((output=(AliAnalysisDataContainer*)next())) {
541 // Do not consider outputs of post event loop tasks
542 isManagedByHandler = kFALSE;
543 if (output->GetProducer() && output->GetProducer()->IsPostEventLoop()) continue;
544 const char *filename = output->GetFileName();
545 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
546 isManagedByHandler = kTRUE;
547 printf("#### Handler output. Extra: %s\n", fExtraFiles.Data());
548 filename = fOutputEventHandler->GetOutputFileName();
550 // Check if data was posted to this container. If not, issue an error.
551 if (!output->GetData() && !isManagedByHandler) {
552 Error("PackOutput", "No data for output container %s. Forgot to PostData ?", output->GetName());
555 if (!output->IsSpecialOutput()) {
557 if (strlen(filename) && !isManagedByHandler) {
558 // Backup current folder
559 TDirectory *opwd = gDirectory;
560 // File resident outputs.
561 // Check first if the file exists.
562 TString openoption = "RECREATE";
563 Bool_t firsttime = kTRUE;
564 if (filestmp.FindObject(output->GetFileName())) {
567 filestmp.Add(new TNamed(output->GetFileName(),""));
569 if (!gSystem->AccessPathName(output->GetFileName()) && !firsttime) openoption = "UPDATE";
570 // TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
571 // Save data to file, then close.
572 if (output->GetData()->InheritsFrom(TCollection::Class())) {
573 // If data is a collection, we set the name of the collection
574 // as the one of the container and we save as a single key.
575 TCollection *coll = (TCollection*)output->GetData();
576 coll->SetName(output->GetName());
577 // coll->Write(output->GetName(), TObject::kSingleKey);
579 if (output->GetData()->InheritsFrom(TTree::Class())) {
580 TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
581 // Save data to file, then close.
582 TTree *tree = (TTree*)output->GetData();
583 // Check if tree is in memory
584 if (tree->GetDirectory()==gROOT) tree->SetDirectory(gDirectory);
588 // output->GetData()->Write();
591 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
593 // printf(" file %s listing content:\n", filename);
596 // Clear file list to release object ownership to user.
599 output->SetFile(NULL);
600 // Restore current directory
601 if (opwd) opwd->cd();
603 // Memory-resident outputs
604 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
606 AliAnalysisDataWrapper *wrap = 0;
607 if (isManagedByHandler) {
608 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
609 wrap->SetName(output->GetName());
611 else wrap =output->ExportData();
612 // Output wrappers must NOT delete data after merging - the user owns them
613 wrap->SetDeleteData(kFALSE);
616 // Special outputs. The file must be opened and connected to the container.
617 TDirectory *opwd = gDirectory;
618 TFile *file = output->GetFile();
620 AliAnalysisTask *producer = output->GetProducer();
622 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
623 output->GetFileName(), output->GetName(), producer->ClassName());
626 TString outFilename = file->GetName();
627 if (fDebug > 1) printf("PackOutput %s: special output\n", output->GetName());
628 if (isManagedByHandler) {
629 // Terminate IO for files managed by the output handler
630 // file->Write() moved to AOD handler (A.G. 11.01.10)
631 // if (file) file->Write();
632 if (file && fDebug > 2) {
633 printf(" handled file %s listing content:\n", file->GetName());
636 fOutputEventHandler->TerminateIO();
639 // Release object ownership to users after writing data to file
640 if (output->GetData()->InheritsFrom(TCollection::Class())) {
641 // If data is a collection, we set the name of the collection
642 // as the one of the container and we save as a single key.
643 TCollection *coll = (TCollection*)output->GetData();
644 coll->SetName(output->GetName());
645 coll->Write(output->GetName(), TObject::kSingleKey);
647 if (output->GetData()->InheritsFrom(TTree::Class())) {
648 TTree *tree = (TTree*)output->GetData();
649 tree->SetDirectory(file);
652 output->GetData()->Write();
656 printf(" file %s listing content:\n", output->GetFileName());
659 // Clear file list to release object ownership to user.
662 output->SetFile(NULL);
664 // Restore current directory
665 if (opwd) opwd->cd();
666 // Check if a special output location was provided or the output files have to be merged
667 if (strlen(fSpecialOutputLocation.Data())) {
668 TString remote = fSpecialOutputLocation;
670 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
671 if (remote.BeginsWith("alien:")) {
672 gROOT->ProcessLine("TGrid::Connect(\"alien:\", gProofServ->GetUser());");
673 remote += outFilename;
674 remote.ReplaceAll(".root", Form("_%d.root", gid));
676 remote += Form("%s_%d_", gSystem->HostName(), gid);
677 remote += outFilename;
680 Info("PackOutput", "Output file for container %s to be copied \n at: %s. No merging.",
681 output->GetName(), remote.Data());
682 TFile::Cp ( outFilename.Data(), remote.Data() );
683 // Copy extra outputs
684 if (fExtraFiles.Length() && isManagedByHandler) {
685 TObjArray *arr = fExtraFiles.Tokenize(" ");
687 TIter nextfilename(arr);
688 while ((os=(TObjString*)nextfilename())) {
689 outFilename = os->GetString();
690 remote = fSpecialOutputLocation;
692 if (remote.BeginsWith("alien://")) {
693 remote += outFilename;
694 remote.ReplaceAll(".root", Form("_%d.root", gid));
696 remote += Form("%s_%d_", gSystem->HostName(), gid);
697 remote += outFilename;
700 Info("PackOutput", "Extra AOD file %s to be copied \n at: %s. No merging.",
701 outFilename.Data(), remote.Data());
702 TFile::Cp ( outFilename.Data(), remote.Data() );
707 // No special location specified-> use TProofOutputFile as merging utility
708 // The file at this output slot must be opened in CreateOutputObjects
709 if (fDebug > 1) printf(" File for container %s to be merged via file merger...\n", output->GetName());
715 if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
718 //______________________________________________________________________________
719 void AliAnalysisManager::ImportWrappers(TList *source)
721 // Import data in output containers from wrappers coming in source.
722 if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
723 TIter next(fOutputs);
724 AliAnalysisDataContainer *cont;
725 AliAnalysisDataWrapper *wrap;
727 Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
728 TDirectory *cdir = gDirectory;
729 while ((cont=(AliAnalysisDataContainer*)next())) {
731 if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
732 if (cont->IsRegisterDataset()) continue;
733 const char *filename = cont->GetFileName();
734 Bool_t isManagedByHandler = kFALSE;
735 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
736 isManagedByHandler = kTRUE;
737 filename = fOutputEventHandler->GetOutputFileName();
739 if (cont->IsSpecialOutput() || inGrid) {
740 if (strlen(fSpecialOutputLocation.Data())) continue;
741 // Copy merged file from PROOF scratch space.
742 // In case of grid the files are already in the current directory.
744 if (isManagedByHandler && fExtraFiles.Length()) {
745 // Copy extra registered dAOD files.
746 TObjArray *arr = fExtraFiles.Tokenize(" ");
748 TIter nextfilename(arr);
749 while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
752 if (!GetFileFromWrapper(filename, source)) continue;
754 // Normally we should connect data from the copied file to the
755 // corresponding output container, but it is not obvious how to do this
756 // automatically if several objects in file...
757 TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
758 if (!f) f = TFile::Open(filename, "READ");
760 Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
764 // Cd to the directory pointed by the container
765 TString folder = cont->GetFolderName();
766 if (!folder.IsNull()) f->cd(folder);
767 // Try to fetch first an object having the container name.
768 obj = gDirectory->Get(cont->GetName());
770 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",
771 cont->GetType()->GetName(), cont->GetName(), filename, cont->GetFolderName(), cont->GetName());
774 wrap = new AliAnalysisDataWrapper(obj);
775 wrap->SetDeleteData(kFALSE);
777 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
779 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
784 printf(" Importing data for container %s\n", cont->GetName());
785 if (strlen(filename)) printf(" -> file %s\n", filename);
788 cont->ImportData(wrap);
790 if (cdir) cdir->cd();
791 if (fDebug > 1) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
794 //______________________________________________________________________________
795 void AliAnalysisManager::UnpackOutput(TList *source)
797 // Called by AliAnalysisSelector::Terminate only on the client.
798 if (fDebug > 1) printf("->AliAnalysisManager::UnpackOutput()\n");
800 Error("UnpackOutput", "No target. Exiting.");
803 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
805 if (fMode == kProofAnalysis) ImportWrappers(source);
807 TIter next(fOutputs);
808 AliAnalysisDataContainer *output;
809 while ((output=(AliAnalysisDataContainer*)next())) {
810 if (!output->GetData()) continue;
811 // Check if there are client tasks that run post event loop
812 if (output->HasConsumers()) {
813 // Disable event loop semaphore
814 output->SetPostEventLoop(kTRUE);
815 TObjArray *list = output->GetConsumers();
816 Int_t ncons = list->GetEntriesFast();
817 for (Int_t i=0; i<ncons; i++) {
818 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
819 task->CheckNotify(kTRUE);
820 // If task is active, execute it
821 if (task->IsPostEventLoop() && task->IsActive()) {
822 if (fDebug > 1) printf("== Executing post event loop task %s\n", task->GetName());
828 if (fDebug > 1) printf("<-AliAnalysisManager::UnpackOutput()\n");
831 //______________________________________________________________________________
832 void AliAnalysisManager::Terminate()
834 // The Terminate() function is the last function to be called during
835 // a query. It always runs on the client, it can be used to present
836 // the results graphically.
837 if (fDebug > 1) printf("->AliAnalysisManager::Terminate()\n");
838 TDirectory *cdir = gDirectory;
840 AliAnalysisTask *task;
841 AliAnalysisDataContainer *output;
844 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
845 // Call Terminate() for tasks
847 while (!IsSkipTerminate() && (task=(AliAnalysisTask*)next())) {
848 // Save all the canvases produced by the Terminate
849 TString pictname = Form("%s_%s", task->GetName(), task->ClassName());
853 AliSysInfo::AddStamp(Form("%s_TERMINATE",task->ClassName()),0, itask, 2);
855 if (TObject::TestBit(kSaveCanvases)) {
856 if (!gROOT->IsBatch()) {
857 if (fDebug>1) printf("Waiting 5 sec for %s::Terminate() to finish drawing ...\n", task->ClassName());
859 while (timer.CpuTime()<5) {
861 gSystem->ProcessEvents();
864 Int_t iend = gROOT->GetListOfCanvases()->GetEntries();
865 if (iend==0) continue;
867 for (Int_t ipict=0; ipict<iend; ipict++) {
868 canvas = (TCanvas*)gROOT->GetListOfCanvases()->At(ipict);
869 if (!canvas) continue;
870 canvas->SaveAs(Form("%s_%02d.gif", pictname.Data(),ipict));
872 gROOT->GetListOfCanvases()->Delete();
876 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
877 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
878 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
880 TObjArray *allOutputs = new TObjArray();
882 for (icont=0; icont<fOutputs->GetEntriesFast(); icont++) allOutputs->Add(fOutputs->At(icont));
883 if (!IsSkipTerminate())
884 for (icont=0; icont<fParamCont->GetEntriesFast(); icont++) allOutputs->Add(fParamCont->At(icont));
885 TIter next1(allOutputs);
886 TString handlerFile = "";
887 TString extraOutputs = "";
888 if (fOutputEventHandler) {
889 handlerFile = fOutputEventHandler->GetOutputFileName();
890 extraOutputs = fOutputEventHandler->GetExtraOutputs();
894 while ((output=(AliAnalysisDataContainer*)next1())) {
895 // Special outputs or grid files have the files already closed and written.
897 if (fMode == kGridAnalysis && icont<=fOutputs->GetEntriesFast()) continue;
898 if (fMode == kProofAnalysis) {
899 if (output->IsSpecialOutput() || output->IsRegisterDataset()) continue;
901 const char *filename = output->GetFileName();
902 TString openoption = "RECREATE";
903 if (!(strcmp(filename, "default"))) continue;
904 if (!strlen(filename)) continue;
905 if (!output->GetData()) continue;
906 TDirectory *opwd = gDirectory;
907 TFile *file = output->GetFile();
908 if (!file) file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
910 //if (handlerFile == filename && !gSystem->AccessPathName(filename)) openoption = "UPDATE";
911 Bool_t firsttime = kTRUE;
912 if (filestmp.FindObject(filename) || extraOutputs.Contains(filename)) {
915 filestmp.Add(new TNamed(filename,""));
917 if (!gSystem->AccessPathName(filename) && !firsttime) openoption = "UPDATE";
918 if (fDebug>1) printf("Opening file: %s option=%s\n",filename, openoption.Data());
919 file = new TFile(filename, openoption);
921 if (fDebug>1) printf("File <%s> already opened with option: <%s> \n", filename, file->GetOption());
922 openoption = file->GetOption();
923 if (openoption == "READ") {
924 if (fDebug>1) printf("...reopening in UPDATE mode\n");
925 file->ReOpen("UPDATE");
928 if (file->IsZombie()) {
929 Error("Terminate", "Cannot open output file %s", filename);
932 output->SetFile(file);
934 // Check for a folder request
935 TString dir = output->GetFolderName();
937 if (!file->GetDirectory(dir)) file->mkdir(dir);
940 if (fDebug > 1) printf("...writing container %s to file %s:%s\n", output->GetName(), file->GetName(), output->GetFolderName());
941 if (output->GetData()->InheritsFrom(TCollection::Class())) {
942 // If data is a collection, we set the name of the collection
943 // as the one of the container and we save as a single key.
944 TCollection *coll = (TCollection*)output->GetData();
945 coll->SetName(output->GetName());
946 coll->Write(output->GetName(), TObject::kSingleKey);
948 if (output->GetData()->InheritsFrom(TTree::Class())) {
949 TTree *tree = (TTree*)output->GetData();
950 tree->SetDirectory(gDirectory);
953 output->GetData()->Write();
956 if (opwd) opwd->cd();
960 while ((output=(AliAnalysisDataContainer*)next1())) {
961 // Close all files at output
962 TDirectory *opwd = gDirectory;
963 if (output->GetFile()) {
964 // Clear file list to release object ownership to user.
965 // output->GetFile()->Clear();
966 output->GetFile()->Close();
967 output->SetFile(NULL);
968 // Copy merged outputs in alien if requested
969 if (fSpecialOutputLocation.Length() &&
970 fSpecialOutputLocation.BeginsWith("alien://")) {
971 Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data());
972 TFile::Cp(output->GetFile()->GetName(),
973 Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
976 if (opwd) opwd->cd();
979 //Write statistics information on the client
980 if (fStatistics) WriteStatisticsMsg(fNcalls);
982 TDirectory *crtdir = gDirectory;
983 TFile f("syswatch.root", "RECREATE");
987 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
988 tree->SetName("syswatch");
989 tree->SetMarkerStyle(kCircle);
990 tree->SetMarkerColor(kBlue);
991 tree->SetMarkerSize(0.5);
992 if (!gROOT->IsBatch()) {
993 tree->SetAlias("event", "id0");
994 tree->SetAlias("task", "id1");
995 tree->SetAlias("stage", "id2");
996 // Already defined aliases
997 // tree->SetAlias("deltaT","stampSec-stampOldSec");
998 // tree->SetAlias("T","stampSec-first");
999 // tree->SetAlias("deltaVM","(pI.fMemVirtual-pIOld.fMemVirtual)");
1000 // tree->SetAlias("VM","pI.fMemVirtual");
1001 TCanvas *canvas = new TCanvas("SysInfo","SysInfo",10,10,1200,800);
1002 Int_t npads = 1 /*COO plot for all tasks*/ +
1003 fTopTasks->GetEntries() /*Exec plot per task*/ +
1004 1 /*Terminate plot for all tasks*/ +
1007 Int_t iopt = (Int_t)TMath::Sqrt((Double_t)npads);
1008 if (npads<iopt*(iopt+1))
1009 canvas->Divide(iopt, iopt+1, 0.01, 0.01);
1011 canvas->Divide(iopt+1, iopt+1, 0.01, 0.01);
1013 // draw the plot of deltaVM for Exec for each task
1014 for (itask=0; itask<fTopTasks->GetEntriesFast(); itask++) {
1015 task = (AliAnalysisTask*)fTopTasks->At(itask);
1017 cut = Form("task==%d && stage==1", itask);
1018 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1019 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1021 hist->SetTitle(Form("%s: Exec dVM[kB]/event", task->GetName()));
1022 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1025 // Draw the plot of deltaVM for CreateOutputObjects for all tasks
1027 tree->SetMarkerStyle(kFullTriangleUp);
1028 tree->SetMarkerColor(kRed);
1029 tree->SetMarkerSize(0.8);
1030 cut = "task>=0 && task<1000 && stage==0";
1031 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1032 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1034 hist->SetTitle("Memory in CreateOutputObjects()");
1035 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1036 hist->GetXaxis()->SetTitle("task");
1038 // draw the plot of deltaVM for Terminate for all tasks
1040 tree->SetMarkerStyle(kOpenSquare);
1041 tree->SetMarkerColor(kMagenta);
1042 cut = "task>=0 && task<1000 && stage==2";
1043 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1044 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1046 hist->SetTitle("Memory in Terminate()");
1047 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1048 hist->GetXaxis()->SetTitle("task");
1052 tree->SetMarkerStyle(kFullCircle);
1053 tree->SetMarkerColor(kGreen);
1054 cut = Form("task==%d && stage==1",fTopTasks->GetEntriesFast()-1);
1055 tree->Draw("VM:event",cut,"", 1234567890, 0);
1056 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1058 hist->SetTitle("Virtual memory");
1059 hist->GetYaxis()->SetTitle("VM [MB]");
1063 tree->SetMarkerStyle(kCircle);
1064 tree->SetMarkerColor(kBlue);
1065 tree->SetMarkerSize(0.5);
1070 if (crtdir) crtdir->cd();
1072 // Validate the output files
1073 if (ValidateOutputFiles() && fIsRemote && fMode!=kProofAnalysis) {
1075 out.open("outputs_valid", ios::out);
1079 if (fDebug > 1) printf("<-AliAnalysisManager::Terminate()\n");
1081 //______________________________________________________________________________
1082 void AliAnalysisManager::ProfileTask(Int_t itop, const char *option) const
1084 // Profiles the task having the itop index in the list of top (first level) tasks.
1085 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->At(itop);
1087 Error("ProfileTask", "There are only %d top tasks in the manager", fTopTasks->GetEntries());
1090 ProfileTask(task->GetName(), option);
1093 //______________________________________________________________________________
1094 void AliAnalysisManager::ProfileTask(const char *name, const char */*option*/) const
1096 // Profile a managed task after the execution of the analysis in case NSysInfo
1098 if (gSystem->AccessPathName("syswatch.root")) {
1099 Error("ProfileTask", "No file syswatch.root found in the current directory");
1102 if (gROOT->IsBatch()) return;
1103 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->FindObject(name);
1105 Error("ProfileTask", "No top task named %s known by the manager.", name);
1108 Int_t itop = fTopTasks->IndexOf(task);
1109 Int_t itask = fTasks->IndexOf(task);
1110 // Create canvas with 2 pads: first draw COO + Terminate, second Exec
1111 TDirectory *cdir = gDirectory;
1112 TFile f("syswatch.root");
1113 TTree *tree = (TTree*)f.Get("syswatch");
1115 Error("ProfileTask", "No tree named <syswatch> found in file syswatch.root");
1118 if (fDebug > 1) printf("=== Profiling task %s (class %s)\n", name, task->ClassName());
1119 TCanvas *canvas = new TCanvas(Form("profile_%d",itop),Form("Profile of task %s (class %s)",name,task->ClassName()),10,10,800,600);
1120 canvas->Divide(2, 2, 0.01, 0.01);
1124 // VM profile for COO and Terminate methods
1126 cut = Form("task==%d && (stage==0 || stage==2)",itask);
1127 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1128 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1130 hist->SetTitle("Alocated VM[kB] for COO and Terminate");
1131 hist->GetYaxis()->SetTitle("deltaVM [kB]");
1132 hist->GetXaxis()->SetTitle("method");
1134 // CPU profile per event
1136 cut = Form("task==%d && stage==1",itop);
1137 tree->Draw("deltaT:event",cut,"", 1234567890, 0);
1138 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1140 hist->SetTitle("Execution time per event");
1141 hist->GetYaxis()->SetTitle("CPU/event [s]");
1143 // VM profile for Exec
1145 cut = Form("task==%d && stage==1",itop);
1146 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1147 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1149 hist->SetTitle("Alocated VM[kB] per event");
1150 hist->GetYaxis()->SetTitle("deltaVM [kB]");
1155 if (cdir) cdir->cd();
1158 //______________________________________________________________________________
1159 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1161 // Adds a user task to the global list of tasks.
1162 if (fTasks->FindObject(task)) {
1163 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1166 task->SetActive(kFALSE);
1170 //______________________________________________________________________________
1171 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1173 // Retreive task by name.
1174 if (!fTasks) return NULL;
1175 return (AliAnalysisTask*)fTasks->FindObject(name);
1178 //______________________________________________________________________________
1179 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
1180 TClass *datatype, EAliAnalysisContType type, const char *filename)
1182 // Create a data container of a certain type. Types can be:
1183 // kExchangeContainer = 0, used to exchange data between tasks
1184 // kInputContainer = 1, used to store input data
1185 // kOutputContainer = 2, used for writing result to a file
1186 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1187 // the output object to a folder inside the output file
1188 if (fContainers->FindObject(name)) {
1189 Error("CreateContainer","A container named %s already defined !",name);
1192 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1193 fContainers->Add(cont);
1195 case kInputContainer:
1198 case kOutputContainer:
1199 fOutputs->Add(cont);
1200 if (filename && strlen(filename)) {
1201 cont->SetFileName(filename);
1202 cont->SetDataOwned(kFALSE); // data owned by the file
1205 case kParamContainer:
1206 fParamCont->Add(cont);
1207 if (filename && strlen(filename)) {
1208 cont->SetFileName(filename);
1209 cont->SetDataOwned(kFALSE); // data owned by the file
1212 case kExchangeContainer:
1218 //______________________________________________________________________________
1219 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1220 AliAnalysisDataContainer *cont)
1222 // Connect input of an existing task to a data container.
1224 Error("ConnectInput", "Task pointer is NULL");
1227 if (!fTasks->FindObject(task)) {
1229 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1231 Bool_t connected = task->ConnectInput(islot, cont);
1235 //______________________________________________________________________________
1236 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1237 AliAnalysisDataContainer *cont)
1239 // Connect output of an existing task to a data container.
1241 Error("ConnectOutput", "Task pointer is NULL");
1244 if (!fTasks->FindObject(task)) {
1246 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1248 Bool_t connected = task->ConnectOutput(islot, cont);
1252 //______________________________________________________________________________
1253 void AliAnalysisManager::CleanContainers()
1255 // Clean data from all containers that have already finished all client tasks.
1256 TIter next(fContainers);
1257 AliAnalysisDataContainer *cont;
1258 while ((cont=(AliAnalysisDataContainer *)next())) {
1259 if (cont->IsOwnedData() &&
1260 cont->IsDataReady() &&
1261 cont->ClientsExecuted()) cont->DeleteData();
1265 //______________________________________________________________________________
1266 Bool_t AliAnalysisManager::InitAnalysis()
1268 // Initialization of analysis chain of tasks. Should be called after all tasks
1269 // and data containers are properly connected
1270 // Reset flag and remove valid_outputs file if exists
1272 if (!gSystem->AccessPathName("outputs_valid"))
1273 gSystem->Unlink("outputs_valid");
1274 // Check for top tasks (depending only on input data containers)
1275 if (!fTasks->First()) {
1276 Error("InitAnalysis", "Analysis has no tasks !");
1280 AliAnalysisTask *task;
1281 AliAnalysisDataContainer *cont;
1284 Bool_t iszombie = kFALSE;
1285 Bool_t istop = kTRUE;
1287 while ((task=(AliAnalysisTask*)next())) {
1290 Int_t ninputs = task->GetNinputs();
1291 for (i=0; i<ninputs; i++) {
1292 cont = task->GetInputSlot(i)->GetContainer();
1296 fZombies->Add(task);
1300 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
1301 i, task->GetName());
1303 if (iszombie) continue;
1304 // Check if cont is an input container
1305 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1306 // Connect to parent task
1310 fTopTasks->Add(task);
1314 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1317 // Check now if there are orphan tasks
1318 for (i=0; i<ntop; i++) {
1319 task = (AliAnalysisTask*)fTopTasks->At(i);
1324 while ((task=(AliAnalysisTask*)next())) {
1325 if (!task->IsUsed()) {
1327 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1330 // Check the task hierarchy (no parent task should depend on data provided
1331 // by a daughter task)
1332 for (i=0; i<ntop; i++) {
1333 task = (AliAnalysisTask*)fTopTasks->At(i);
1334 if (task->CheckCircularDeps()) {
1335 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1340 // Check that all containers feeding post-event loop tasks are in the outputs list
1341 TIter nextcont(fContainers); // loop over all containers
1342 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1343 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1344 if (cont->HasConsumers()) {
1345 // Check if one of the consumers is post event loop
1346 TIter nextconsumer(cont->GetConsumers());
1347 while ((task=(AliAnalysisTask*)nextconsumer())) {
1348 if (task->IsPostEventLoop()) {
1349 fOutputs->Add(cont);
1356 // Check if all special output containers have a file name provided
1357 TIter nextout(fOutputs);
1358 while ((cont=(AliAnalysisDataContainer*)nextout())) {
1359 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1360 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1364 // Initialize requested branch list if needed
1365 if (!fAutoBranchHandling) {
1367 while ((task=(AliAnalysisTask*)next())) {
1368 if (!task->HasBranches()) {
1369 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\"",
1370 task->GetName(), task->ClassName());
1373 if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1374 Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1377 TString taskbranches;
1378 task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1379 if (taskbranches.IsNull()) {
1380 Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1381 task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1384 AddBranches(taskbranches);
1391 //______________________________________________________________________________
1392 void AliAnalysisManager::AddBranches(const char *branches)
1394 // Add branches to the existing fRequestedBranches.
1395 TString br(branches);
1396 TObjArray *arr = br.Tokenize(",");
1399 while ((obj=next())) {
1400 if (!fRequestedBranches.Contains(obj->GetName())) {
1401 if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1402 fRequestedBranches += obj->GetName();
1408 //______________________________________________________________________________
1409 void AliAnalysisManager::CheckBranches(Bool_t load)
1411 // The method checks the input branches to be loaded during the analysis.
1412 if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;
1413 TObjArray *arr = fRequestedBranches.Tokenize(",");
1416 while ((obj=next())) {
1417 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1419 br = fTree->GetBranch(obj->GetName());
1421 Error("CheckBranches", "Could not find branch %s",obj->GetName());
1426 if (load && br->GetReadEntry()!=GetCurrentEntry()) br->GetEntry(GetCurrentEntry());
1431 //______________________________________________________________________________
1432 void AliAnalysisManager::PrintStatus(Option_t *option) const
1434 // Print task hierarchy.
1436 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1439 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1441 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1442 TIter next(fTopTasks);
1443 AliAnalysisTask *task;
1444 while ((task=(AliAnalysisTask*)next()))
1445 task->PrintTask(option);
1447 if (!fAutoBranchHandling && !fRequestedBranches.IsNull())
1448 printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1450 TString sopt(option);
1453 if (sopt.Contains("ALL"))
1455 if ( fOutputEventHandler )
1457 cout << TString('_',78) << endl;
1458 cout << "OutputEventHandler:" << endl;
1459 fOutputEventHandler->Print(" ");
1464 //______________________________________________________________________________
1465 void AliAnalysisManager::ResetAnalysis()
1467 // Reset all execution flags and clean containers.
1471 //______________________________________________________________________________
1472 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1474 // Start analysis having a grid handler.
1475 if (!fGridHandler) {
1476 Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1477 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1481 return StartAnalysis(type, tree, nentries, firstentry);
1484 //______________________________________________________________________________
1485 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1487 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1488 // MIX. Process nentries starting from firstentry
1490 // Backup current directory and make sure gDirectory points to gROOT
1491 TDirectory *cdir = gDirectory;
1494 Error("StartAnalysis","Analysis manager was not initialized !");
1498 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1499 fMaxEntries = nentries;
1501 TString anaType = type;
1503 fMode = kLocalAnalysis;
1504 Bool_t runlocalinit = kTRUE;
1505 if (anaType.Contains("file")) {
1506 runlocalinit = kFALSE;
1509 if (anaType.Contains("proof")) fMode = kProofAnalysis;
1510 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1511 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
1513 if (fMode == kGridAnalysis) {
1515 if (!anaType.Contains("terminate")) {
1516 if (!fGridHandler) {
1517 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1518 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1522 // Write analysis manager in the analysis file
1523 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1524 // run local task configuration
1525 TIter nextTask(fTasks);
1526 AliAnalysisTask *task;
1527 while ((task=(AliAnalysisTask*)nextTask())) {
1531 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1532 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1537 // Terminate grid analysis
1538 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1539 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1540 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1541 if (!fGridHandler->MergeOutputs()) {
1542 // Return if outputs could not be merged or if it alien handler
1543 // was configured for offline mode or local testing.
1548 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1549 ImportWrappers(NULL);
1555 SetEventLoop(kFALSE);
1556 // Enable event loop mode if a tree was provided
1557 if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1560 TString ttype = "TTree";
1561 if (tree && tree->IsA() == TChain::Class()) {
1562 chain = (TChain*)tree;
1563 if (!chain || !chain->GetListOfFiles()->First()) {
1564 Error("StartAnalysis", "Cannot process null or empty chain...");
1571 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1572 if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1573 // Initialize locally all tasks (happens for all modes)
1575 AliAnalysisTask *task;
1577 while ((task=(AliAnalysisTask*)next())) {
1581 if (getsysInfo) AliSysInfo::AddStamp("LocalInit_all", 0);
1585 case kLocalAnalysis:
1586 if (!tree && !fGridHandler) {
1587 TIter nextT(fTasks);
1588 // Call CreateOutputObjects for all tasks
1590 Bool_t dirStatus = TH1::AddDirectoryStatus();
1591 while ((task=(AliAnalysisTask*)nextT())) {
1592 TH1::AddDirectory(kFALSE);
1593 task->CreateOutputObjects();
1594 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1598 TH1::AddDirectory(dirStatus);
1599 if (IsExternalLoop()) {
1600 Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1601 \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1608 fSelector = new AliAnalysisSelector(this);
1609 // Check if a plugin handler is used
1611 // Get the chain from the plugin
1612 TString dataType = "esdTree";
1613 if (fInputEventHandler) {
1614 dataType = fInputEventHandler->GetDataType();
1618 chain = fGridHandler->GetChainForTestMode(dataType);
1620 Error("StartAnalysis", "No chain for test mode. Aborting.");
1623 cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1624 retv = chain->Process(fSelector, "", nentries, firstentry);
1627 // Run tree-based analysis via AliAnalysisSelector
1628 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1629 retv = tree->Process(fSelector, "", nentries, firstentry);
1631 case kProofAnalysis:
1633 // Check if the plugin is used
1635 return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1637 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1638 Error("StartAnalysis", "No PROOF!!! Exiting.");
1642 line = Form("gProof->AddInput((TObject*)%p);", this);
1643 gROOT->ProcessLine(line);
1646 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1647 retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1649 Error("StartAnalysis", "No chain!!! Exiting.");
1656 if (!anaType.Contains("terminate")) {
1657 if (!fGridHandler) {
1658 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1659 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1663 // Write analysis manager in the analysis file
1664 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1665 // Start the analysis via the handler
1666 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1667 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1672 // Terminate grid analysis
1673 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1674 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1675 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1676 if (!fGridHandler->MergeOutputs()) {
1677 // Return if outputs could not be merged or if it alien handler
1678 // was configured for offline mode or local testing.
1683 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1684 ImportWrappers(NULL);
1688 case kMixingAnalysis:
1689 // Run event mixing analysis
1691 Error("StartAnalysis", "Cannot run event mixing without event pool");
1695 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1696 fSelector = new AliAnalysisSelector(this);
1697 while ((chain=fEventPool->GetNextChain())) {
1699 // Call NotifyBinChange for all tasks
1700 while ((task=(AliAnalysisTask*)next()))
1701 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1702 retv = chain->Process(fSelector);
1704 Error("StartAnalysis", "Mixing analysis failed");
1709 PackOutput(fSelector->GetOutputList());
1716 //______________________________________________________________________________
1717 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1719 // Start analysis for this manager on a given dataset. Analysis task can be:
1720 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1722 Error("StartAnalysis","Analysis manager was not initialized !");
1726 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1727 TString anaType = type;
1729 if (!anaType.Contains("proof")) {
1730 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1733 fMode = kProofAnalysis;
1735 SetEventLoop(kTRUE);
1736 // Set the dataset flag
1737 TObject::SetBit(kUseDataSet);
1740 // Start proof analysis using the grid handler
1741 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1742 Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
1745 // Check if the plugin is in test mode
1746 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
1747 dataset = "test_collection";
1749 dataset = fGridHandler->GetProofDataSet();
1753 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1754 Error("StartAnalysis", "No PROOF!!! Exiting.");
1758 // Initialize locally all tasks
1760 AliAnalysisTask *task;
1761 while ((task=(AliAnalysisTask*)next())) {
1765 line = Form("gProof->AddInput((TObject*)%p);", this);
1766 gROOT->ProcessLine(line);
1768 line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1769 dataset, nentries, firstentry);
1770 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1771 retv = (Long_t)gROOT->ProcessLine(line);
1775 //______________________________________________________________________________
1776 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
1778 // Opens according the option the file specified by cont->GetFileName() and changes
1779 // current directory to cont->GetFolderName(). If the file was already opened, it
1780 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
1781 // be optionally ignored.
1782 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1783 TString filename = cont->GetFileName();
1785 if (filename.IsNull()) {
1786 ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
1789 if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
1791 f = mgr->OpenProofFile(cont,option);
1793 // Check first if the file is already opened
1794 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1796 // Check if option "UPDATE" was preserved
1797 TString opt(option);
1799 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1800 ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1802 f = TFile::Open(filename, option);
1805 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1809 // Check for a folder request
1810 TString dir = cont->GetFolderName();
1811 if (!dir.IsNull()) {
1812 if (!f->GetDirectory(dir)) f->mkdir(dir);
1817 ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
1818 cont->SetFile(NULL);
1822 //______________________________________________________________________________
1823 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
1825 // Opens a special output file used in PROOF.
1827 TString filename = cont->GetFileName();
1828 if (cont == fCommonOutput) {
1829 if (fOutputEventHandler) {
1830 if (strlen(extaod)) filename = extaod;
1831 filename = fOutputEventHandler->GetOutputFileName();
1833 else Fatal("OpenProofFile","No output container. Exiting.");
1836 if (fMode!=kProofAnalysis || !fSelector) {
1837 Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
1840 if (fSpecialOutputLocation.Length()) {
1841 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1843 // Check if option "UPDATE" was preserved
1844 TString opt(option);
1846 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1847 ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1849 f = new TFile(filename, option);
1851 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1855 // Check for a folder request
1856 TString dir = cont->GetFolderName();
1858 if (!f->GetDirectory(dir)) f->mkdir(dir);
1863 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1864 cont->SetFile(NULL);
1867 // Check if there is already a proof output file in the output list
1868 TObject *pof = fSelector->GetOutputList()->FindObject(filename);
1870 // Get the actual file
1871 line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
1872 filename = (const char*)gROOT->ProcessLine(line);
1874 printf("File: %s already booked via TProofOutputFile\n", filename.Data());
1876 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1878 Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
1881 // Check if option "UPDATE" was preserved
1882 TString opt(option);
1884 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1885 Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
1887 if (cont->IsRegisterDataset()) {
1888 TString dsetName = filename;
1889 dsetName.ReplaceAll(".root", cont->GetTitle());
1890 dsetName.ReplaceAll(":","_");
1891 if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
1892 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
1894 if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
1895 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
1897 if (fDebug > 1) printf("=== %s\n", line.Data());
1898 gROOT->ProcessLine(line);
1899 line = Form("pf->OpenFile(\"%s\");", option);
1900 gROOT->ProcessLine(line);
1903 gROOT->ProcessLine("pf->Print()");
1904 printf(" == proof file name: %s", f->GetName());
1906 // Add to proof output list
1907 line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
1908 if (fDebug > 1) printf("=== %s\n", line.Data());
1909 gROOT->ProcessLine(line);
1911 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1915 // Check for a folder request
1916 TString dir = cont->GetFolderName();
1917 if (!dir.IsNull()) {
1918 if (!f->GetDirectory(dir)) f->mkdir(dir);
1923 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1924 cont->SetFile(NULL);
1928 //______________________________________________________________________________
1929 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1931 // Execute analysis.
1932 static Long64_t nentries = 0;
1933 static TTree *lastTree = 0;
1934 static TStopwatch *timer = new TStopwatch();
1935 if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
1937 if (fTree && (fTree != lastTree)) {
1938 nentries += fTree->GetEntries();
1941 if (!fNcalls) timer->Start();
1942 if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, TMath::Min(fMaxEntries,nentries), timer, kFALSE);
1945 TDirectory *cdir = gDirectory;
1946 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1947 if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
1949 Error("ExecAnalysis", "Analysis manager was not initialized !");
1954 AliAnalysisTask *task;
1955 // Check if the top tree is active.
1957 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1958 AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
1960 // De-activate all tasks
1961 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1962 AliAnalysisDataContainer *cont = fCommonInput;
1963 if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
1965 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1969 cont->SetData(fTree); // This will notify all consumers
1970 Long64_t entry = fTree->GetTree()->GetReadEntry();
1972 // Call BeginEvent() for optional input/output and MC services
1973 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1974 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1975 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1977 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1978 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
1980 // Execute the tasks
1981 // TIter next1(cont->GetConsumers());
1982 TIter next1(fTopTasks);
1984 while ((task=(AliAnalysisTask*)next1())) {
1986 cout << " Executing task " << task->GetName() << endl;
1988 task->ExecuteTask(option);
1990 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1991 AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
1995 // Call FinishEvent() for optional output and MC services
1996 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1997 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1998 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1999 // Gather system information if requested
2000 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2001 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
2005 // The event loop is not controlled by TSelector
2007 // Call BeginEvent() for optional input/output and MC services
2008 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
2009 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
2010 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
2012 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2013 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2014 TIter next2(fTopTasks);
2015 while ((task=(AliAnalysisTask*)next2())) {
2016 task->SetActive(kTRUE);
2018 cout << " Executing task " << task->GetName() << endl;
2020 task->ExecuteTask(option);
2024 // Call FinishEvent() for optional output and MC services
2025 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2026 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2027 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2028 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2029 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2033 //______________________________________________________________________________
2034 Bool_t AliAnalysisManager::IsPipe(std::ostream &out)
2036 // Check if the stdout is connected to a pipe (C.Holm)
2037 Bool_t ispipe = kFALSE;
2038 out.seekp(0, std::ios_base::cur);
2041 if (errno == ESPIPE) ispipe = kTRUE;
2046 //______________________________________________________________________________
2047 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2049 // Set the input event handler and create a container for it.
2050 fInputEventHandler = handler;
2051 fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2054 //______________________________________________________________________________
2055 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2057 // Set the input event handler and create a container for it.
2058 fOutputEventHandler = handler;
2059 fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2060 fCommonOutput->SetSpecialOutput();
2063 //______________________________________________________________________________
2064 void AliAnalysisManager::SetDebugLevel(UInt_t level)
2066 // Set verbosity of the analysis manager. If the progress bar is used, the call is ignored
2067 if (TObject::TestBit(kUseProgressBar)) {
2068 Info("SetDebugLevel","Ignored. Disable the progress bar first.");
2074 //______________________________________________________________________________
2075 void AliAnalysisManager::SetUseProgressBar(Bool_t flag, Int_t freq)
2077 // Enable a text mode progress bar. Resets debug level to 0.
2078 Info("SetUseProgressBar", "Progress bar enabled, updated every %d events.\n ### NOTE: Debug level reset to 0 ###", freq);
2079 TObject::SetBit(kUseProgressBar,flag);
2080 fPBUpdateFreq = freq;
2084 //______________________________________________________________________________
2085 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2087 // This method is used externally to register output files which are not
2088 // connected to any output container, so that the manager can properly register,
2089 // retrieve or merge them when running in distributed mode. The file names are
2090 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2091 if (fExtraFiles.Contains(fname)) return;
2092 if (fExtraFiles.Length()) fExtraFiles += " ";
2093 fExtraFiles += fname;
2096 //______________________________________________________________________________
2097 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2099 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2103 TObject *pof = source->FindObject(filename);
2104 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2105 Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2108 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2109 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2110 TString clientUrl(chUrl);
2111 TString fullPath_str(fullPath);
2112 if (clientUrl.Contains("localhost")){
2113 TObjArray* array = fullPath_str.Tokenize ( "//" );
2114 TObjString *strobj = ( TObjString *)array->At(1);
2115 TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2116 TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2117 fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2118 fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2119 if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2123 else if (clientUrl.Contains("__lite__")) {
2124 // Special case for ProofLite environement - get file info and copy.
2125 gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2126 fullPath_str = Form("%s/%s", tmp, fullPath);
2129 Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2130 Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename);
2132 Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2136 //______________________________________________________________________________
2137 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2139 // Fill analysis type in the provided string.
2141 case kLocalAnalysis:
2144 case kProofAnalysis:
2150 case kMixingAnalysis:
2155 //______________________________________________________________________________
2156 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2158 // Validate all output files.
2159 TIter next(fOutputs);
2160 AliAnalysisDataContainer *output;
2161 TDirectory *cdir = gDirectory;
2162 TString openedFiles;
2163 while ((output=(AliAnalysisDataContainer*)next())) {
2164 if (output->IsRegisterDataset()) continue;
2165 TString filename = output->GetFileName();
2166 if (filename == "default") {
2167 if (!fOutputEventHandler) continue;
2168 filename = fOutputEventHandler->GetOutputFileName();
2169 // Main AOD may not be there
2170 if (gSystem->AccessPathName(filename)) continue;
2172 // Check if the file is closed
2173 if (openedFiles.Contains(filename)) continue;;
2174 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2176 Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2177 // Clear file list to release object ownership to user.
2181 file = TFile::Open(filename);
2182 if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2183 Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2188 openedFiles += filename;
2195 //______________________________________________________________________________
2196 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2198 // Implements a nice text mode progress bar.
2199 static Long64_t icount = 0;
2200 static TString oname;
2201 static TString nname;
2202 static Long64_t ocurrent = 0;
2203 static Long64_t osize = 0;
2204 static Int_t oseconds = 0;
2205 static TStopwatch *owatch = 0;
2206 static Bool_t oneoftwo = kFALSE;
2207 static Int_t nrefresh = 0;
2208 static Int_t nchecks = 0;
2209 static char lastChar = 0;
2210 const char symbol[4] = {'-','\\','|','/'};
2212 if (!lastChar) lastChar = (IsPipe(std::cerr))?'\r':'\n';
2218 ocurrent = TMath::Abs(current);
2219 osize = TMath::Abs(size);
2220 if (ocurrent > osize) ocurrent=osize;
2225 if ((current % fPBUpdateFreq) != 0) return;
2227 char progress[11] = " ";
2228 Int_t ichar = icount%4;
2233 if (owatch && !last) {
2235 time = owatch->RealTime();
2236 seconds = int(time) % 60;
2237 minutes = (int(time) / 60) % 60;
2238 hours = (int(time) / 60 / 60);
2240 if (oseconds==seconds) {
2244 oneoftwo = !oneoftwo;
2248 if (refresh && oneoftwo) {
2250 if (nchecks <= 0) nchecks = nrefresh+1;
2251 Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2252 oname = Form(" == %d%% ==", pctdone);
2254 Double_t percent = 100.0*ocurrent/osize;
2255 Int_t nchar = Int_t(percent/10);
2256 if (nchar>10) nchar=10;
2258 for (i=0; i<nchar; i++) progress[i] = '=';
2259 progress[nchar] = symbol[ichar];
2260 for (i=nchar+1; i<10; i++) progress[i] = ' ';
2261 progress[10] = '\0';
2264 if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2265 else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2266 else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2268 Int_t full = Int_t(ocurrent > 0 ?
2269 time * (float(osize)/ocurrent) + .5 :
2271 Int_t remain = full - time;
2272 Int_t rsec = remain % 60;
2273 Int_t rmin = (remain / 60) % 60;
2274 Int_t rhour = (remain / 60 / 60);
2275 fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d ETA %.2d:%.2d:%.2d%c",
2276 percent, hours, minutes, seconds, rhour, rmin, rsec, lastChar);
2278 else fprintf(stderr, "[%6.2f %%]%c", percent, lastChar);
2279 if (refresh && oneoftwo) oname = nname;
2280 if (owatch) owatch->Continue();
2289 fprintf(stderr, "\n");
2293 //______________________________________________________________________________
2294 void AliAnalysisManager::DoLoadBranch(const char *name)
2296 // Get tree and load branch if needed.
2297 static Long64_t crtEntry = -100;
2299 if (fAutoBranchHandling || !fTree)
2302 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2304 br = fTree->GetBranch(name);
2306 Error("DoLoadBranch", "Could not find branch %s",name);
2311 if (br->GetReadEntry()==fCurrentEntry) return;
2312 Int_t ret = br->GetEntry(GetCurrentEntry());
2314 Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2315 if (crtEntry != fCurrentEntry) {
2316 CountEvent(1,0,1,0);
2317 crtEntry = fCurrentEntry;
2320 if (crtEntry != fCurrentEntry) {
2321 CountEvent(1,1,0,0);
2322 crtEntry = fCurrentEntry;
2327 //______________________________________________________________________________
2328 void AliAnalysisManager::AddStatisticsTask(UInt_t offlineMask)
2330 // Add the statistics task to the manager.
2332 Info("AddStatisticsTask", "Already added");
2335 TString line = Form("AliAnalysisTaskStat::AddToManager(%u);", offlineMask);
2336 gROOT->ProcessLine(line);
2339 //______________________________________________________________________________
2340 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2342 // Bookkeep current event;
2343 if (!fStatistics) return;
2344 fStatistics->AddInput(ninput);
2345 fStatistics->AddProcessed(nprocessed);
2346 fStatistics->AddFailed(nfailed);
2347 fStatistics->AddAccepted(naccepted);
2350 //______________________________________________________________________________
2351 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2353 // Add a line in the statistics message. If available, the statistics message is written
2354 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2356 if (!strlen(line)) return;
2357 if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2358 fStatisticsMsg += line;
2361 //______________________________________________________________________________
2362 void AliAnalysisManager::WriteStatisticsMsg(Int_t)
2364 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2365 static Bool_t done = kFALSE;
2368 if (!fStatistics) return;
2370 AddStatisticsMsg(Form("Number of input events: %lld",fStatistics->GetNinput()));
2371 AddStatisticsMsg(Form("Number of processed events: %lld",fStatistics->GetNprocessed()));
2372 AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2373 AddStatisticsMsg(Form("Number of accepted events for mask %s: %lld", AliAnalysisStatistics::GetMaskAsString(fStatistics->GetOfflineMask()), fStatistics->GetNaccepted()));
2374 out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2375 fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2376 fStatistics->GetNaccepted()), ios::out);
2377 out << fStatisticsMsg << endl;
2381 //______________________________________________________________________________
2382 const char* AliAnalysisManager::GetOADBPath()
2384 // returns the path of the OADB
2385 // this static function just depends on environment variables
2387 static TString oadbPath;
2389 if (gSystem->Getenv("OADB_PATH"))
2390 oadbPath = gSystem->Getenv("OADB_PATH");
2391 else if (gSystem->Getenv("ALICE_ROOT"))
2392 oadbPath.Form("%s/OADB", gSystem->Getenv("ALICE_ROOT"));
2394 ::Fatal("AliAnalysisManager::GetOADBPath", "Cannot figure out AODB path. Define ALICE_ROOT or OADB_PATH!");