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>
45 #include "AliAnalysisSelector.h"
46 #include "AliAnalysisGrid.h"
47 #include "AliAnalysisTask.h"
48 #include "AliAnalysisDataContainer.h"
49 #include "AliAnalysisDataSlot.h"
50 #include "AliVEventHandler.h"
51 #include "AliVEventPool.h"
52 #include "AliSysInfo.h"
53 #include "AliAnalysisStatistics.h"
55 ClassImp(AliAnalysisManager)
57 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
58 TString AliAnalysisManager::fgCommonFileName = "";
59 Int_t AliAnalysisManager::fPBUpdateFreq = 1;
61 //______________________________________________________________________________
62 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
65 fInputEventHandler(NULL),
66 fOutputEventHandler(NULL),
67 fMCtruthEventHandler(NULL),
71 fMode(kLocalAnalysis),
75 fSpecialOutputLocation(""),
88 fAutoBranchHandling(kTRUE),
97 // Default constructor.
98 fgAnalysisManager = this;
99 fgCommonFileName = "AnalysisResults.root";
100 fTasks = new TObjArray();
101 fTopTasks = new TObjArray();
102 fZombies = new TObjArray();
103 fContainers = new TObjArray();
104 fInputs = new TObjArray();
105 fOutputs = new TObjArray();
106 fParamCont = new TObjArray();
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;
158 //______________________________________________________________________________
159 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
162 if (&other != this) {
163 TNamed::operator=(other);
164 fInputEventHandler = other.fInputEventHandler;
165 fOutputEventHandler = other.fOutputEventHandler;
166 fMCtruthEventHandler = other.fMCtruthEventHandler;
167 fEventPool = other.fEventPool;
170 fNSysInfo = other.fNSysInfo;
172 fInitOK = other.fInitOK;
173 fIsRemote = other.fIsRemote;
174 fDebug = other.fDebug;
175 fTasks = new TObjArray(*other.fTasks);
176 fTopTasks = new TObjArray(*other.fTopTasks);
177 fZombies = new TObjArray(*other.fZombies);
178 fContainers = new TObjArray(*other.fContainers);
179 fInputs = new TObjArray(*other.fInputs);
180 fOutputs = new TObjArray(*other.fOutputs);
181 fParamCont = new TObjArray(*other.fParamCont);
183 fCommonOutput = NULL;
186 fExtraFiles = other.fExtraFiles;
187 fgCommonFileName = "AnalysisResults.root";
188 fgAnalysisManager = this;
189 fAutoBranchHandling = other.fAutoBranchHandling;
190 fTable.Clear("nodelete");
191 fRunFromPath = other.fRunFromPath;
192 fNcalls = other. fNcalls;
193 fMaxEntries = other.fMaxEntries;
194 fStatisticsMsg = other.fStatisticsMsg;
195 fRequestedBranches = other.fRequestedBranches;
196 fStatistics = other.fStatistics;
201 //______________________________________________________________________________
202 AliAnalysisManager::~AliAnalysisManager()
205 if (fTasks) {fTasks->Delete(); delete fTasks;}
206 if (fTopTasks) delete fTopTasks;
207 if (fZombies) delete fZombies;
208 if (fContainers) {fContainers->Delete(); delete fContainers;}
209 if (fInputs) delete fInputs;
210 if (fOutputs) delete fOutputs;
211 if (fParamCont) delete fParamCont;
212 if (fGridHandler) delete fGridHandler;
213 if (fInputEventHandler) delete fInputEventHandler;
214 if (fOutputEventHandler) delete fOutputEventHandler;
215 if (fMCtruthEventHandler) delete fMCtruthEventHandler;
216 if (fEventPool) delete fEventPool;
217 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
220 //______________________________________________________________________________
221 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
223 // Read one entry of the tree or a whole branch.
224 fCurrentEntry = entry;
225 if (!fAutoBranchHandling)
227 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : -1;
230 //______________________________________________________________________________
231 Int_t AliAnalysisManager::GetRunFromAlienPath(const char *path)
233 // Attempt to extract run number from input data path. Works only for paths to
234 // alice data in alien.
235 // sim: /alice/sim/<production>/run_no/...
236 // data: /alice/data/year/period/000run_no/... (ESD or AOD)
237 TString type = "unknown";
239 if (s.Contains("/alice/data")) type = "real";
240 else if (s.Contains("/alice/sim")) type = "simulated";
243 ind1 = s.Index("/00");
245 ind2 = s.Index("/",ind1+1);
246 if (ind2>0) srun = s(ind1+1, ind2-ind1-1);
249 ind1 = s.Index("/LHC");
251 ind1 = s.Index("/",ind1+1);
253 ind2 = s.Index("/",ind1+1);
254 if (ind2>0) srun = s(ind1+1, ind2-ind1-1);
258 Int_t run = srun.Atoi();
259 if (run>0) printf("=== GetRunFromAlienPath: run %d of %s data ===\n", run, type.Data());
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;
435 if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) return kFALSE;
437 fTable.Clear("nodelete"); // clearing the hash table may not be needed -> C.L.
438 if (fMode == kProofAnalysis) fIsRemote = kTRUE;
440 TFile *curfile = fTree->GetCurrentFile();
442 Error("Notify","No current file");
446 if (fDebug > 1) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
447 Int_t run = AliAnalysisManager::GetRunFromAlienPath(curfile->GetName());
448 if (run && (run != fRunFromPath)) {
450 if (fDebug > 1) printf(" ### run found from path: %d\n", run);
453 AliAnalysisTask *task;
455 // Call Notify of the event handlers
456 if (fInputEventHandler) {
457 fInputEventHandler->Notify(curfile->GetName());
460 if (fOutputEventHandler) {
461 fOutputEventHandler->Notify(curfile->GetName());
464 if (fMCtruthEventHandler) {
465 fMCtruthEventHandler->Notify(curfile->GetName());
468 // Call Notify for all tasks
469 while ((task=(AliAnalysisTask*)next()))
472 if (fDebug > 1) printf("<-AliAnalysisManager::Notify()\n");
476 //______________________________________________________________________________
477 Bool_t AliAnalysisManager::Process(Long64_t)
479 // The Process() function is called for each entry in the tree (or possibly
480 // keyed object in the case of PROOF) to be processed. The entry argument
481 // specifies which entry in the currently loaded tree is to be processed.
482 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
483 // to read either all or the required parts of the data. When processing
484 // keyed objects with PROOF, the object is already loaded and is available
485 // via the fObject pointer.
487 // This function should contain the "body" of the analysis. It can contain
488 // simple or elaborate selection criteria, run algorithms on the data
489 // of the event and typically fill histograms.
491 // WARNING when a selector is used with a TChain, you must use
492 // the pointer to the current TTree to call GetEntry(entry).
493 // The entry is always the local entry number in the current tree.
494 // Assuming that fChain is the pointer to the TChain being processed,
495 // use fChain->GetTree()->GetEntry(entry).
497 // This method is obsolete. ExecAnalysis is called instead.
501 //______________________________________________________________________________
502 void AliAnalysisManager::PackOutput(TList *target)
504 // Pack all output data containers in the output list. Called at SlaveTerminate
505 // stage in PROOF case for each slave.
506 if (fDebug > 1) printf("->AliAnalysisManager::PackOutput()\n");
508 Error("PackOutput", "No target. Exiting.");
511 TDirectory *cdir = gDirectory;
513 if (fInputEventHandler) fInputEventHandler ->Terminate();
514 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
515 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
518 // Call FinishTaskOutput() for each event loop task (not called for
519 // post-event loop tasks - use Terminate() fo those)
520 TIter nexttask(fTasks);
521 AliAnalysisTask *task;
522 while ((task=(AliAnalysisTask*)nexttask())) {
523 if (!task->IsPostEventLoop()) {
524 if (fDebug > 1) printf("->FinishTaskOutput: task %s\n", task->GetName());
525 task->FinishTaskOutput();
527 if (fDebug > 1) printf("<-FinishTaskOutput: task %s\n", task->GetName());
530 // Write statistics message on the workers.
531 if (fStatistics) WriteStatisticsMsg(fNcalls);
533 if (fMode == kProofAnalysis) {
534 TIter next(fOutputs);
535 AliAnalysisDataContainer *output;
536 Bool_t isManagedByHandler = kFALSE;
539 while ((output=(AliAnalysisDataContainer*)next())) {
540 // Do not consider outputs of post event loop tasks
541 isManagedByHandler = kFALSE;
542 if (output->GetProducer() && output->GetProducer()->IsPostEventLoop()) continue;
543 const char *filename = output->GetFileName();
544 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
545 isManagedByHandler = kTRUE;
546 printf("#### Handler output. Extra: %s\n", fExtraFiles.Data());
547 filename = fOutputEventHandler->GetOutputFileName();
549 // Check if data was posted to this container. If not, issue an error.
550 if (!output->GetData() && !isManagedByHandler) {
551 Error("PackOutput", "No data for output container %s. Forgot to PostData ?", output->GetName());
554 if (!output->IsSpecialOutput()) {
556 if (strlen(filename) && !isManagedByHandler) {
557 // Backup current folder
558 TDirectory *opwd = gDirectory;
559 // File resident outputs.
560 // Check first if the file exists.
561 TString openoption = "RECREATE";
562 Bool_t firsttime = kTRUE;
563 if (filestmp.FindObject(output->GetFileName())) {
566 filestmp.Add(new TNamed(output->GetFileName(),""));
568 if (!gSystem->AccessPathName(output->GetFileName()) && !firsttime) openoption = "UPDATE";
569 // TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
570 // Save data to file, then close.
571 if (output->GetData()->InheritsFrom(TCollection::Class())) {
572 // If data is a collection, we set the name of the collection
573 // as the one of the container and we save as a single key.
574 TCollection *coll = (TCollection*)output->GetData();
575 coll->SetName(output->GetName());
576 // coll->Write(output->GetName(), TObject::kSingleKey);
578 if (output->GetData()->InheritsFrom(TTree::Class())) {
579 TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
580 // Save data to file, then close.
581 TTree *tree = (TTree*)output->GetData();
582 // Check if tree is in memory
583 if (tree->GetDirectory()==gROOT) tree->SetDirectory(gDirectory);
587 // output->GetData()->Write();
590 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
592 // printf(" file %s listing content:\n", filename);
595 // Clear file list to release object ownership to user.
598 output->SetFile(NULL);
599 // Restore current directory
600 if (opwd) opwd->cd();
602 // Memory-resident outputs
603 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
605 AliAnalysisDataWrapper *wrap = 0;
606 if (isManagedByHandler) {
607 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
608 wrap->SetName(output->GetName());
610 else wrap =output->ExportData();
611 // Output wrappers must NOT delete data after merging - the user owns them
612 wrap->SetDeleteData(kFALSE);
615 // Special outputs. The file must be opened and connected to the container.
616 TDirectory *opwd = gDirectory;
617 TFile *file = output->GetFile();
619 AliAnalysisTask *producer = output->GetProducer();
621 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
622 output->GetFileName(), output->GetName(), producer->ClassName());
625 TString outFilename = file->GetName();
626 if (fDebug > 1) printf("PackOutput %s: special output\n", output->GetName());
627 if (isManagedByHandler) {
628 // Terminate IO for files managed by the output handler
629 // file->Write() moved to AOD handler (A.G. 11.01.10)
630 // if (file) file->Write();
631 if (file && fDebug > 2) {
632 printf(" handled file %s listing content:\n", file->GetName());
635 fOutputEventHandler->TerminateIO();
638 // Release object ownership to users after writing data to file
639 if (output->GetData()->InheritsFrom(TCollection::Class())) {
640 // If data is a collection, we set the name of the collection
641 // as the one of the container and we save as a single key.
642 TCollection *coll = (TCollection*)output->GetData();
643 coll->SetName(output->GetName());
644 coll->Write(output->GetName(), TObject::kSingleKey);
646 if (output->GetData()->InheritsFrom(TTree::Class())) {
647 TTree *tree = (TTree*)output->GetData();
648 tree->SetDirectory(file);
651 output->GetData()->Write();
655 printf(" file %s listing content:\n", output->GetFileName());
658 // Clear file list to release object ownership to user.
661 output->SetFile(NULL);
663 // Restore current directory
664 if (opwd) opwd->cd();
665 // Check if a special output location was provided or the output files have to be merged
666 if (strlen(fSpecialOutputLocation.Data())) {
667 TString remote = fSpecialOutputLocation;
669 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
670 if (remote.BeginsWith("alien:")) {
671 gROOT->ProcessLine("TGrid::Connect(\"alien:\", gProofServ->GetUser());");
672 remote += outFilename;
673 remote.ReplaceAll(".root", Form("_%d.root", gid));
675 remote += Form("%s_%d_", gSystem->HostName(), gid);
676 remote += outFilename;
679 Info("PackOutput", "Output file for container %s to be copied \n at: %s. No merging.",
680 output->GetName(), remote.Data());
681 TFile::Cp ( outFilename.Data(), remote.Data() );
682 // Copy extra outputs
683 if (fExtraFiles.Length() && isManagedByHandler) {
684 TObjArray *arr = fExtraFiles.Tokenize(" ");
686 TIter nextfilename(arr);
687 while ((os=(TObjString*)nextfilename())) {
688 outFilename = os->GetString();
689 remote = fSpecialOutputLocation;
691 if (remote.BeginsWith("alien://")) {
692 remote += outFilename;
693 remote.ReplaceAll(".root", Form("_%d.root", gid));
695 remote += Form("%s_%d_", gSystem->HostName(), gid);
696 remote += outFilename;
699 Info("PackOutput", "Extra AOD file %s to be copied \n at: %s. No merging.",
700 outFilename.Data(), remote.Data());
701 TFile::Cp ( outFilename.Data(), remote.Data() );
706 // No special location specified-> use TProofOutputFile as merging utility
707 // The file at this output slot must be opened in CreateOutputObjects
708 if (fDebug > 1) printf(" File for container %s to be merged via file merger...\n", output->GetName());
714 if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
717 //______________________________________________________________________________
718 void AliAnalysisManager::ImportWrappers(TList *source)
720 // Import data in output containers from wrappers coming in source.
721 if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
722 TIter next(fOutputs);
723 AliAnalysisDataContainer *cont;
724 AliAnalysisDataWrapper *wrap;
726 Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
727 TDirectory *cdir = gDirectory;
728 while ((cont=(AliAnalysisDataContainer*)next())) {
730 if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
731 if (cont->IsRegisterDataset()) continue;
732 const char *filename = cont->GetFileName();
733 Bool_t isManagedByHandler = kFALSE;
734 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
735 isManagedByHandler = kTRUE;
736 filename = fOutputEventHandler->GetOutputFileName();
738 if (cont->IsSpecialOutput() || inGrid) {
739 if (strlen(fSpecialOutputLocation.Data())) continue;
740 // Copy merged file from PROOF scratch space.
741 // In case of grid the files are already in the current directory.
743 if (isManagedByHandler && fExtraFiles.Length()) {
744 // Copy extra registered dAOD files.
745 TObjArray *arr = fExtraFiles.Tokenize(" ");
747 TIter nextfilename(arr);
748 while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
751 if (!GetFileFromWrapper(filename, source)) continue;
753 // Normally we should connect data from the copied file to the
754 // corresponding output container, but it is not obvious how to do this
755 // automatically if several objects in file...
756 TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
757 if (!f) f = TFile::Open(filename, "READ");
759 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[MB]/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[MB] for COO and Terminate");
1131 hist->GetYaxis()->SetTitle("deltaVM [MB]");
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[MB] per event");
1150 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1155 if (cdir) cdir->cd();
1158 //______________________________________________________________________________
1159 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1161 // Adds a user task to the global list of tasks.
1163 Error("AddTask", "Cannot add task %s since InitAnalysis was already called", task->GetName());
1167 if (fTasks->FindObject(task)) {
1168 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1171 task->SetActive(kFALSE);
1175 //______________________________________________________________________________
1176 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1178 // Retreive task by name.
1179 if (!fTasks) return NULL;
1180 return (AliAnalysisTask*)fTasks->FindObject(name);
1183 //______________________________________________________________________________
1184 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
1185 TClass *datatype, EAliAnalysisContType type, const char *filename)
1187 // Create a data container of a certain type. Types can be:
1188 // kExchangeContainer = 0, used to exchange data between tasks
1189 // kInputContainer = 1, used to store input data
1190 // kOutputContainer = 2, used for writing result to a file
1191 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1192 // the output object to a folder inside the output file
1193 if (fContainers->FindObject(name)) {
1194 Error("CreateContainer","A container named %s already defined !",name);
1197 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1198 fContainers->Add(cont);
1200 case kInputContainer:
1203 case kOutputContainer:
1204 fOutputs->Add(cont);
1205 if (filename && strlen(filename)) {
1206 cont->SetFileName(filename);
1207 cont->SetDataOwned(kFALSE); // data owned by the file
1210 case kParamContainer:
1211 fParamCont->Add(cont);
1212 if (filename && strlen(filename)) {
1213 cont->SetFileName(filename);
1214 cont->SetDataOwned(kFALSE); // data owned by the file
1217 case kExchangeContainer:
1223 //______________________________________________________________________________
1224 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1225 AliAnalysisDataContainer *cont)
1227 // Connect input of an existing task to a data container.
1229 Error("ConnectInput", "Task pointer is NULL");
1232 if (!fTasks->FindObject(task)) {
1234 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1236 Bool_t connected = task->ConnectInput(islot, cont);
1240 //______________________________________________________________________________
1241 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1242 AliAnalysisDataContainer *cont)
1244 // Connect output of an existing task to a data container.
1246 Error("ConnectOutput", "Task pointer is NULL");
1249 if (!fTasks->FindObject(task)) {
1251 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1253 Bool_t connected = task->ConnectOutput(islot, cont);
1257 //______________________________________________________________________________
1258 void AliAnalysisManager::CleanContainers()
1260 // Clean data from all containers that have already finished all client tasks.
1261 TIter next(fContainers);
1262 AliAnalysisDataContainer *cont;
1263 while ((cont=(AliAnalysisDataContainer *)next())) {
1264 if (cont->IsOwnedData() &&
1265 cont->IsDataReady() &&
1266 cont->ClientsExecuted()) cont->DeleteData();
1270 //______________________________________________________________________________
1271 Bool_t AliAnalysisManager::InitAnalysis()
1273 // Initialization of analysis chain of tasks. Should be called after all tasks
1274 // and data containers are properly connected
1275 // Reset flag and remove valid_outputs file if exists
1276 if (fInitOK) return kTRUE;
1277 if (!gSystem->AccessPathName("outputs_valid"))
1278 gSystem->Unlink("outputs_valid");
1279 // Check for top tasks (depending only on input data containers)
1280 if (!fTasks->First()) {
1281 Error("InitAnalysis", "Analysis has no tasks !");
1285 AliAnalysisTask *task;
1286 AliAnalysisDataContainer *cont;
1289 Bool_t iszombie = kFALSE;
1290 Bool_t istop = kTRUE;
1292 while ((task=(AliAnalysisTask*)next())) {
1295 Int_t ninputs = task->GetNinputs();
1296 for (i=0; i<ninputs; i++) {
1297 cont = task->GetInputSlot(i)->GetContainer();
1301 fZombies->Add(task);
1305 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
1306 i, task->GetName());
1308 if (iszombie) continue;
1309 // Check if cont is an input container
1310 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1311 // Connect to parent task
1315 fTopTasks->Add(task);
1319 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1322 // Check now if there are orphan tasks
1323 for (i=0; i<ntop; i++) {
1324 task = (AliAnalysisTask*)fTopTasks->At(i);
1329 while ((task=(AliAnalysisTask*)next())) {
1330 if (!task->IsUsed()) {
1332 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1335 // Check the task hierarchy (no parent task should depend on data provided
1336 // by a daughter task)
1337 for (i=0; i<ntop; i++) {
1338 task = (AliAnalysisTask*)fTopTasks->At(i);
1339 if (task->CheckCircularDeps()) {
1340 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1345 // Check that all containers feeding post-event loop tasks are in the outputs list
1346 TIter nextcont(fContainers); // loop over all containers
1347 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1348 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1349 if (cont->HasConsumers()) {
1350 // Check if one of the consumers is post event loop
1351 TIter nextconsumer(cont->GetConsumers());
1352 while ((task=(AliAnalysisTask*)nextconsumer())) {
1353 if (task->IsPostEventLoop()) {
1354 fOutputs->Add(cont);
1361 // Check if all special output containers have a file name provided
1362 TIter nextout(fOutputs);
1363 while ((cont=(AliAnalysisDataContainer*)nextout())) {
1364 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1365 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1369 // Initialize requested branch list if needed
1370 if (!fAutoBranchHandling) {
1372 while ((task=(AliAnalysisTask*)next())) {
1373 if (!task->HasBranches()) {
1374 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\"",
1375 task->GetName(), task->ClassName());
1378 if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1379 Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1382 TString taskbranches;
1383 task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1384 if (taskbranches.IsNull()) {
1385 Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1386 task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1389 AddBranches(taskbranches);
1396 //______________________________________________________________________________
1397 void AliAnalysisManager::AddBranches(const char *branches)
1399 // Add branches to the existing fRequestedBranches.
1400 TString br(branches);
1401 TObjArray *arr = br.Tokenize(",");
1404 while ((obj=next())) {
1405 if (!fRequestedBranches.Contains(obj->GetName())) {
1406 if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1407 fRequestedBranches += obj->GetName();
1413 //______________________________________________________________________________
1414 void AliAnalysisManager::CheckBranches(Bool_t load)
1416 // The method checks the input branches to be loaded during the analysis.
1417 if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;
1418 TObjArray *arr = fRequestedBranches.Tokenize(",");
1421 while ((obj=next())) {
1422 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1424 br = fTree->GetBranch(obj->GetName());
1426 Error("CheckBranches", "Could not find branch %s",obj->GetName());
1431 if (load && br->GetReadEntry()!=GetCurrentEntry()) br->GetEntry(GetCurrentEntry());
1436 //______________________________________________________________________________
1437 Bool_t AliAnalysisManager::CheckTasks() const
1439 // Check consistency of tasks.
1440 Int_t ntasks = fTasks->GetEntries();
1442 Error("CheckTasks", "No tasks connected to the manager. This may be due to forgetting to compile the task or to load their library.");
1445 // Get the pointer to AliAnalysisTaskSE::Class()
1446 TClass *badptr = (TClass*)gROOT->ProcessLine("AliAnalysisTaskSE::Class()");
1447 // Loop all tasks to check if their corresponding library was loaded
1450 while ((obj=next())) {
1451 if (obj->IsA() == badptr) {
1452 Error("CheckTasks", "##################\n \
1453 Class for task %s NOT loaded. You probably forgot to load the library for this task (or compile it dynamically).\n###########################\n",obj->GetName());
1460 //______________________________________________________________________________
1461 void AliAnalysisManager::PrintStatus(Option_t *option) const
1463 // Print task hierarchy.
1465 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1468 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1470 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1471 TIter next(fTopTasks);
1472 AliAnalysisTask *task;
1473 while ((task=(AliAnalysisTask*)next()))
1474 task->PrintTask(option);
1476 if (!fAutoBranchHandling && !fRequestedBranches.IsNull())
1477 printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1479 TString sopt(option);
1482 if (sopt.Contains("ALL"))
1484 if ( fOutputEventHandler )
1486 cout << TString('_',78) << endl;
1487 cout << "OutputEventHandler:" << endl;
1488 fOutputEventHandler->Print(" ");
1493 //______________________________________________________________________________
1494 void AliAnalysisManager::ResetAnalysis()
1496 // Reset all execution flags and clean containers.
1500 //______________________________________________________________________________
1501 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1503 // Start analysis having a grid handler.
1504 if (!fGridHandler) {
1505 Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1506 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1510 return StartAnalysis(type, tree, nentries, firstentry);
1513 //______________________________________________________________________________
1514 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1516 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1517 // MIX. Process nentries starting from firstentry
1519 // Backup current directory and make sure gDirectory points to gROOT
1520 TDirectory *cdir = gDirectory;
1523 Error("StartAnalysis","Analysis manager was not initialized !");
1527 if (!CheckTasks()) Fatal("StartAnalysis", "Not all needed libraries were loaded");
1529 printf("StartAnalysis %s\n",GetName());
1530 AliLog::SetGlobalLogLevel(AliLog::kInfo);
1532 fMaxEntries = nentries;
1534 TString anaType = type;
1536 fMode = kLocalAnalysis;
1537 Bool_t runlocalinit = kTRUE;
1538 if (anaType.Contains("file")) {
1539 runlocalinit = kFALSE;
1542 if (anaType.Contains("proof")) fMode = kProofAnalysis;
1543 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1544 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
1546 if (fMode == kGridAnalysis) {
1548 if (!anaType.Contains("terminate")) {
1549 if (!fGridHandler) {
1550 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1551 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1555 // Write analysis manager in the analysis file
1556 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1557 // run local task configuration
1558 TIter nextTask(fTasks);
1559 AliAnalysisTask *task;
1560 while ((task=(AliAnalysisTask*)nextTask())) {
1564 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1565 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1570 // Terminate grid analysis
1571 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1572 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1573 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1574 if (!fGridHandler->MergeOutputs()) {
1575 // Return if outputs could not be merged or if it alien handler
1576 // was configured for offline mode or local testing.
1581 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1582 ImportWrappers(NULL);
1588 SetEventLoop(kFALSE);
1589 // Enable event loop mode if a tree was provided
1590 if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1593 TString ttype = "TTree";
1594 if (tree && tree->IsA() == TChain::Class()) {
1595 chain = (TChain*)tree;
1596 if (!chain || !chain->GetListOfFiles()->First()) {
1597 Error("StartAnalysis", "Cannot process null or empty chain...");
1604 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1605 if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1606 // Initialize locally all tasks (happens for all modes)
1608 AliAnalysisTask *task;
1610 while ((task=(AliAnalysisTask*)next())) {
1614 if (getsysInfo) AliSysInfo::AddStamp("LocalInit_all", 0);
1618 case kLocalAnalysis:
1619 if (!tree && !fGridHandler) {
1620 TIter nextT(fTasks);
1621 // Call CreateOutputObjects for all tasks
1623 Bool_t dirStatus = TH1::AddDirectoryStatus();
1624 while ((task=(AliAnalysisTask*)nextT())) {
1625 TH1::AddDirectory(kFALSE);
1626 task->CreateOutputObjects();
1627 if (!task->CheckPostData()) {
1628 Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
1629 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
1630 ########### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ###########", task->GetName(), task->ClassName());
1632 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1636 TH1::AddDirectory(dirStatus);
1637 if (IsExternalLoop()) {
1638 Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1639 \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1646 fSelector = new AliAnalysisSelector(this);
1647 // Check if a plugin handler is used
1649 // Get the chain from the plugin
1650 TString dataType = "esdTree";
1651 if (fInputEventHandler) {
1652 dataType = fInputEventHandler->GetDataType();
1656 chain = fGridHandler->GetChainForTestMode(dataType);
1658 Error("StartAnalysis", "No chain for test mode. Aborting.");
1661 cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1662 retv = chain->Process(fSelector, "", nentries, firstentry);
1665 // Run tree-based analysis via AliAnalysisSelector
1666 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1667 retv = tree->Process(fSelector, "", nentries, firstentry);
1669 case kProofAnalysis:
1671 // Check if the plugin is used
1673 return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1675 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1676 Error("StartAnalysis", "No PROOF!!! Exiting.");
1680 line = Form("gProof->AddInput((TObject*)%p);", this);
1681 gROOT->ProcessLine(line);
1684 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1685 retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1687 Error("StartAnalysis", "No chain!!! Exiting.");
1694 if (!anaType.Contains("terminate")) {
1695 if (!fGridHandler) {
1696 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1697 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1701 // Write analysis manager in the analysis file
1702 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1703 // Start the analysis via the handler
1704 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1705 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1710 // Terminate grid analysis
1711 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1712 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1713 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1714 if (!fGridHandler->MergeOutputs()) {
1715 // Return if outputs could not be merged or if it alien handler
1716 // was configured for offline mode or local testing.
1721 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1722 ImportWrappers(NULL);
1726 case kMixingAnalysis:
1727 // Run event mixing analysis
1729 Error("StartAnalysis", "Cannot run event mixing without event pool");
1733 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1734 fSelector = new AliAnalysisSelector(this);
1735 while ((chain=fEventPool->GetNextChain())) {
1737 // Call NotifyBinChange for all tasks
1738 while ((task=(AliAnalysisTask*)next()))
1739 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1740 retv = chain->Process(fSelector);
1742 Error("StartAnalysis", "Mixing analysis failed");
1747 PackOutput(fSelector->GetOutputList());
1754 //______________________________________________________________________________
1755 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1757 // Start analysis for this manager on a given dataset. Analysis task can be:
1758 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1760 Error("StartAnalysis","Analysis manager was not initialized !");
1764 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1765 TString anaType = type;
1767 if (!anaType.Contains("proof")) {
1768 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1771 fMode = kProofAnalysis;
1773 SetEventLoop(kTRUE);
1774 // Set the dataset flag
1775 TObject::SetBit(kUseDataSet);
1778 // Start proof analysis using the grid handler
1779 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1780 Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
1783 // Check if the plugin is in test mode
1784 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
1785 dataset = "test_collection";
1787 dataset = fGridHandler->GetProofDataSet();
1791 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1792 Error("StartAnalysis", "No PROOF!!! Exiting.");
1796 // Initialize locally all tasks
1798 AliAnalysisTask *task;
1799 while ((task=(AliAnalysisTask*)next())) {
1803 line = Form("gProof->AddInput((TObject*)%p);", this);
1804 gROOT->ProcessLine(line);
1806 line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1807 dataset, nentries, firstentry);
1808 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1809 retv = (Long_t)gROOT->ProcessLine(line);
1813 //______________________________________________________________________________
1814 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
1816 // Opens according the option the file specified by cont->GetFileName() and changes
1817 // current directory to cont->GetFolderName(). If the file was already opened, it
1818 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
1819 // be optionally ignored.
1820 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1821 TString filename = cont->GetFileName();
1823 if (filename.IsNull()) {
1824 ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
1827 if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
1829 f = mgr->OpenProofFile(cont,option);
1831 // Check first if the file is already opened
1832 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1834 // Check if option "UPDATE" was preserved
1835 TString opt(option);
1837 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1838 ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1840 f = TFile::Open(filename, option);
1843 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1847 // Check for a folder request
1848 TString dir = cont->GetFolderName();
1849 if (!dir.IsNull()) {
1850 if (!f->GetDirectory(dir)) f->mkdir(dir);
1855 ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
1856 cont->SetFile(NULL);
1860 //______________________________________________________________________________
1861 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
1863 // Opens a special output file used in PROOF.
1865 TString filename = cont->GetFileName();
1866 if (cont == fCommonOutput) {
1867 if (fOutputEventHandler) {
1868 if (strlen(extaod)) filename = extaod;
1869 filename = fOutputEventHandler->GetOutputFileName();
1871 else Fatal("OpenProofFile","No output container. Exiting.");
1874 if (fMode!=kProofAnalysis || !fSelector) {
1875 Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
1878 if (fSpecialOutputLocation.Length()) {
1879 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1881 // Check if option "UPDATE" was preserved
1882 TString opt(option);
1884 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1885 ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1887 f = new TFile(filename, option);
1889 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1893 // Check for a folder request
1894 TString dir = cont->GetFolderName();
1896 if (!f->GetDirectory(dir)) f->mkdir(dir);
1901 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1902 cont->SetFile(NULL);
1905 // Check if there is already a proof output file in the output list
1906 TObject *pof = fSelector->GetOutputList()->FindObject(filename);
1908 // Get the actual file
1909 line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
1910 filename = (const char*)gROOT->ProcessLine(line);
1912 printf("File: %s already booked via TProofOutputFile\n", filename.Data());
1914 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1916 Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
1919 // Check if option "UPDATE" was preserved
1920 TString opt(option);
1922 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1923 Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
1925 if (cont->IsRegisterDataset()) {
1926 TString dsetName = filename;
1927 dsetName.ReplaceAll(".root", cont->GetTitle());
1928 dsetName.ReplaceAll(":","_");
1929 if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
1930 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
1932 if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
1933 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
1935 if (fDebug > 1) printf("=== %s\n", line.Data());
1936 gROOT->ProcessLine(line);
1937 line = Form("pf->OpenFile(\"%s\");", option);
1938 gROOT->ProcessLine(line);
1941 gROOT->ProcessLine("pf->Print()");
1942 printf(" == proof file name: %s", f->GetName());
1944 // Add to proof output list
1945 line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
1946 if (fDebug > 1) printf("=== %s\n", line.Data());
1947 gROOT->ProcessLine(line);
1949 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1953 // Check for a folder request
1954 TString dir = cont->GetFolderName();
1955 if (!dir.IsNull()) {
1956 if (!f->GetDirectory(dir)) f->mkdir(dir);
1961 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1962 cont->SetFile(NULL);
1966 //______________________________________________________________________________
1967 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1969 // Execute analysis.
1970 static Long64_t nentries = 0;
1971 static TTree *lastTree = 0;
1972 static TStopwatch *timer = new TStopwatch();
1973 // Only the first call to Process will trigger a true Notify. Other Notify
1974 // coming before is ignored.
1975 if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) {
1976 TObject::SetBit(AliAnalysisManager::kTrueNotify);
1979 if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
1981 if (fTree && (fTree != lastTree)) {
1982 nentries += fTree->GetEntries();
1985 if (!fNcalls) timer->Start();
1986 if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, TMath::Min(fMaxEntries,nentries), timer, kFALSE);
1989 TDirectory *cdir = gDirectory;
1990 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1991 if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
1993 Error("ExecAnalysis", "Analysis manager was not initialized !");
1998 AliAnalysisTask *task;
1999 // Check if the top tree is active.
2001 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2002 AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
2004 // De-activate all tasks
2005 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
2006 AliAnalysisDataContainer *cont = fCommonInput;
2007 if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
2009 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
2013 cont->SetData(fTree); // This will notify all consumers
2014 Long64_t entry = fTree->GetTree()->GetReadEntry();
2016 // Call BeginEvent() for optional input/output and MC services
2017 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
2018 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
2019 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
2021 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2022 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2024 // Execute the tasks
2025 // TIter next1(cont->GetConsumers());
2026 TIter next1(fTopTasks);
2028 while ((task=(AliAnalysisTask*)next1())) {
2030 cout << " Executing task " << task->GetName() << endl;
2032 task->ExecuteTask(option);
2034 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2035 AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
2039 // Call FinishEvent() for optional output and MC services
2040 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2041 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2042 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2043 // Gather system information if requested
2044 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2045 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
2049 // The event loop is not controlled by TSelector
2051 // Call BeginEvent() for optional input/output and MC services
2052 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
2053 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
2054 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
2056 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2057 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2058 TIter next2(fTopTasks);
2059 while ((task=(AliAnalysisTask*)next2())) {
2060 task->SetActive(kTRUE);
2062 cout << " Executing task " << task->GetName() << endl;
2064 task->ExecuteTask(option);
2068 // Call FinishEvent() for optional output and MC services
2069 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2070 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2071 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2072 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2073 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2077 //______________________________________________________________________________
2078 Bool_t AliAnalysisManager::IsPipe(std::ostream &out)
2080 // Check if the stdout is connected to a pipe (C.Holm)
2081 Bool_t ispipe = kFALSE;
2082 out.seekp(0, std::ios_base::cur);
2085 if (errno == ESPIPE) ispipe = kTRUE;
2090 //______________________________________________________________________________
2091 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2093 // Set the input event handler and create a container for it.
2094 fInputEventHandler = handler;
2095 if (!fCommonInput) fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2098 //______________________________________________________________________________
2099 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2101 // Set the input event handler and create a container for it.
2102 fOutputEventHandler = handler;
2103 if (!fCommonOutput) fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2104 fCommonOutput->SetSpecialOutput();
2107 //______________________________________________________________________________
2108 void AliAnalysisManager::SetDebugLevel(UInt_t level)
2110 // Set verbosity of the analysis manager. If the progress bar is used, the call is ignored
2111 if (TObject::TestBit(kUseProgressBar)) {
2112 Info("SetDebugLevel","Ignored. Disable the progress bar first.");
2118 //______________________________________________________________________________
2119 void AliAnalysisManager::SetUseProgressBar(Bool_t flag, Int_t freq)
2121 // Enable a text mode progress bar. Resets debug level to 0.
2122 Info("SetUseProgressBar", "Progress bar enabled, updated every %d events.\n ### NOTE: Debug level reset to 0 ###", freq);
2123 TObject::SetBit(kUseProgressBar,flag);
2124 fPBUpdateFreq = freq;
2128 //______________________________________________________________________________
2129 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2131 // This method is used externally to register output files which are not
2132 // connected to any output container, so that the manager can properly register,
2133 // retrieve or merge them when running in distributed mode. The file names are
2134 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2135 if (fExtraFiles.Contains(fname)) return;
2136 if (fExtraFiles.Length()) fExtraFiles += " ";
2137 fExtraFiles += fname;
2140 //______________________________________________________________________________
2141 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2143 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2147 TObject *pof = source->FindObject(filename);
2148 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2149 Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2152 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2153 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2154 TString clientUrl(chUrl);
2155 TString fullPath_str(fullPath);
2156 if (clientUrl.Contains("localhost")){
2157 TObjArray* array = fullPath_str.Tokenize ( "//" );
2158 TObjString *strobj = ( TObjString *)array->At(1);
2159 TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2160 TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2161 fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2162 fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2163 if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2167 else if (clientUrl.Contains("__lite__")) {
2168 // Special case for ProofLite environement - get file info and copy.
2169 gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2170 fullPath_str = Form("%s/%s", tmp, fullPath);
2173 Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2174 Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename);
2176 Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2180 //______________________________________________________________________________
2181 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2183 // Fill analysis type in the provided string.
2185 case kLocalAnalysis:
2188 case kProofAnalysis:
2194 case kMixingAnalysis:
2199 //______________________________________________________________________________
2200 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2202 // Validate all output files.
2203 TIter next(fOutputs);
2204 AliAnalysisDataContainer *output;
2205 TDirectory *cdir = gDirectory;
2206 TString openedFiles;
2207 while ((output=(AliAnalysisDataContainer*)next())) {
2208 if (output->IsRegisterDataset()) continue;
2209 TString filename = output->GetFileName();
2210 if (filename == "default") {
2211 if (!fOutputEventHandler) continue;
2212 filename = fOutputEventHandler->GetOutputFileName();
2213 // Main AOD may not be there
2214 if (gSystem->AccessPathName(filename)) continue;
2216 // Check if the file is closed
2217 if (openedFiles.Contains(filename)) continue;;
2218 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2220 Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2221 // Clear file list to release object ownership to user.
2225 file = TFile::Open(filename);
2226 if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2227 Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2232 openedFiles += filename;
2239 //______________________________________________________________________________
2240 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2242 // Implements a nice text mode progress bar.
2243 static Long64_t icount = 0;
2244 static TString oname;
2245 static TString nname;
2246 static Long64_t ocurrent = 0;
2247 static Long64_t osize = 0;
2248 static Int_t oseconds = 0;
2249 static TStopwatch *owatch = 0;
2250 static Bool_t oneoftwo = kFALSE;
2251 static Int_t nrefresh = 0;
2252 static Int_t nchecks = 0;
2253 static char lastChar = 0;
2254 const char symbol[4] = {'-','\\','|','/'};
2256 if (!lastChar) lastChar = (IsPipe(std::cerr))?'\r':'\n';
2262 ocurrent = TMath::Abs(current);
2263 osize = TMath::Abs(size);
2264 if (ocurrent > osize) ocurrent=osize;
2269 if ((current % fPBUpdateFreq) != 0) return;
2271 char progress[11] = " ";
2272 Int_t ichar = icount%4;
2277 if (owatch && !last) {
2279 time = owatch->RealTime();
2280 seconds = int(time) % 60;
2281 minutes = (int(time) / 60) % 60;
2282 hours = (int(time) / 60 / 60);
2284 if (oseconds==seconds) {
2288 oneoftwo = !oneoftwo;
2292 if (refresh && oneoftwo) {
2294 if (nchecks <= 0) nchecks = nrefresh+1;
2295 Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2296 oname = Form(" == %d%% ==", pctdone);
2298 Double_t percent = 100.0*ocurrent/osize;
2299 Int_t nchar = Int_t(percent/10);
2300 if (nchar>10) nchar=10;
2302 for (i=0; i<nchar; i++) progress[i] = '=';
2303 progress[nchar] = symbol[ichar];
2304 for (i=nchar+1; i<10; i++) progress[i] = ' ';
2305 progress[10] = '\0';
2308 if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2309 else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2310 else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2312 Int_t full = Int_t(ocurrent > 0 ?
2313 time * (float(osize)/ocurrent) + .5 :
2315 Int_t remain = Int_t(full - time);
2316 Int_t rsec = remain % 60;
2317 Int_t rmin = (remain / 60) % 60;
2318 Int_t rhour = (remain / 60 / 60);
2319 fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d ETA %.2d:%.2d:%.2d%c",
2320 percent, hours, minutes, seconds, rhour, rmin, rsec, lastChar);
2322 else fprintf(stderr, "[%6.2f %%]%c", percent, lastChar);
2323 if (refresh && oneoftwo) oname = nname;
2324 if (owatch) owatch->Continue();
2333 fprintf(stderr, "\n");
2337 //______________________________________________________________________________
2338 void AliAnalysisManager::DoLoadBranch(const char *name)
2340 // Get tree and load branch if needed.
2341 static Long64_t crtEntry = -100;
2343 if (fAutoBranchHandling || !fTree)
2346 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2348 br = fTree->GetBranch(name);
2350 Error("DoLoadBranch", "Could not find branch %s",name);
2355 if (br->GetReadEntry()==fCurrentEntry) return;
2356 Int_t ret = br->GetEntry(GetCurrentEntry());
2358 Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2359 if (crtEntry != fCurrentEntry) {
2360 CountEvent(1,0,1,0);
2361 crtEntry = fCurrentEntry;
2364 if (crtEntry != fCurrentEntry) {
2365 CountEvent(1,1,0,0);
2366 crtEntry = fCurrentEntry;
2371 //______________________________________________________________________________
2372 void AliAnalysisManager::AddStatisticsTask(UInt_t offlineMask)
2374 // Add the statistics task to the manager.
2376 Info("AddStatisticsTask", "Already added");
2379 TString line = Form("AliAnalysisTaskStat::AddToManager(%u);", offlineMask);
2380 gROOT->ProcessLine(line);
2383 //______________________________________________________________________________
2384 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2386 // Bookkeep current event;
2387 if (!fStatistics) return;
2388 fStatistics->AddInput(ninput);
2389 fStatistics->AddProcessed(nprocessed);
2390 fStatistics->AddFailed(nfailed);
2391 fStatistics->AddAccepted(naccepted);
2394 //______________________________________________________________________________
2395 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2397 // Add a line in the statistics message. If available, the statistics message is written
2398 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2400 if (!strlen(line)) return;
2401 if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2402 fStatisticsMsg += line;
2405 //______________________________________________________________________________
2406 void AliAnalysisManager::WriteStatisticsMsg(Int_t)
2408 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2409 static Bool_t done = kFALSE;
2412 if (!fStatistics) return;
2414 AddStatisticsMsg(Form("Number of input events: %lld",fStatistics->GetNinput()));
2415 AddStatisticsMsg(Form("Number of processed events: %lld",fStatistics->GetNprocessed()));
2416 AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2417 AddStatisticsMsg(Form("Number of accepted events for mask %s: %lld", AliAnalysisStatistics::GetMaskAsString(fStatistics->GetOfflineMask()), fStatistics->GetNaccepted()));
2418 out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2419 fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2420 fStatistics->GetNaccepted()), ios::out);
2421 out << fStatisticsMsg << endl;
2425 //______________________________________________________________________________
2426 const char* AliAnalysisManager::GetOADBPath()
2428 // returns the path of the OADB
2429 // this static function just depends on environment variables
2431 static TString oadbPath;
2433 if (gSystem->Getenv("OADB_PATH"))
2434 oadbPath = gSystem->Getenv("OADB_PATH");
2435 else if (gSystem->Getenv("ALICE_ROOT"))
2436 oadbPath.Form("%s/OADB", gSystem->Getenv("ALICE_ROOT"));
2438 ::Fatal("AliAnalysisManager::GetOADBPath", "Cannot figure out AODB path. Define ALICE_ROOT or OADB_PATH!");