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>
38 #include <TMethodCall.h>
43 #include <TStopwatch.h>
46 #include "AliAnalysisSelector.h"
47 #include "AliAnalysisGrid.h"
48 #include "AliAnalysisTask.h"
49 #include "AliAnalysisDataContainer.h"
50 #include "AliAnalysisDataSlot.h"
51 #include "AliVEventHandler.h"
52 #include "AliVEventPool.h"
53 #include "AliSysInfo.h"
54 #include "AliAnalysisStatistics.h"
56 ClassImp(AliAnalysisManager)
58 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
59 TString AliAnalysisManager::fgCommonFileName = "";
60 Int_t AliAnalysisManager::fPBUpdateFreq = 1;
62 //______________________________________________________________________________
63 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
66 fInputEventHandler(NULL),
67 fOutputEventHandler(NULL),
68 fMCtruthEventHandler(NULL),
72 fMode(kLocalAnalysis),
76 fSpecialOutputLocation(""),
89 fAutoBranchHandling(kTRUE),
99 // Default constructor.
100 fgAnalysisManager = this;
101 fgCommonFileName = "AnalysisResults.root";
102 if (TClass::IsCallingNew() != TClass::kDummyNew) {
103 fTasks = new TObjArray();
104 fTopTasks = new TObjArray();
105 fZombies = new TObjArray();
106 fContainers = new TObjArray();
107 fInputs = new TObjArray();
108 fOutputs = new TObjArray();
109 fParamCont = new TObjArray();
110 fGlobals = new TMap();
115 //______________________________________________________________________________
116 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
119 fInputEventHandler(NULL),
120 fOutputEventHandler(NULL),
121 fMCtruthEventHandler(NULL),
126 fInitOK(other.fInitOK),
127 fIsRemote(other.fIsRemote),
128 fDebug(other.fDebug),
129 fSpecialOutputLocation(""),
142 fAutoBranchHandling(other.fAutoBranchHandling),
145 fNcalls(other.fNcalls),
146 fMaxEntries(other.fMaxEntries),
147 fStatisticsMsg(other.fStatisticsMsg),
148 fRequestedBranches(other.fRequestedBranches),
149 fStatistics(other.fStatistics),
150 fGlobals(other.fGlobals)
153 fTasks = new TObjArray(*other.fTasks);
154 fTopTasks = new TObjArray(*other.fTopTasks);
155 fZombies = new TObjArray(*other.fZombies);
156 fContainers = new TObjArray(*other.fContainers);
157 fInputs = new TObjArray(*other.fInputs);
158 fOutputs = new TObjArray(*other.fOutputs);
159 fParamCont = new TObjArray(*other.fParamCont);
160 fgCommonFileName = "AnalysisResults.root";
161 fgAnalysisManager = this;
164 //______________________________________________________________________________
165 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
168 if (&other != this) {
169 TNamed::operator=(other);
170 fInputEventHandler = other.fInputEventHandler;
171 fOutputEventHandler = other.fOutputEventHandler;
172 fMCtruthEventHandler = other.fMCtruthEventHandler;
173 fEventPool = other.fEventPool;
176 fNSysInfo = other.fNSysInfo;
178 fInitOK = other.fInitOK;
179 fIsRemote = other.fIsRemote;
180 fDebug = other.fDebug;
181 fTasks = new TObjArray(*other.fTasks);
182 fTopTasks = new TObjArray(*other.fTopTasks);
183 fZombies = new TObjArray(*other.fZombies);
184 fContainers = new TObjArray(*other.fContainers);
185 fInputs = new TObjArray(*other.fInputs);
186 fOutputs = new TObjArray(*other.fOutputs);
187 fParamCont = new TObjArray(*other.fParamCont);
189 fCommonOutput = NULL;
192 fExtraFiles = other.fExtraFiles;
193 fgCommonFileName = "AnalysisResults.root";
194 fgAnalysisManager = this;
195 fAutoBranchHandling = other.fAutoBranchHandling;
196 fTable.Clear("nodelete");
197 fRunFromPath = other.fRunFromPath;
198 fNcalls = other. fNcalls;
199 fMaxEntries = other.fMaxEntries;
200 fStatisticsMsg = other.fStatisticsMsg;
201 fRequestedBranches = other.fRequestedBranches;
202 fStatistics = other.fStatistics;
203 fGlobals = new TMap();
208 //______________________________________________________________________________
209 AliAnalysisManager::~AliAnalysisManager()
212 if (fTasks) {fTasks->Delete(); delete fTasks;}
213 if (fTopTasks) delete fTopTasks;
214 if (fZombies) delete fZombies;
215 if (fContainers) {fContainers->Delete(); delete fContainers;}
216 if (fInputs) delete fInputs;
217 if (fOutputs) delete fOutputs;
218 if (fParamCont) delete fParamCont;
219 if (fGridHandler) delete fGridHandler;
220 if (fInputEventHandler) delete fInputEventHandler;
221 if (fOutputEventHandler) delete fOutputEventHandler;
222 if (fMCtruthEventHandler) delete fMCtruthEventHandler;
223 if (fEventPool) delete fEventPool;
224 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
225 if (fGlobals) {fGlobals->DeleteAll(); delete fGlobals;}
228 //______________________________________________________________________________
229 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
231 // Read one entry of the tree or a whole branch.
232 fCurrentEntry = entry;
233 if (!fAutoBranchHandling)
235 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : -1;
238 //______________________________________________________________________________
239 Int_t AliAnalysisManager::GetRunFromAlienPath(const char *path)
241 // Attempt to extract run number from input data path. Works only for paths to
242 // alice data in alien.
243 // sim: /alice/sim/<production>/run_no/...
244 // data: /alice/data/year/period/000run_no/... (ESD or AOD)
245 TString type = "unknown";
247 if (s.Contains("/alice/data")) type = "real";
248 else if (s.Contains("/alice/sim")) type = "simulated";
251 ind1 = s.Index("/00");
253 ind2 = s.Index("/",ind1+1);
254 if (ind2-ind1>8) srun = s(ind1+1, ind2-ind1-1);
257 ind1 = s.Index("/LHC");
259 ind1 = s.Index("/",ind1+1);
261 ind2 = s.Index("/",ind1+1);
262 if (ind2>0) srun = s(ind1+1, ind2-ind1-1);
266 Int_t run = srun.Atoi();
267 if (run>0) printf("=== GetRunFromAlienPath: run %d of %s data ===\n", run, type.Data());
271 //______________________________________________________________________________
272 Bool_t AliAnalysisManager::Init(TTree *tree)
274 // The Init() function is called when the selector needs to initialize
275 // a new tree or chain. Typically here the branch addresses of the tree
276 // will be set. It is normaly not necessary to make changes to the
277 // generated code, but the routine can be extended by the user if needed.
278 // Init() will be called many times when running with PROOF.
279 Bool_t init = kFALSE;
280 if (!tree) return kFALSE; // Should not happen - protected in selector caller
282 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
284 // Call InitTree of EventHandler
285 if (fOutputEventHandler) {
286 if (fMode == kProofAnalysis) {
287 init = fOutputEventHandler->Init(0x0, "proof");
289 init = fOutputEventHandler->Init(0x0, "local");
292 Error("Init", "Output event handler failed to initialize");
297 if (fInputEventHandler) {
298 if (fMode == kProofAnalysis) {
299 init = fInputEventHandler->Init(tree, "proof");
301 init = fInputEventHandler->Init(tree, "local");
304 Error("Init", "Input event handler failed to initialize tree");
308 // If no input event handler we need to get the tree once
310 if(!tree->GetTree()) {
311 Long64_t readEntry = tree->LoadTree(0);
312 if (readEntry == -2) {
313 Error("Init", "Input tree has no entry. Exiting");
319 if (fMCtruthEventHandler) {
320 if (fMode == kProofAnalysis) {
321 init = fMCtruthEventHandler->Init(0x0, "proof");
323 init = fMCtruthEventHandler->Init(0x0, "local");
326 Error("Init", "MC event handler failed to initialize");
331 if (!fInitOK) InitAnalysis();
332 if (!fInitOK) return kFALSE;
335 AliAnalysisDataContainer *top = fCommonInput;
336 if (!top) top = (AliAnalysisDataContainer*)fInputs->At(0);
338 Error("Init","No top input container !");
342 CheckBranches(kFALSE);
344 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
349 //______________________________________________________________________________
350 void AliAnalysisManager::SlaveBegin(TTree *tree)
352 // The SlaveBegin() function is called after the Begin() function.
353 // When running with PROOF SlaveBegin() is called on each slave server.
354 // The tree argument is deprecated (on PROOF 0 is passed).
355 if (fDebug > 1) printf("->AliAnalysisManager::SlaveBegin()\n");
356 if (!CheckTasks()) Fatal("SlaveBegin", "Not all needed libraries were loaded");
357 static Bool_t isCalled = kFALSE;
358 Bool_t init = kFALSE;
359 Bool_t initOK = kTRUE;
361 TDirectory *curdir = gDirectory;
362 // Call SlaveBegin only once in case of mixing
363 if (isCalled && fMode==kMixingAnalysis) return;
365 // Call Init of EventHandler
366 if (fOutputEventHandler) {
367 if (fMode == kProofAnalysis) {
368 // Merging AOD's in PROOF via TProofOutputFile
369 if (fDebug > 1) printf(" Initializing AOD output file %s...\n", fOutputEventHandler->GetOutputFileName());
370 init = fOutputEventHandler->Init("proof");
371 if (!init) msg = "Failed to initialize output handler on worker";
373 init = fOutputEventHandler->Init("local");
374 if (!init) msg = "Failed to initialize output handler";
377 if (!fSelector) Error("SlaveBegin", "Selector not set");
378 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
381 if (fInputEventHandler) {
382 fInputEventHandler->SetInputTree(tree);
383 if (fMode == kProofAnalysis) {
384 init = fInputEventHandler->Init("proof");
385 if (!init) msg = "Failed to initialize input handler on worker";
387 init = fInputEventHandler->Init("local");
388 if (!init) msg = "Failed to initialize input handler";
391 if (!fSelector) Error("SlaveBegin", "Selector not set");
392 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
395 if (fMCtruthEventHandler) {
396 if (fMode == kProofAnalysis) {
397 init = fMCtruthEventHandler->Init("proof");
398 if (!init) msg = "Failed to initialize MC handler on worker";
400 init = fMCtruthEventHandler->Init("local");
401 if (!init) msg = "Failed to initialize MC handler";
404 if (!fSelector) Error("SlaveBegin", "Selector not set");
405 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
407 if (curdir) curdir->cd();
411 AliAnalysisTask *task;
412 // Call CreateOutputObjects for all tasks
413 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
414 Bool_t dirStatus = TH1::AddDirectoryStatus();
416 while ((task=(AliAnalysisTask*)next())) {
418 // Start with memory as current dir and make sure by default histograms do not get attached to files.
419 TH1::AddDirectory(kFALSE);
420 task->CreateOutputObjects();
421 if (!task->CheckPostData()) {
422 Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
423 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
424 ####### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ##########", task->GetName(), task->ClassName());
426 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
429 TH1::AddDirectory(dirStatus);
430 if (curdir) curdir->cd();
431 if (fDebug > 1) printf("<-AliAnalysisManager::SlaveBegin()\n");
434 //______________________________________________________________________________
435 Bool_t AliAnalysisManager::Notify()
437 // The Notify() function is called when a new file is opened. This
438 // can be either for a new TTree in a TChain or when when a new TTree
439 // is started when using PROOF. It is normaly not necessary to make changes
440 // to the generated code, but the routine can be extended by the
441 // user if needed. The return value is currently not used.
442 if (!fTree) return kFALSE;
443 if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) return kFALSE;
445 fTable.Clear("nodelete"); // clearing the hash table may not be needed -> C.L.
446 if (fMode == kProofAnalysis) fIsRemote = kTRUE;
448 TFile *curfile = fTree->GetCurrentFile();
450 Error("Notify","No current file");
454 if (fDebug > 1) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
455 Int_t run = AliAnalysisManager::GetRunFromAlienPath(curfile->GetName());
456 if (run && (run != fRunFromPath)) {
458 if (fDebug > 1) printf(" ### run found from path: %d\n", run);
461 AliAnalysisTask *task;
463 // Call Notify of the event handlers
464 if (fInputEventHandler) {
465 fInputEventHandler->Notify(curfile->GetName());
468 if (fOutputEventHandler) {
469 fOutputEventHandler->Notify(curfile->GetName());
472 if (fMCtruthEventHandler) {
473 fMCtruthEventHandler->Notify(curfile->GetName());
476 // Call Notify for all tasks
477 while ((task=(AliAnalysisTask*)next()))
480 if (fDebug > 1) printf("<-AliAnalysisManager::Notify()\n");
484 //______________________________________________________________________________
485 Bool_t AliAnalysisManager::Process(Long64_t)
487 // The Process() function is called for each entry in the tree (or possibly
488 // keyed object in the case of PROOF) to be processed. The entry argument
489 // specifies which entry in the currently loaded tree is to be processed.
490 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
491 // to read either all or the required parts of the data. When processing
492 // keyed objects with PROOF, the object is already loaded and is available
493 // via the fObject pointer.
495 // This function should contain the "body" of the analysis. It can contain
496 // simple or elaborate selection criteria, run algorithms on the data
497 // of the event and typically fill histograms.
499 // WARNING when a selector is used with a TChain, you must use
500 // the pointer to the current TTree to call GetEntry(entry).
501 // The entry is always the local entry number in the current tree.
502 // Assuming that fChain is the pointer to the TChain being processed,
503 // use fChain->GetTree()->GetEntry(entry).
505 // This method is obsolete. ExecAnalysis is called instead.
509 //______________________________________________________________________________
510 void AliAnalysisManager::PackOutput(TList *target)
512 // Pack all output data containers in the output list. Called at SlaveTerminate
513 // stage in PROOF case for each slave.
514 if (fDebug > 1) printf("->AliAnalysisManager::PackOutput()\n");
516 Error("PackOutput", "No target. Exiting.");
519 TDirectory *cdir = gDirectory;
521 if (fInputEventHandler) fInputEventHandler ->Terminate();
522 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
523 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
526 // Call FinishTaskOutput() for each event loop task (not called for
527 // post-event loop tasks - use Terminate() fo those)
528 TIter nexttask(fTasks);
529 AliAnalysisTask *task;
530 while ((task=(AliAnalysisTask*)nexttask())) {
531 if (!task->IsPostEventLoop()) {
532 if (fDebug > 1) printf("->FinishTaskOutput: task %s\n", task->GetName());
533 task->FinishTaskOutput();
535 if (fDebug > 1) printf("<-FinishTaskOutput: task %s\n", task->GetName());
538 // Write statistics message on the workers.
539 if (fStatistics) WriteStatisticsMsg(fNcalls);
541 if (fMode == kProofAnalysis) {
542 TIter next(fOutputs);
543 AliAnalysisDataContainer *output;
544 Bool_t isManagedByHandler = kFALSE;
547 while ((output=(AliAnalysisDataContainer*)next())) {
548 // Do not consider outputs of post event loop tasks
549 isManagedByHandler = kFALSE;
550 if (output->GetProducer() && output->GetProducer()->IsPostEventLoop()) continue;
551 const char *filename = output->GetFileName();
552 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
553 isManagedByHandler = kTRUE;
554 printf("#### Handler output. Extra: %s\n", fExtraFiles.Data());
555 filename = fOutputEventHandler->GetOutputFileName();
557 // Check if data was posted to this container. If not, issue an error.
558 if (!output->GetData() && !isManagedByHandler) {
559 Error("PackOutput", "No data for output container %s. Forgot to PostData ?", output->GetName());
562 if (!output->IsSpecialOutput()) {
564 if (strlen(filename) && !isManagedByHandler) {
565 // Backup current folder
566 TDirectory *opwd = gDirectory;
567 // File resident outputs.
568 // Check first if the file exists.
569 TString openoption = "RECREATE";
570 Bool_t firsttime = kTRUE;
571 if (filestmp.FindObject(output->GetFileName())) {
574 filestmp.Add(new TNamed(output->GetFileName(),""));
576 if (!gSystem->AccessPathName(output->GetFileName()) && !firsttime) openoption = "UPDATE";
577 // TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
578 // Save data to file, then close.
579 if (output->GetData()->InheritsFrom(TCollection::Class())) {
580 // If data is a collection, we set the name of the collection
581 // as the one of the container and we save as a single key.
582 TCollection *coll = (TCollection*)output->GetData();
583 coll->SetName(output->GetName());
584 // coll->Write(output->GetName(), TObject::kSingleKey);
586 if (output->GetData()->InheritsFrom(TTree::Class())) {
587 TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
588 // Save data to file, then close.
589 TTree *tree = (TTree*)output->GetData();
590 // Check if tree is in memory
591 if (tree->GetDirectory()==gROOT) tree->SetDirectory(gDirectory);
595 // output->GetData()->Write();
598 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
600 // printf(" file %s listing content:\n", filename);
603 // Clear file list to release object ownership to user.
606 output->SetFile(NULL);
607 // Restore current directory
608 if (opwd) opwd->cd();
610 // Memory-resident outputs
611 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
613 AliAnalysisDataWrapper *wrap = 0;
614 if (isManagedByHandler) {
615 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
616 wrap->SetName(output->GetName());
618 else wrap =output->ExportData();
619 // Output wrappers must NOT delete data after merging - the user owns them
620 wrap->SetDeleteData(kFALSE);
623 // Special outputs. The file must be opened and connected to the container.
624 TDirectory *opwd = gDirectory;
625 TFile *file = output->GetFile();
627 AliAnalysisTask *producer = output->GetProducer();
629 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
630 output->GetFileName(), output->GetName(), producer->ClassName());
633 TString outFilename = file->GetName();
634 if (fDebug > 1) printf("PackOutput %s: special output\n", output->GetName());
635 if (isManagedByHandler) {
636 // Terminate IO for files managed by the output handler
637 // file->Write() moved to AOD handler (A.G. 11.01.10)
638 // if (file) file->Write();
639 if (file && fDebug > 2) {
640 printf(" handled file %s listing content:\n", file->GetName());
643 fOutputEventHandler->TerminateIO();
646 // Release object ownership to users after writing data to file
647 if (output->GetData()->InheritsFrom(TCollection::Class())) {
648 // If data is a collection, we set the name of the collection
649 // as the one of the container and we save as a single key.
650 TCollection *coll = (TCollection*)output->GetData();
651 coll->SetName(output->GetName());
652 coll->Write(output->GetName(), TObject::kSingleKey);
654 if (output->GetData()->InheritsFrom(TTree::Class())) {
655 TTree *tree = (TTree*)output->GetData();
656 tree->SetDirectory(file);
659 output->GetData()->Write();
663 printf(" file %s listing content:\n", output->GetFileName());
666 // Clear file list to release object ownership to user.
669 output->SetFile(NULL);
671 // Restore current directory
672 if (opwd) opwd->cd();
673 // Check if a special output location was provided or the output files have to be merged
674 if (strlen(fSpecialOutputLocation.Data())) {
675 TString remote = fSpecialOutputLocation;
677 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
678 if (remote.BeginsWith("alien:")) {
679 gROOT->ProcessLine("TGrid::Connect(\"alien:\", gProofServ->GetUser());");
680 remote += outFilename;
681 remote.ReplaceAll(".root", Form("_%d.root", gid));
683 remote += Form("%s_%d_", gSystem->HostName(), gid);
684 remote += outFilename;
687 Info("PackOutput", "Output file for container %s to be copied \n at: %s. No merging.",
688 output->GetName(), remote.Data());
689 TFile::Cp ( outFilename.Data(), remote.Data() );
690 // Copy extra outputs
691 if (fExtraFiles.Length() && isManagedByHandler) {
692 TObjArray *arr = fExtraFiles.Tokenize(" ");
694 TIter nextfilename(arr);
695 while ((os=(TObjString*)nextfilename())) {
696 outFilename = os->GetString();
697 remote = fSpecialOutputLocation;
699 if (remote.BeginsWith("alien://")) {
700 remote += outFilename;
701 remote.ReplaceAll(".root", Form("_%d.root", gid));
703 remote += Form("%s_%d_", gSystem->HostName(), gid);
704 remote += outFilename;
707 Info("PackOutput", "Extra AOD file %s to be copied \n at: %s. No merging.",
708 outFilename.Data(), remote.Data());
709 TFile::Cp ( outFilename.Data(), remote.Data() );
714 // No special location specified-> use TProofOutputFile as merging utility
715 // The file at this output slot must be opened in CreateOutputObjects
716 if (fDebug > 1) printf(" File for container %s to be merged via file merger...\n", output->GetName());
722 if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
725 //______________________________________________________________________________
726 void AliAnalysisManager::ImportWrappers(TList *source)
728 // Import data in output containers from wrappers coming in source.
729 if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
730 TIter next(fOutputs);
731 AliAnalysisDataContainer *cont;
732 AliAnalysisDataWrapper *wrap;
734 Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
735 TDirectory *cdir = gDirectory;
736 while ((cont=(AliAnalysisDataContainer*)next())) {
738 if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
739 if (cont->IsRegisterDataset()) continue;
740 const char *filename = cont->GetFileName();
741 Bool_t isManagedByHandler = kFALSE;
742 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
743 isManagedByHandler = kTRUE;
744 filename = fOutputEventHandler->GetOutputFileName();
746 if (cont->IsSpecialOutput() || inGrid) {
747 if (strlen(fSpecialOutputLocation.Data())) continue;
748 // Copy merged file from PROOF scratch space.
749 // In case of grid the files are already in the current directory.
751 if (isManagedByHandler && fExtraFiles.Length()) {
752 // Copy extra registered dAOD files.
753 TObjArray *arr = fExtraFiles.Tokenize(" ");
755 TIter nextfilename(arr);
756 while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
759 if (!GetFileFromWrapper(filename, source)) continue;
761 // Normally we should connect data from the copied file to the
762 // corresponding output container, but it is not obvious how to do this
763 // automatically if several objects in file...
764 TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
765 if (!f) f = TFile::Open(filename, "READ");
767 Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
772 // Cd to the directory pointed by the container
773 TString folder = cont->GetFolderName();
774 if (!folder.IsNull()) f->cd(folder);
775 // Try to fetch first an object having the container name.
776 obj = gDirectory->Get(cont->GetName());
778 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",
779 cont->GetType()->GetName(), cont->GetName(), filename, cont->GetFolderName(), cont->GetName());
782 wrap = new AliAnalysisDataWrapper(obj);
783 wrap->SetDeleteData(kFALSE);
785 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
787 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
792 printf(" Importing data for container %s\n", cont->GetName());
793 if (strlen(filename)) printf(" -> file %s\n", filename);
796 cont->ImportData(wrap);
798 if (cdir) cdir->cd();
799 if (fDebug > 1) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
802 //______________________________________________________________________________
803 void AliAnalysisManager::UnpackOutput(TList *source)
805 // Called by AliAnalysisSelector::Terminate only on the client.
806 if (fDebug > 1) printf("->AliAnalysisManager::UnpackOutput()\n");
808 Error("UnpackOutput", "No target. Exiting.");
811 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
813 if (fMode == kProofAnalysis) ImportWrappers(source);
815 TIter next(fOutputs);
816 AliAnalysisDataContainer *output;
817 while ((output=(AliAnalysisDataContainer*)next())) {
818 if (!output->GetData()) continue;
819 // Check if there are client tasks that run post event loop
820 if (output->HasConsumers()) {
821 // Disable event loop semaphore
822 output->SetPostEventLoop(kTRUE);
823 TObjArray *list = output->GetConsumers();
824 Int_t ncons = list->GetEntriesFast();
825 for (Int_t i=0; i<ncons; i++) {
826 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
827 task->CheckNotify(kTRUE);
828 // If task is active, execute it
829 if (task->IsPostEventLoop() && task->IsActive()) {
830 if (fDebug > 1) printf("== Executing post event loop task %s\n", task->GetName());
836 if (fDebug > 1) printf("<-AliAnalysisManager::UnpackOutput()\n");
839 //______________________________________________________________________________
840 void AliAnalysisManager::Terminate()
842 // The Terminate() function is the last function to be called during
843 // a query. It always runs on the client, it can be used to present
844 // the results graphically.
845 if (fDebug > 1) printf("->AliAnalysisManager::Terminate()\n");
846 TDirectory *cdir = gDirectory;
848 AliAnalysisTask *task;
849 AliAnalysisDataContainer *output;
852 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
853 // Call Terminate() for tasks
855 while (!IsSkipTerminate() && (task=(AliAnalysisTask*)next())) {
856 // Save all the canvases produced by the Terminate
857 TString pictname = Form("%s_%s", task->GetName(), task->ClassName());
861 AliSysInfo::AddStamp(Form("%s_TERMINATE",task->ClassName()),0, itask, 2);
863 if (TObject::TestBit(kSaveCanvases)) {
864 if (!gROOT->IsBatch()) {
865 if (fDebug>1) printf("Waiting 5 sec for %s::Terminate() to finish drawing ...\n", task->ClassName());
867 while (timer.CpuTime()<5) {
869 gSystem->ProcessEvents();
872 Int_t iend = gROOT->GetListOfCanvases()->GetEntries();
873 if (iend==0) continue;
875 for (Int_t ipict=0; ipict<iend; ipict++) {
876 canvas = (TCanvas*)gROOT->GetListOfCanvases()->At(ipict);
877 if (!canvas) continue;
878 canvas->SaveAs(Form("%s_%02d.gif", pictname.Data(),ipict));
880 gROOT->GetListOfCanvases()->Delete();
884 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
885 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
886 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
888 TObjArray *allOutputs = new TObjArray();
890 for (icont=0; icont<fOutputs->GetEntriesFast(); icont++) allOutputs->Add(fOutputs->At(icont));
891 if (!IsSkipTerminate())
892 for (icont=0; icont<fParamCont->GetEntriesFast(); icont++) allOutputs->Add(fParamCont->At(icont));
893 TIter next1(allOutputs);
894 TString handlerFile = "";
895 TString extraOutputs = "";
896 if (fOutputEventHandler) {
897 handlerFile = fOutputEventHandler->GetOutputFileName();
898 extraOutputs = fOutputEventHandler->GetExtraOutputs();
902 while ((output=(AliAnalysisDataContainer*)next1())) {
903 // Special outputs or grid files have the files already closed and written.
905 if (fMode == kGridAnalysis && icont<=fOutputs->GetEntriesFast()) continue;
906 if (fMode == kProofAnalysis) {
907 if (output->IsSpecialOutput() || output->IsRegisterDataset()) continue;
909 const char *filename = output->GetFileName();
910 TString openoption = "RECREATE";
911 if (!(strcmp(filename, "default"))) continue;
912 if (!strlen(filename)) continue;
913 if (!output->GetData()) continue;
914 TDirectory *opwd = gDirectory;
915 TFile *file = output->GetFile();
916 if (!file) file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
918 //if (handlerFile == filename && !gSystem->AccessPathName(filename)) openoption = "UPDATE";
919 Bool_t firsttime = kTRUE;
920 if (filestmp.FindObject(filename) || extraOutputs.Contains(filename)) {
923 filestmp.Add(new TNamed(filename,""));
925 if (!gSystem->AccessPathName(filename) && !firsttime) openoption = "UPDATE";
926 if (fDebug>1) printf("Opening file: %s option=%s\n",filename, openoption.Data());
927 file = new TFile(filename, openoption);
929 if (fDebug>1) printf("File <%s> already opened with option: <%s> \n", filename, file->GetOption());
930 openoption = file->GetOption();
931 if (openoption == "READ") {
932 if (fDebug>1) printf("...reopening in UPDATE mode\n");
933 file->ReOpen("UPDATE");
936 if (file->IsZombie()) {
937 Error("Terminate", "Cannot open output file %s", filename);
940 output->SetFile(file);
942 // Check for a folder request
943 TString dir = output->GetFolderName();
945 if (!file->GetDirectory(dir)) file->mkdir(dir);
948 if (fDebug > 1) printf("...writing container %s to file %s:%s\n", output->GetName(), file->GetName(), output->GetFolderName());
949 if (output->GetData()->InheritsFrom(TCollection::Class())) {
950 // If data is a collection, we set the name of the collection
951 // as the one of the container and we save as a single key.
952 TCollection *coll = (TCollection*)output->GetData();
953 coll->SetName(output->GetName());
954 coll->Write(output->GetName(), TObject::kSingleKey);
956 if (output->GetData()->InheritsFrom(TTree::Class())) {
957 TTree *tree = (TTree*)output->GetData();
958 tree->SetDirectory(gDirectory);
961 output->GetData()->Write();
964 if (opwd) opwd->cd();
969 while ((output=(AliAnalysisDataContainer*)next1())) {
970 // Close all files at output
971 TDirectory *opwd = gDirectory;
972 if (output->GetFile()) {
973 // Clear file list to release object ownership to user.
974 // output->GetFile()->Clear();
975 output->GetFile()->Close();
976 // Copy merged outputs in alien if requested
977 if (fSpecialOutputLocation.BeginsWith("alien://")) {
978 if (copiedFiles.Contains(output->GetFile()->GetName())) {
979 if (opwd) opwd->cd();
980 output->SetFile(NULL);
983 Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data());
984 gROOT->ProcessLine("if (!gGrid) TGrid::Connect(\"alien:\");");
985 TFile::Cp(output->GetFile()->GetName(),
986 Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
987 copiedFiles += output->GetFile()->GetName();
989 output->SetFile(NULL);
991 if (opwd) opwd->cd();
994 //Write statistics information on the client
995 if (fStatistics) WriteStatisticsMsg(fNcalls);
997 TDirectory *crtdir = gDirectory;
998 TFile f("syswatch.root", "RECREATE");
1001 if (!f.IsZombie()) {
1002 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
1003 tree->SetName("syswatch");
1004 tree->SetMarkerStyle(kCircle);
1005 tree->SetMarkerColor(kBlue);
1006 tree->SetMarkerSize(0.5);
1007 if (!gROOT->IsBatch()) {
1008 tree->SetAlias("event", "id0");
1009 tree->SetAlias("task", "id1");
1010 tree->SetAlias("stage", "id2");
1011 // Already defined aliases
1012 // tree->SetAlias("deltaT","stampSec-stampOldSec");
1013 // tree->SetAlias("T","stampSec-first");
1014 // tree->SetAlias("deltaVM","(pI.fMemVirtual-pIOld.fMemVirtual)");
1015 // tree->SetAlias("VM","pI.fMemVirtual");
1016 TCanvas *canvas = new TCanvas("SysInfo","SysInfo",10,10,1200,800);
1017 Int_t npads = 1 /*COO plot for all tasks*/ +
1018 fTopTasks->GetEntries() /*Exec plot per task*/ +
1019 1 /*Terminate plot for all tasks*/ +
1022 Int_t iopt = (Int_t)TMath::Sqrt((Double_t)npads);
1023 if (npads<iopt*(iopt+1))
1024 canvas->Divide(iopt, iopt+1, 0.01, 0.01);
1026 canvas->Divide(iopt+1, iopt+1, 0.01, 0.01);
1028 // draw the plot of deltaVM for Exec for each task
1029 for (itask=0; itask<fTopTasks->GetEntriesFast(); itask++) {
1030 task = (AliAnalysisTask*)fTopTasks->At(itask);
1032 cut = Form("task==%d && stage==1", itask);
1033 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1034 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1036 hist->SetTitle(Form("%s: Exec dVM[MB]/event", task->GetName()));
1037 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1040 // Draw the plot of deltaVM for CreateOutputObjects for all tasks
1042 tree->SetMarkerStyle(kFullTriangleUp);
1043 tree->SetMarkerColor(kRed);
1044 tree->SetMarkerSize(0.8);
1045 cut = "task>=0 && task<1000 && stage==0";
1046 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1047 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1049 hist->SetTitle("Memory in CreateOutputObjects()");
1050 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1051 hist->GetXaxis()->SetTitle("task");
1053 // draw the plot of deltaVM for Terminate for all tasks
1055 tree->SetMarkerStyle(kOpenSquare);
1056 tree->SetMarkerColor(kMagenta);
1057 cut = "task>=0 && task<1000 && stage==2";
1058 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1059 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1061 hist->SetTitle("Memory in Terminate()");
1062 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1063 hist->GetXaxis()->SetTitle("task");
1067 tree->SetMarkerStyle(kFullCircle);
1068 tree->SetMarkerColor(kGreen);
1069 cut = Form("task==%d && stage==1",fTopTasks->GetEntriesFast()-1);
1070 tree->Draw("VM:event",cut,"", 1234567890, 0);
1071 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1073 hist->SetTitle("Virtual memory");
1074 hist->GetYaxis()->SetTitle("VM [MB]");
1078 tree->SetMarkerStyle(kCircle);
1079 tree->SetMarkerColor(kBlue);
1080 tree->SetMarkerSize(0.5);
1085 if (crtdir) crtdir->cd();
1087 // Validate the output files
1088 if (ValidateOutputFiles() && fIsRemote && fMode!=kProofAnalysis) {
1090 out.open("outputs_valid", ios::out);
1094 if (fDebug > 1) printf("<-AliAnalysisManager::Terminate()\n");
1096 //______________________________________________________________________________
1097 void AliAnalysisManager::ProfileTask(Int_t itop, const char *option) const
1099 // Profiles the task having the itop index in the list of top (first level) tasks.
1100 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->At(itop);
1102 Error("ProfileTask", "There are only %d top tasks in the manager", fTopTasks->GetEntries());
1105 ProfileTask(task->GetName(), option);
1108 //______________________________________________________________________________
1109 void AliAnalysisManager::ProfileTask(const char *name, const char */*option*/) const
1111 // Profile a managed task after the execution of the analysis in case NSysInfo
1113 if (gSystem->AccessPathName("syswatch.root")) {
1114 Error("ProfileTask", "No file syswatch.root found in the current directory");
1117 if (gROOT->IsBatch()) return;
1118 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->FindObject(name);
1120 Error("ProfileTask", "No top task named %s known by the manager.", name);
1123 Int_t itop = fTopTasks->IndexOf(task);
1124 Int_t itask = fTasks->IndexOf(task);
1125 // Create canvas with 2 pads: first draw COO + Terminate, second Exec
1126 TDirectory *cdir = gDirectory;
1127 TFile f("syswatch.root");
1128 TTree *tree = (TTree*)f.Get("syswatch");
1130 Error("ProfileTask", "No tree named <syswatch> found in file syswatch.root");
1133 if (fDebug > 1) printf("=== Profiling task %s (class %s)\n", name, task->ClassName());
1134 TCanvas *canvas = new TCanvas(Form("profile_%d",itop),Form("Profile of task %s (class %s)",name,task->ClassName()),10,10,800,600);
1135 canvas->Divide(2, 2, 0.01, 0.01);
1139 // VM profile for COO and Terminate methods
1141 cut = Form("task==%d && (stage==0 || stage==2)",itask);
1142 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1143 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1145 hist->SetTitle("Alocated VM[MB] for COO and Terminate");
1146 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1147 hist->GetXaxis()->SetTitle("method");
1149 // CPU profile per event
1151 cut = Form("task==%d && stage==1",itop);
1152 tree->Draw("deltaT:event",cut,"", 1234567890, 0);
1153 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1155 hist->SetTitle("Execution time per event");
1156 hist->GetYaxis()->SetTitle("CPU/event [s]");
1158 // VM profile for Exec
1160 cut = Form("task==%d && stage==1",itop);
1161 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1162 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1164 hist->SetTitle("Alocated VM[MB] per event");
1165 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1170 if (cdir) cdir->cd();
1173 //______________________________________________________________________________
1174 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1176 // Adds a user task to the global list of tasks.
1178 Error("AddTask", "Cannot add task %s since InitAnalysis was already called", task->GetName());
1182 if (fTasks->FindObject(task)) {
1183 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1186 task->SetActive(kFALSE);
1190 //______________________________________________________________________________
1191 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1193 // Retreive task by name.
1194 if (!fTasks) return NULL;
1195 return (AliAnalysisTask*)fTasks->FindObject(name);
1198 //______________________________________________________________________________
1199 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
1200 TClass *datatype, EAliAnalysisContType type, const char *filename)
1202 // Create a data container of a certain type. Types can be:
1203 // kExchangeContainer = 0, used to exchange data between tasks
1204 // kInputContainer = 1, used to store input data
1205 // kOutputContainer = 2, used for writing result to a file
1206 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1207 // the output object to a folder inside the output file
1208 if (fContainers->FindObject(name)) {
1209 Error("CreateContainer","A container named %s already defined !",name);
1212 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1213 fContainers->Add(cont);
1215 case kInputContainer:
1218 case kOutputContainer:
1219 fOutputs->Add(cont);
1220 if (filename && strlen(filename)) {
1221 cont->SetFileName(filename);
1222 cont->SetDataOwned(kFALSE); // data owned by the file
1225 case kParamContainer:
1226 fParamCont->Add(cont);
1227 if (filename && strlen(filename)) {
1228 cont->SetFileName(filename);
1229 cont->SetDataOwned(kFALSE); // data owned by the file
1232 case kExchangeContainer:
1238 //______________________________________________________________________________
1239 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1240 AliAnalysisDataContainer *cont)
1242 // Connect input of an existing task to a data container.
1244 Error("ConnectInput", "Task pointer is NULL");
1247 if (!fTasks->FindObject(task)) {
1249 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1251 Bool_t connected = task->ConnectInput(islot, cont);
1255 //______________________________________________________________________________
1256 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1257 AliAnalysisDataContainer *cont)
1259 // Connect output of an existing task to a data container.
1261 Error("ConnectOutput", "Task pointer is NULL");
1264 if (!fTasks->FindObject(task)) {
1266 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1268 Bool_t connected = task->ConnectOutput(islot, cont);
1272 //______________________________________________________________________________
1273 void AliAnalysisManager::CleanContainers()
1275 // Clean data from all containers that have already finished all client tasks.
1276 TIter next(fContainers);
1277 AliAnalysisDataContainer *cont;
1278 while ((cont=(AliAnalysisDataContainer *)next())) {
1279 if (cont->IsOwnedData() &&
1280 cont->IsDataReady() &&
1281 cont->ClientsExecuted()) cont->DeleteData();
1285 //______________________________________________________________________________
1286 Bool_t AliAnalysisManager::InitAnalysis()
1288 // Initialization of analysis chain of tasks. Should be called after all tasks
1289 // and data containers are properly connected
1290 // Reset flag and remove valid_outputs file if exists
1291 if (fInitOK) return kTRUE;
1292 if (!gSystem->AccessPathName("outputs_valid"))
1293 gSystem->Unlink("outputs_valid");
1294 // Check for top tasks (depending only on input data containers)
1295 if (!fTasks->First()) {
1296 Error("InitAnalysis", "Analysis has no tasks !");
1300 AliAnalysisTask *task;
1301 AliAnalysisDataContainer *cont;
1304 Bool_t iszombie = kFALSE;
1305 Bool_t istop = kTRUE;
1307 while ((task=(AliAnalysisTask*)next())) {
1310 Int_t ninputs = task->GetNinputs();
1311 for (i=0; i<ninputs; i++) {
1312 cont = task->GetInputSlot(i)->GetContainer();
1316 fZombies->Add(task);
1320 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
1321 i, task->GetName());
1323 if (iszombie) continue;
1324 // Check if cont is an input container
1325 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1326 // Connect to parent task
1330 fTopTasks->Add(task);
1334 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1337 // Check now if there are orphan tasks
1338 for (i=0; i<ntop; i++) {
1339 task = (AliAnalysisTask*)fTopTasks->At(i);
1344 while ((task=(AliAnalysisTask*)next())) {
1345 if (!task->IsUsed()) {
1347 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1350 // Check the task hierarchy (no parent task should depend on data provided
1351 // by a daughter task)
1352 for (i=0; i<ntop; i++) {
1353 task = (AliAnalysisTask*)fTopTasks->At(i);
1354 if (task->CheckCircularDeps()) {
1355 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1360 // Check that all containers feeding post-event loop tasks are in the outputs list
1361 TIter nextcont(fContainers); // loop over all containers
1362 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1363 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1364 if (cont->HasConsumers()) {
1365 // Check if one of the consumers is post event loop
1366 TIter nextconsumer(cont->GetConsumers());
1367 while ((task=(AliAnalysisTask*)nextconsumer())) {
1368 if (task->IsPostEventLoop()) {
1369 fOutputs->Add(cont);
1376 // Check if all special output containers have a file name provided
1377 TIter nextout(fOutputs);
1378 while ((cont=(AliAnalysisDataContainer*)nextout())) {
1379 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1380 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1384 // Initialize requested branch list if needed
1385 if (!fAutoBranchHandling) {
1387 while ((task=(AliAnalysisTask*)next())) {
1388 if (!task->HasBranches()) {
1389 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\"",
1390 task->GetName(), task->ClassName());
1393 if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1394 Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1397 TString taskbranches;
1398 task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1399 if (taskbranches.IsNull()) {
1400 Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1401 task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1404 AddBranches(taskbranches);
1411 //______________________________________________________________________________
1412 void AliAnalysisManager::AddBranches(const char *branches)
1414 // Add branches to the existing fRequestedBranches.
1415 TString br(branches);
1416 TObjArray *arr = br.Tokenize(",");
1419 while ((obj=next())) {
1420 if (!fRequestedBranches.Contains(obj->GetName())) {
1421 if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1422 fRequestedBranches += obj->GetName();
1428 //______________________________________________________________________________
1429 void AliAnalysisManager::CheckBranches(Bool_t load)
1431 // The method checks the input branches to be loaded during the analysis.
1432 if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;
1433 TObjArray *arr = fRequestedBranches.Tokenize(",");
1436 while ((obj=next())) {
1437 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1439 br = fTree->GetBranch(obj->GetName());
1441 Error("CheckBranches", "Could not find branch %s",obj->GetName());
1446 if (load && br->GetReadEntry()!=GetCurrentEntry()) br->GetEntry(GetCurrentEntry());
1451 //______________________________________________________________________________
1452 Bool_t AliAnalysisManager::CheckTasks() const
1454 // Check consistency of tasks.
1455 Int_t ntasks = fTasks->GetEntries();
1457 Error("CheckTasks", "No tasks connected to the manager. This may be due to forgetting to compile the task or to load their library.");
1460 // Get the pointer to AliAnalysisTaskSE::Class()
1461 TClass *badptr = (TClass*)gROOT->ProcessLine("AliAnalysisTaskSE::Class()");
1462 // Loop all tasks to check if their corresponding library was loaded
1465 while ((obj=next())) {
1466 if (obj->IsA() == badptr) {
1467 Error("CheckTasks", "##################\n \
1468 Class for task %s NOT loaded. You probably forgot to load the library for this task (or compile it dynamically).\n###########################\n",obj->GetName());
1475 //______________________________________________________________________________
1476 void AliAnalysisManager::PrintStatus(Option_t *option) const
1478 // Print task hierarchy.
1480 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1483 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1485 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1486 TIter next(fTopTasks);
1487 AliAnalysisTask *task;
1488 while ((task=(AliAnalysisTask*)next()))
1489 task->PrintTask(option);
1491 if (!fAutoBranchHandling && !fRequestedBranches.IsNull())
1492 printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1494 TString sopt(option);
1497 if (sopt.Contains("ALL"))
1499 if ( fOutputEventHandler )
1501 cout << TString('_',78) << endl;
1502 cout << "OutputEventHandler:" << endl;
1503 fOutputEventHandler->Print(" ");
1508 //______________________________________________________________________________
1509 void AliAnalysisManager::ResetAnalysis()
1511 // Reset all execution flags and clean containers.
1515 //______________________________________________________________________________
1516 void AliAnalysisManager::RunLocalInit()
1518 // Run LocalInit method for all tasks.
1519 TDirectory *cdir = gDirectory;
1520 if (IsTrainInitialized()) return;
1521 TIter nextTask(fTasks);
1522 AliAnalysisTask *task;
1523 while ((task=(AliAnalysisTask*)nextTask())) {
1528 TObject::SetBit(kTasksInitialized, kTRUE);
1531 //______________________________________________________________________________
1532 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1534 // Start analysis having a grid handler.
1535 if (!fGridHandler) {
1536 Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1537 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1541 return StartAnalysis(type, tree, nentries, firstentry);
1544 //______________________________________________________________________________
1545 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1547 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1548 // MIX. Process nentries starting from firstentry
1550 // Backup current directory and make sure gDirectory points to gROOT
1551 TDirectory *cdir = gDirectory;
1554 Error("StartAnalysis","Analysis manager was not initialized !");
1558 if (!CheckTasks()) Fatal("StartAnalysis", "Not all needed libraries were loaded");
1560 printf("StartAnalysis %s\n",GetName());
1561 AliLog::SetGlobalLogLevel(AliLog::kInfo);
1563 fMaxEntries = nentries;
1565 TString anaType = type;
1567 fMode = kLocalAnalysis;
1568 if (anaType.Contains("file")) fIsRemote = kTRUE;
1569 if (anaType.Contains("proof")) fMode = kProofAnalysis;
1570 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1571 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
1573 if (fMode == kGridAnalysis) {
1575 if (!anaType.Contains("terminate")) {
1576 if (!fGridHandler) {
1577 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1578 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1582 // Write analysis manager in the analysis file
1583 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1584 // run local task configuration
1586 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1587 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1592 // Terminate grid analysis
1593 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1594 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1595 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1596 if (!fGridHandler->MergeOutputs()) {
1597 // Return if outputs could not be merged or if it alien handler
1598 // was configured for offline mode or local testing.
1603 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1604 ImportWrappers(NULL);
1610 SetEventLoop(kFALSE);
1611 // Enable event loop mode if a tree was provided
1612 if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1615 TString ttype = "TTree";
1616 if (tree && tree->IsA() == TChain::Class()) {
1617 chain = (TChain*)tree;
1618 if (!chain || !chain->GetListOfFiles()->First()) {
1619 Error("StartAnalysis", "Cannot process null or empty chain...");
1626 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1627 if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1628 // Initialize locally all tasks (happens for all modes)
1630 AliAnalysisTask *task;
1634 case kLocalAnalysis:
1635 if (!tree && !fGridHandler) {
1636 TIter nextT(fTasks);
1637 // Call CreateOutputObjects for all tasks
1639 Bool_t dirStatus = TH1::AddDirectoryStatus();
1640 while ((task=(AliAnalysisTask*)nextT())) {
1641 TH1::AddDirectory(kFALSE);
1642 task->CreateOutputObjects();
1643 if (!task->CheckPostData()) {
1644 Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
1645 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
1646 ########### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ###########", task->GetName(), task->ClassName());
1648 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1652 TH1::AddDirectory(dirStatus);
1653 if (IsExternalLoop()) {
1654 Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1655 \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1662 fSelector = new AliAnalysisSelector(this);
1663 // Check if a plugin handler is used
1665 // Get the chain from the plugin
1666 TString dataType = "esdTree";
1667 if (fInputEventHandler) {
1668 dataType = fInputEventHandler->GetDataType();
1672 chain = fGridHandler->GetChainForTestMode(dataType);
1674 Error("StartAnalysis", "No chain for test mode. Aborting.");
1677 cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1678 retv = chain->Process(fSelector, "", nentries, firstentry);
1681 // Run tree-based analysis via AliAnalysisSelector
1682 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1683 retv = tree->Process(fSelector, "", nentries, firstentry);
1685 case kProofAnalysis:
1687 // Check if the plugin is used
1689 return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1691 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1692 Error("StartAnalysis", "No PROOF!!! Exiting.");
1696 line = Form("gProof->AddInput((TObject*)%p);", this);
1697 gROOT->ProcessLine(line);
1700 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1701 retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1703 Error("StartAnalysis", "No chain!!! Exiting.");
1710 if (!anaType.Contains("terminate")) {
1711 if (!fGridHandler) {
1712 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1713 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1717 // Write analysis manager in the analysis file
1718 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1719 // Start the analysis via the handler
1720 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1721 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1726 // Terminate grid analysis
1727 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1728 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1729 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1730 if (!fGridHandler->MergeOutputs()) {
1731 // Return if outputs could not be merged or if it alien handler
1732 // was configured for offline mode or local testing.
1737 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1738 ImportWrappers(NULL);
1742 case kMixingAnalysis:
1743 // Run event mixing analysis
1745 Error("StartAnalysis", "Cannot run event mixing without event pool");
1749 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1750 fSelector = new AliAnalysisSelector(this);
1751 while ((chain=fEventPool->GetNextChain())) {
1753 // Call NotifyBinChange for all tasks
1754 while ((task=(AliAnalysisTask*)next()))
1755 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1756 retv = chain->Process(fSelector);
1758 Error("StartAnalysis", "Mixing analysis failed");
1763 PackOutput(fSelector->GetOutputList());
1770 //______________________________________________________________________________
1771 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1773 // Start analysis for this manager on a given dataset. Analysis task can be:
1774 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1776 Error("StartAnalysis","Analysis manager was not initialized !");
1780 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1781 TString anaType = type;
1783 if (!anaType.Contains("proof")) {
1784 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1787 fMode = kProofAnalysis;
1789 SetEventLoop(kTRUE);
1790 // Set the dataset flag
1791 TObject::SetBit(kUseDataSet);
1794 // Start proof analysis using the grid handler
1795 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1796 Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
1799 // Check if the plugin is in test mode
1800 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
1801 dataset = "test_collection";
1803 dataset = fGridHandler->GetProofDataSet();
1807 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1808 Error("StartAnalysis", "No PROOF!!! Exiting.");
1812 // Initialize locally all tasks
1815 line = Form("gProof->AddInput((TObject*)%p);", this);
1816 gROOT->ProcessLine(line);
1818 line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1819 dataset, nentries, firstentry);
1820 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1821 retv = (Long_t)gROOT->ProcessLine(line);
1825 //______________________________________________________________________________
1826 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
1828 // Opens according the option the file specified by cont->GetFileName() and changes
1829 // current directory to cont->GetFolderName(). If the file was already opened, it
1830 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
1831 // be optionally ignored.
1832 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1833 TString filename = cont->GetFileName();
1835 if (filename.IsNull()) {
1836 ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
1839 if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
1841 f = mgr->OpenProofFile(cont,option);
1843 // Check first if the file is already opened
1844 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1846 // Check if option "UPDATE" was preserved
1847 TString opt(option);
1849 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1850 ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1852 f = TFile::Open(filename, option);
1855 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1859 // Check for a folder request
1860 TString dir = cont->GetFolderName();
1861 if (!dir.IsNull()) {
1862 if (!f->GetDirectory(dir)) f->mkdir(dir);
1867 ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
1868 cont->SetFile(NULL);
1872 //______________________________________________________________________________
1873 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
1875 // Opens a special output file used in PROOF.
1877 TString filename = cont->GetFileName();
1878 if (cont == fCommonOutput) {
1879 if (fOutputEventHandler) {
1880 if (strlen(extaod)) filename = extaod;
1881 filename = fOutputEventHandler->GetOutputFileName();
1883 else Fatal("OpenProofFile","No output container. Exiting.");
1886 if (fMode!=kProofAnalysis || !fSelector) {
1887 Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
1890 if (fSpecialOutputLocation.Length()) {
1891 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1893 // Check if option "UPDATE" was preserved
1894 TString opt(option);
1896 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1897 ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1899 f = new TFile(filename, option);
1901 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1905 // Check for a folder request
1906 TString dir = cont->GetFolderName();
1908 if (!f->GetDirectory(dir)) f->mkdir(dir);
1913 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1914 cont->SetFile(NULL);
1917 // Check if there is already a proof output file in the output list
1918 TObject *pof = fSelector->GetOutputList()->FindObject(filename);
1920 // Get the actual file
1921 line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
1922 filename = (const char*)gROOT->ProcessLine(line);
1924 printf("File: %s already booked via TProofOutputFile\n", filename.Data());
1926 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1928 Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
1931 // Check if option "UPDATE" was preserved
1932 TString opt(option);
1934 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1935 Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
1937 if (cont->IsRegisterDataset()) {
1938 TString dsetName = filename;
1939 dsetName.ReplaceAll(".root", cont->GetTitle());
1940 dsetName.ReplaceAll(":","_");
1941 if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
1942 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
1944 if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
1945 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
1947 if (fDebug > 1) printf("=== %s\n", line.Data());
1948 gROOT->ProcessLine(line);
1949 line = Form("pf->OpenFile(\"%s\");", option);
1950 gROOT->ProcessLine(line);
1953 gROOT->ProcessLine("pf->Print()");
1954 printf(" == proof file name: %s", f->GetName());
1956 // Add to proof output list
1957 line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
1958 if (fDebug > 1) printf("=== %s\n", line.Data());
1959 gROOT->ProcessLine(line);
1961 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1965 // Check for a folder request
1966 TString dir = cont->GetFolderName();
1967 if (!dir.IsNull()) {
1968 if (!f->GetDirectory(dir)) f->mkdir(dir);
1973 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1974 cont->SetFile(NULL);
1978 //______________________________________________________________________________
1979 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1981 // Execute analysis.
1982 static Long64_t nentries = 0;
1983 static TTree *lastTree = 0;
1984 static TStopwatch *timer = new TStopwatch();
1985 // Only the first call to Process will trigger a true Notify. Other Notify
1986 // coming before is ignored.
1987 if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) {
1988 TObject::SetBit(AliAnalysisManager::kTrueNotify);
1991 if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
1993 if (fTree && (fTree != lastTree)) {
1994 nentries += fTree->GetEntries();
1997 if (!fNcalls) timer->Start();
1998 if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, TMath::Min(fMaxEntries,nentries), timer, kFALSE);
2001 TDirectory *cdir = gDirectory;
2002 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
2003 if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
2005 Error("ExecAnalysis", "Analysis manager was not initialized !");
2010 AliAnalysisTask *task;
2011 // Check if the top tree is active.
2013 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2014 AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
2016 // De-activate all tasks
2017 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
2018 AliAnalysisDataContainer *cont = fCommonInput;
2019 if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
2021 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
2025 cont->SetData(fTree); // This will notify all consumers
2026 Long64_t entry = fTree->GetTree()->GetReadEntry();
2028 // Call BeginEvent() for optional input/output and MC services
2029 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
2030 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
2031 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
2033 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2034 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2036 // Execute the tasks
2037 // TIter next1(cont->GetConsumers());
2038 TIter next1(fTopTasks);
2040 while ((task=(AliAnalysisTask*)next1())) {
2042 cout << " Executing task " << task->GetName() << endl;
2044 task->ExecuteTask(option);
2046 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2047 AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
2051 // Call FinishEvent() for optional output and MC services
2052 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2053 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2054 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2055 // Gather system information if requested
2056 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2057 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
2061 // The event loop is not controlled by TSelector
2063 // Call BeginEvent() for optional input/output and MC services
2064 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
2065 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
2066 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
2068 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2069 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2070 TIter next2(fTopTasks);
2071 while ((task=(AliAnalysisTask*)next2())) {
2072 task->SetActive(kTRUE);
2074 cout << " Executing task " << task->GetName() << endl;
2076 task->ExecuteTask(option);
2080 // Call FinishEvent() for optional output and MC services
2081 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2082 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2083 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2084 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2085 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2089 //______________________________________________________________________________
2090 Bool_t AliAnalysisManager::IsPipe(std::ostream &out)
2092 // Check if the stdout is connected to a pipe (C.Holm)
2093 Bool_t ispipe = kFALSE;
2094 out.seekp(0, std::ios_base::cur);
2097 if (errno == ESPIPE) ispipe = kTRUE;
2102 //______________________________________________________________________________
2103 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2105 // Set the input event handler and create a container for it.
2106 fInputEventHandler = handler;
2107 if (!fCommonInput) fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2110 //______________________________________________________________________________
2111 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2113 // Set the input event handler and create a container for it.
2114 fOutputEventHandler = handler;
2115 if (!fCommonOutput) fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2116 fCommonOutput->SetSpecialOutput();
2119 //______________________________________________________________________________
2120 void AliAnalysisManager::SetDebugLevel(UInt_t level)
2122 // Set verbosity of the analysis manager. If the progress bar is used, the call is ignored
2123 if (TObject::TestBit(kUseProgressBar)) {
2124 Info("SetDebugLevel","Ignored. Disable the progress bar first.");
2130 //______________________________________________________________________________
2131 void AliAnalysisManager::SetUseProgressBar(Bool_t flag, Int_t freq)
2133 // Enable a text mode progress bar. Resets debug level to 0.
2134 Info("SetUseProgressBar", "Progress bar enabled, updated every %d events.\n ### NOTE: Debug level reset to 0 ###", freq);
2135 TObject::SetBit(kUseProgressBar,flag);
2136 fPBUpdateFreq = freq;
2140 //______________________________________________________________________________
2141 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2143 // This method is used externally to register output files which are not
2144 // connected to any output container, so that the manager can properly register,
2145 // retrieve or merge them when running in distributed mode. The file names are
2146 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2147 if (fExtraFiles.Contains(fname)) return;
2148 if (fExtraFiles.Length()) fExtraFiles += " ";
2149 fExtraFiles += fname;
2152 //______________________________________________________________________________
2153 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2155 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2159 TObject *pof = source->FindObject(filename);
2160 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2161 Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2164 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2165 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2166 TString clientUrl(chUrl);
2167 TString fullPath_str(fullPath);
2168 if (clientUrl.Contains("localhost")){
2169 TObjArray* array = fullPath_str.Tokenize ( "//" );
2170 TObjString *strobj = ( TObjString *)array->At(1);
2171 TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2172 TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2173 fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2174 fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2175 if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2179 else if (clientUrl.Contains("__lite__")) {
2180 // Special case for ProofLite environement - get file info and copy.
2181 gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2182 fullPath_str = Form("%s/%s", tmp, fullPath);
2185 Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2186 Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename);
2188 Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2192 //______________________________________________________________________________
2193 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2195 // Fill analysis type in the provided string.
2197 case kLocalAnalysis:
2200 case kProofAnalysis:
2206 case kMixingAnalysis:
2211 //______________________________________________________________________________
2212 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2214 // Validate all output files.
2215 TIter next(fOutputs);
2216 AliAnalysisDataContainer *output;
2217 TDirectory *cdir = gDirectory;
2218 TString openedFiles;
2219 while ((output=(AliAnalysisDataContainer*)next())) {
2220 if (output->IsRegisterDataset()) continue;
2221 TString filename = output->GetFileName();
2222 if (filename == "default") {
2223 if (!fOutputEventHandler) continue;
2224 filename = fOutputEventHandler->GetOutputFileName();
2225 // Main AOD may not be there
2226 if (gSystem->AccessPathName(filename)) continue;
2228 // Check if the file is closed
2229 if (openedFiles.Contains(filename)) continue;;
2230 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2232 Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2233 // Clear file list to release object ownership to user.
2237 file = TFile::Open(filename);
2238 if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2239 Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2244 openedFiles += filename;
2251 //______________________________________________________________________________
2252 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2254 // Implements a nice text mode progress bar.
2255 static Long64_t icount = 0;
2256 static TString oname;
2257 static TString nname;
2258 static Long64_t ocurrent = 0;
2259 static Long64_t osize = 0;
2260 static Int_t oseconds = 0;
2261 static TStopwatch *owatch = 0;
2262 static Bool_t oneoftwo = kFALSE;
2263 static Int_t nrefresh = 0;
2264 static Int_t nchecks = 0;
2265 static char lastChar = 0;
2266 const char symbol[4] = {'-','\\','|','/'};
2268 if (!lastChar) lastChar = (IsPipe(std::cerr))?'\r':'\n';
2274 ocurrent = TMath::Abs(current);
2275 osize = TMath::Abs(size);
2276 if (ocurrent > osize) ocurrent=osize;
2281 if ((current % fPBUpdateFreq) != 0) return;
2283 char progress[11] = " ";
2284 Int_t ichar = icount%4;
2289 if (owatch && !last) {
2291 time = owatch->RealTime();
2292 seconds = int(time) % 60;
2293 minutes = (int(time) / 60) % 60;
2294 hours = (int(time) / 60 / 60);
2296 if (oseconds==seconds) {
2300 oneoftwo = !oneoftwo;
2304 if (refresh && oneoftwo) {
2306 if (nchecks <= 0) nchecks = nrefresh+1;
2307 Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2308 oname = Form(" == %d%% ==", pctdone);
2310 Double_t percent = 100.0*ocurrent/osize;
2311 Int_t nchar = Int_t(percent/10);
2312 if (nchar>10) nchar=10;
2314 for (i=0; i<nchar; i++) progress[i] = '=';
2315 progress[nchar] = symbol[ichar];
2316 for (i=nchar+1; i<10; i++) progress[i] = ' ';
2317 progress[10] = '\0';
2320 if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2321 else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2322 else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2324 Int_t full = Int_t(ocurrent > 0 ?
2325 time * (float(osize)/ocurrent) + .5 :
2327 Int_t remain = Int_t(full - time);
2328 Int_t rsec = remain % 60;
2329 Int_t rmin = (remain / 60) % 60;
2330 Int_t rhour = (remain / 60 / 60);
2331 fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d ETA %.2d:%.2d:%.2d%c",
2332 percent, hours, minutes, seconds, rhour, rmin, rsec, lastChar);
2334 else fprintf(stderr, "[%6.2f %%]%c", percent, lastChar);
2335 if (refresh && oneoftwo) oname = nname;
2336 if (owatch) owatch->Continue();
2345 fprintf(stderr, "\n");
2349 //______________________________________________________________________________
2350 void AliAnalysisManager::DoLoadBranch(const char *name)
2352 // Get tree and load branch if needed.
2353 static Long64_t crtEntry = -100;
2355 if (fAutoBranchHandling || !fTree)
2358 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2360 br = fTree->GetBranch(name);
2362 Error("DoLoadBranch", "Could not find branch %s",name);
2367 if (br->GetReadEntry()==fCurrentEntry) return;
2368 Int_t ret = br->GetEntry(GetCurrentEntry());
2370 Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2371 if (crtEntry != fCurrentEntry) {
2372 CountEvent(1,0,1,0);
2373 crtEntry = fCurrentEntry;
2376 if (crtEntry != fCurrentEntry) {
2377 CountEvent(1,1,0,0);
2378 crtEntry = fCurrentEntry;
2383 //______________________________________________________________________________
2384 void AliAnalysisManager::AddStatisticsTask(UInt_t offlineMask)
2386 // Add the statistics task to the manager.
2388 Info("AddStatisticsTask", "Already added");
2391 TString line = Form("AliAnalysisTaskStat::AddToManager(%u);", offlineMask);
2392 gROOT->ProcessLine(line);
2395 //______________________________________________________________________________
2396 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2398 // Bookkeep current event;
2399 if (!fStatistics) return;
2400 fStatistics->AddInput(ninput);
2401 fStatistics->AddProcessed(nprocessed);
2402 fStatistics->AddFailed(nfailed);
2403 fStatistics->AddAccepted(naccepted);
2406 //______________________________________________________________________________
2407 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2409 // Add a line in the statistics message. If available, the statistics message is written
2410 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2412 if (!strlen(line)) return;
2413 if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2414 fStatisticsMsg += line;
2417 //______________________________________________________________________________
2418 void AliAnalysisManager::WriteStatisticsMsg(Int_t)
2420 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2421 static Bool_t done = kFALSE;
2424 if (!fStatistics) return;
2426 AddStatisticsMsg(Form("Number of input events: %lld",fStatistics->GetNinput()));
2427 AddStatisticsMsg(Form("Number of processed events: %lld",fStatistics->GetNprocessed()));
2428 AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2429 AddStatisticsMsg(Form("Number of accepted events for mask %s: %lld", AliAnalysisStatistics::GetMaskAsString(fStatistics->GetOfflineMask()), fStatistics->GetNaccepted()));
2430 out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2431 fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2432 fStatistics->GetNaccepted()), ios::out);
2433 out << fStatisticsMsg << endl;
2437 //______________________________________________________________________________
2438 const char* AliAnalysisManager::GetOADBPath()
2440 // returns the path of the OADB
2441 // this static function just depends on environment variables
2443 static TString oadbPath;
2445 if (gSystem->Getenv("OADB_PATH"))
2446 oadbPath = gSystem->Getenv("OADB_PATH");
2447 else if (gSystem->Getenv("ALICE_ROOT"))
2448 oadbPath.Form("%s/OADB", gSystem->Getenv("ALICE_ROOT"));
2450 ::Fatal("AliAnalysisManager::GetOADBPath", "Cannot figure out AODB path. Define ALICE_ROOT or OADB_PATH!");
2455 //______________________________________________________________________________
2456 void AliAnalysisManager::SetGlobalStr(const char *key, const char *value)
2458 // Define a custom string variable mapped to a global unique name. The variable
2459 // can be then retrieved by a given analysis macro via GetGlobalStr(key).
2460 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2462 ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2465 Bool_t valid = kFALSE;
2466 TString existing = AliAnalysisManager::GetGlobalStr(key, valid);
2468 ::Error("AliAnalysisManager::SetGlobalStr", "Global %s = %s already defined.", key, existing.Data());
2471 mgr->GetGlobals()->Add(new TObjString(key), new TObjString(value));
2474 //______________________________________________________________________________
2475 const char *AliAnalysisManager::GetGlobalStr(const char *key, Bool_t &valid)
2477 // Static method to retrieve a global variable defined via SetGlobalStr.
2479 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2481 TObject *value = mgr->GetGlobals()->GetValue(key);
2482 if (!value) return 0;
2484 return value->GetName();
2487 //______________________________________________________________________________
2488 void AliAnalysisManager::SetGlobalInt(const char *key, Int_t value)
2490 // Define a custom integer variable mapped to a global unique name. The variable
2491 // can be then retrieved by a given analysis macro via GetGlobalInt(key).
2492 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2494 ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2497 Bool_t valid = kFALSE;
2498 Int_t existing = AliAnalysisManager::GetGlobalInt(key, valid);
2500 ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %i already defined.", key, existing);
2503 mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%i",value)));
2506 //______________________________________________________________________________
2507 Int_t AliAnalysisManager::GetGlobalInt(const char *key, Bool_t &valid)
2509 // Static method to retrieve a global variable defined via SetGlobalInt.
2511 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2513 TObject *value = mgr->GetGlobals()->GetValue(key);
2514 if (!value) return 0;
2516 TString s = value->GetName();
2520 //______________________________________________________________________________
2521 void AliAnalysisManager::SetGlobalDbl(const char *key, Double_t value)
2523 // Define a custom double precision variable mapped to a global unique name. The variable
2524 // can be then retrieved by a given analysis macro via GetGlobalInt(key).
2525 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2527 ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2530 Bool_t valid = kFALSE;
2531 Double_t existing = AliAnalysisManager::GetGlobalDbl(key, valid);
2533 ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %g already defined.", key, existing);
2536 mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%f.16",value)));
2539 //______________________________________________________________________________
2540 Double_t AliAnalysisManager::GetGlobalDbl(const char *key, Bool_t &valid)
2542 // Static method to retrieve a global variable defined via SetGlobalDbl.
2544 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2546 TObject *value = mgr->GetGlobals()->GetValue(key);
2547 if (!value) return 0;
2549 TString s = value->GetName();