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)
240 Int_t index = s.Index("/alice/sim");
242 for (Int_t i=0; i<3; i++) {
243 index = s.Index("/", index+1);
244 if (index<0) return 0;
249 index = s.Index("/alice/data");
251 for (Int_t i=0; i<4; i++) {
252 index = s.Index("/", index+1);
253 if (index<0) return 0;
261 //______________________________________________________________________________
262 Bool_t AliAnalysisManager::Init(TTree *tree)
264 // The Init() function is called when the selector needs to initialize
265 // a new tree or chain. Typically here the branch addresses of the tree
266 // will be set. It is normaly not necessary to make changes to the
267 // generated code, but the routine can be extended by the user if needed.
268 // Init() will be called many times when running with PROOF.
269 Bool_t init = kFALSE;
270 if (!tree) return kFALSE; // Should not happen - protected in selector caller
272 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
274 // Call InitTree of EventHandler
275 if (fOutputEventHandler) {
276 if (fMode == kProofAnalysis) {
277 init = fOutputEventHandler->Init(0x0, "proof");
279 init = fOutputEventHandler->Init(0x0, "local");
282 Error("Init", "Output event handler failed to initialize");
287 if (fInputEventHandler) {
288 if (fMode == kProofAnalysis) {
289 init = fInputEventHandler->Init(tree, "proof");
291 init = fInputEventHandler->Init(tree, "local");
294 Error("Init", "Input event handler failed to initialize tree");
298 // If no input event handler we need to get the tree once
300 if(!tree->GetTree()) {
301 Long64_t readEntry = tree->LoadTree(0);
302 if (readEntry == -2) {
303 Error("Init", "Input tree has no entry. Exiting");
309 if (fMCtruthEventHandler) {
310 if (fMode == kProofAnalysis) {
311 init = fMCtruthEventHandler->Init(0x0, "proof");
313 init = fMCtruthEventHandler->Init(0x0, "local");
316 Error("Init", "MC event handler failed to initialize");
321 if (!fInitOK) InitAnalysis();
322 if (!fInitOK) return kFALSE;
325 AliAnalysisDataContainer *top = fCommonInput;
326 if (!top) top = (AliAnalysisDataContainer*)fInputs->At(0);
328 Error("Init","No top input container !");
332 CheckBranches(kFALSE);
334 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
339 //______________________________________________________________________________
340 void AliAnalysisManager::SlaveBegin(TTree *tree)
342 // The SlaveBegin() function is called after the Begin() function.
343 // When running with PROOF SlaveBegin() is called on each slave server.
344 // The tree argument is deprecated (on PROOF 0 is passed).
345 if (fDebug > 1) printf("->AliAnalysisManager::SlaveBegin()\n");
346 if (!CheckTasks()) Fatal("SlaveBegin", "Not all needed libraries were loaded");
347 static Bool_t isCalled = kFALSE;
348 Bool_t init = kFALSE;
349 Bool_t initOK = kTRUE;
351 TDirectory *curdir = gDirectory;
352 // Call SlaveBegin only once in case of mixing
353 if (isCalled && fMode==kMixingAnalysis) return;
355 // Call Init of EventHandler
356 if (fOutputEventHandler) {
357 if (fMode == kProofAnalysis) {
358 // Merging AOD's in PROOF via TProofOutputFile
359 if (fDebug > 1) printf(" Initializing AOD output file %s...\n", fOutputEventHandler->GetOutputFileName());
360 init = fOutputEventHandler->Init("proof");
361 if (!init) msg = "Failed to initialize output handler on worker";
363 init = fOutputEventHandler->Init("local");
364 if (!init) msg = "Failed to initialize output handler";
367 if (!fSelector) Error("SlaveBegin", "Selector not set");
368 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
371 if (fInputEventHandler) {
372 fInputEventHandler->SetInputTree(tree);
373 if (fMode == kProofAnalysis) {
374 init = fInputEventHandler->Init("proof");
375 if (!init) msg = "Failed to initialize input handler on worker";
377 init = fInputEventHandler->Init("local");
378 if (!init) msg = "Failed to initialize input handler";
381 if (!fSelector) Error("SlaveBegin", "Selector not set");
382 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
385 if (fMCtruthEventHandler) {
386 if (fMode == kProofAnalysis) {
387 init = fMCtruthEventHandler->Init("proof");
388 if (!init) msg = "Failed to initialize MC handler on worker";
390 init = fMCtruthEventHandler->Init("local");
391 if (!init) msg = "Failed to initialize MC handler";
394 if (!fSelector) Error("SlaveBegin", "Selector not set");
395 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
397 if (curdir) curdir->cd();
401 AliAnalysisTask *task;
402 // Call CreateOutputObjects for all tasks
403 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
404 Bool_t dirStatus = TH1::AddDirectoryStatus();
406 while ((task=(AliAnalysisTask*)next())) {
408 // Start with memory as current dir and make sure by default histograms do not get attached to files.
409 TH1::AddDirectory(kFALSE);
410 task->CreateOutputObjects();
411 if (!task->CheckPostData()) {
412 Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
413 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
414 ####### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ##########", task->GetName(), task->ClassName());
416 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
419 TH1::AddDirectory(dirStatus);
420 if (curdir) curdir->cd();
421 if (fDebug > 1) printf("<-AliAnalysisManager::SlaveBegin()\n");
424 //______________________________________________________________________________
425 Bool_t AliAnalysisManager::Notify()
427 // The Notify() function is called when a new file is opened. This
428 // can be either for a new TTree in a TChain or when when a new TTree
429 // is started when using PROOF. It is normaly not necessary to make changes
430 // to the generated code, but the routine can be extended by the
431 // user if needed. The return value is currently not used.
432 if (!fTree) return kFALSE;
433 if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) return kFALSE;
435 fTable.Clear("nodelete"); // clearing the hash table may not be needed -> C.L.
436 if (fMode == kProofAnalysis) fIsRemote = kTRUE;
438 TFile *curfile = fTree->GetCurrentFile();
440 Error("Notify","No current file");
444 if (fDebug > 1) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
445 Int_t run = AliAnalysisManager::GetRunFromAlienPath(curfile->GetName());
446 if (run && (run != fRunFromPath)) {
448 if (fDebug > 1) printf(" ### run found from path: %d\n", run);
451 AliAnalysisTask *task;
453 // Call Notify of the event handlers
454 if (fInputEventHandler) {
455 fInputEventHandler->Notify(curfile->GetName());
458 if (fOutputEventHandler) {
459 fOutputEventHandler->Notify(curfile->GetName());
462 if (fMCtruthEventHandler) {
463 fMCtruthEventHandler->Notify(curfile->GetName());
466 // Call Notify for all tasks
467 while ((task=(AliAnalysisTask*)next()))
470 if (fDebug > 1) printf("<-AliAnalysisManager::Notify()\n");
474 //______________________________________________________________________________
475 Bool_t AliAnalysisManager::Process(Long64_t)
477 // The Process() function is called for each entry in the tree (or possibly
478 // keyed object in the case of PROOF) to be processed. The entry argument
479 // specifies which entry in the currently loaded tree is to be processed.
480 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
481 // to read either all or the required parts of the data. When processing
482 // keyed objects with PROOF, the object is already loaded and is available
483 // via the fObject pointer.
485 // This function should contain the "body" of the analysis. It can contain
486 // simple or elaborate selection criteria, run algorithms on the data
487 // of the event and typically fill histograms.
489 // WARNING when a selector is used with a TChain, you must use
490 // the pointer to the current TTree to call GetEntry(entry).
491 // The entry is always the local entry number in the current tree.
492 // Assuming that fChain is the pointer to the TChain being processed,
493 // use fChain->GetTree()->GetEntry(entry).
495 // This method is obsolete. ExecAnalysis is called instead.
499 //______________________________________________________________________________
500 void AliAnalysisManager::PackOutput(TList *target)
502 // Pack all output data containers in the output list. Called at SlaveTerminate
503 // stage in PROOF case for each slave.
504 if (fDebug > 1) printf("->AliAnalysisManager::PackOutput()\n");
506 Error("PackOutput", "No target. Exiting.");
509 TDirectory *cdir = gDirectory;
511 if (fInputEventHandler) fInputEventHandler ->Terminate();
512 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
513 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
516 // Call FinishTaskOutput() for each event loop task (not called for
517 // post-event loop tasks - use Terminate() fo those)
518 TIter nexttask(fTasks);
519 AliAnalysisTask *task;
520 while ((task=(AliAnalysisTask*)nexttask())) {
521 if (!task->IsPostEventLoop()) {
522 if (fDebug > 1) printf("->FinishTaskOutput: task %s\n", task->GetName());
523 task->FinishTaskOutput();
525 if (fDebug > 1) printf("<-FinishTaskOutput: task %s\n", task->GetName());
528 // Write statistics message on the workers.
529 if (fStatistics) WriteStatisticsMsg(fNcalls);
531 if (fMode == kProofAnalysis) {
532 TIter next(fOutputs);
533 AliAnalysisDataContainer *output;
534 Bool_t isManagedByHandler = kFALSE;
537 while ((output=(AliAnalysisDataContainer*)next())) {
538 // Do not consider outputs of post event loop tasks
539 isManagedByHandler = kFALSE;
540 if (output->GetProducer() && output->GetProducer()->IsPostEventLoop()) continue;
541 const char *filename = output->GetFileName();
542 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
543 isManagedByHandler = kTRUE;
544 printf("#### Handler output. Extra: %s\n", fExtraFiles.Data());
545 filename = fOutputEventHandler->GetOutputFileName();
547 // Check if data was posted to this container. If not, issue an error.
548 if (!output->GetData() && !isManagedByHandler) {
549 Error("PackOutput", "No data for output container %s. Forgot to PostData ?", output->GetName());
552 if (!output->IsSpecialOutput()) {
554 if (strlen(filename) && !isManagedByHandler) {
555 // Backup current folder
556 TDirectory *opwd = gDirectory;
557 // File resident outputs.
558 // Check first if the file exists.
559 TString openoption = "RECREATE";
560 Bool_t firsttime = kTRUE;
561 if (filestmp.FindObject(output->GetFileName())) {
564 filestmp.Add(new TNamed(output->GetFileName(),""));
566 if (!gSystem->AccessPathName(output->GetFileName()) && !firsttime) openoption = "UPDATE";
567 // TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
568 // Save data to file, then close.
569 if (output->GetData()->InheritsFrom(TCollection::Class())) {
570 // If data is a collection, we set the name of the collection
571 // as the one of the container and we save as a single key.
572 TCollection *coll = (TCollection*)output->GetData();
573 coll->SetName(output->GetName());
574 // coll->Write(output->GetName(), TObject::kSingleKey);
576 if (output->GetData()->InheritsFrom(TTree::Class())) {
577 TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
578 // Save data to file, then close.
579 TTree *tree = (TTree*)output->GetData();
580 // Check if tree is in memory
581 if (tree->GetDirectory()==gROOT) tree->SetDirectory(gDirectory);
585 // output->GetData()->Write();
588 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
590 // printf(" file %s listing content:\n", filename);
593 // Clear file list to release object ownership to user.
596 output->SetFile(NULL);
597 // Restore current directory
598 if (opwd) opwd->cd();
600 // Memory-resident outputs
601 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
603 AliAnalysisDataWrapper *wrap = 0;
604 if (isManagedByHandler) {
605 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
606 wrap->SetName(output->GetName());
608 else wrap =output->ExportData();
609 // Output wrappers must NOT delete data after merging - the user owns them
610 wrap->SetDeleteData(kFALSE);
613 // Special outputs. The file must be opened and connected to the container.
614 TDirectory *opwd = gDirectory;
615 TFile *file = output->GetFile();
617 AliAnalysisTask *producer = output->GetProducer();
619 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
620 output->GetFileName(), output->GetName(), producer->ClassName());
623 TString outFilename = file->GetName();
624 if (fDebug > 1) printf("PackOutput %s: special output\n", output->GetName());
625 if (isManagedByHandler) {
626 // Terminate IO for files managed by the output handler
627 // file->Write() moved to AOD handler (A.G. 11.01.10)
628 // if (file) file->Write();
629 if (file && fDebug > 2) {
630 printf(" handled file %s listing content:\n", file->GetName());
633 fOutputEventHandler->TerminateIO();
636 // Release object ownership to users after writing data to file
637 if (output->GetData()->InheritsFrom(TCollection::Class())) {
638 // If data is a collection, we set the name of the collection
639 // as the one of the container and we save as a single key.
640 TCollection *coll = (TCollection*)output->GetData();
641 coll->SetName(output->GetName());
642 coll->Write(output->GetName(), TObject::kSingleKey);
644 if (output->GetData()->InheritsFrom(TTree::Class())) {
645 TTree *tree = (TTree*)output->GetData();
646 tree->SetDirectory(file);
649 output->GetData()->Write();
653 printf(" file %s listing content:\n", output->GetFileName());
656 // Clear file list to release object ownership to user.
659 output->SetFile(NULL);
661 // Restore current directory
662 if (opwd) opwd->cd();
663 // Check if a special output location was provided or the output files have to be merged
664 if (strlen(fSpecialOutputLocation.Data())) {
665 TString remote = fSpecialOutputLocation;
667 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
668 if (remote.BeginsWith("alien:")) {
669 gROOT->ProcessLine("TGrid::Connect(\"alien:\", gProofServ->GetUser());");
670 remote += outFilename;
671 remote.ReplaceAll(".root", Form("_%d.root", gid));
673 remote += Form("%s_%d_", gSystem->HostName(), gid);
674 remote += outFilename;
677 Info("PackOutput", "Output file for container %s to be copied \n at: %s. No merging.",
678 output->GetName(), remote.Data());
679 TFile::Cp ( outFilename.Data(), remote.Data() );
680 // Copy extra outputs
681 if (fExtraFiles.Length() && isManagedByHandler) {
682 TObjArray *arr = fExtraFiles.Tokenize(" ");
684 TIter nextfilename(arr);
685 while ((os=(TObjString*)nextfilename())) {
686 outFilename = os->GetString();
687 remote = fSpecialOutputLocation;
689 if (remote.BeginsWith("alien://")) {
690 remote += outFilename;
691 remote.ReplaceAll(".root", Form("_%d.root", gid));
693 remote += Form("%s_%d_", gSystem->HostName(), gid);
694 remote += outFilename;
697 Info("PackOutput", "Extra AOD file %s to be copied \n at: %s. No merging.",
698 outFilename.Data(), remote.Data());
699 TFile::Cp ( outFilename.Data(), remote.Data() );
704 // No special location specified-> use TProofOutputFile as merging utility
705 // The file at this output slot must be opened in CreateOutputObjects
706 if (fDebug > 1) printf(" File for container %s to be merged via file merger...\n", output->GetName());
712 if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
715 //______________________________________________________________________________
716 void AliAnalysisManager::ImportWrappers(TList *source)
718 // Import data in output containers from wrappers coming in source.
719 if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
720 TIter next(fOutputs);
721 AliAnalysisDataContainer *cont;
722 AliAnalysisDataWrapper *wrap;
724 Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
725 TDirectory *cdir = gDirectory;
726 while ((cont=(AliAnalysisDataContainer*)next())) {
728 if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
729 if (cont->IsRegisterDataset()) continue;
730 const char *filename = cont->GetFileName();
731 Bool_t isManagedByHandler = kFALSE;
732 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
733 isManagedByHandler = kTRUE;
734 filename = fOutputEventHandler->GetOutputFileName();
736 if (cont->IsSpecialOutput() || inGrid) {
737 if (strlen(fSpecialOutputLocation.Data())) continue;
738 // Copy merged file from PROOF scratch space.
739 // In case of grid the files are already in the current directory.
741 if (isManagedByHandler && fExtraFiles.Length()) {
742 // Copy extra registered dAOD files.
743 TObjArray *arr = fExtraFiles.Tokenize(" ");
745 TIter nextfilename(arr);
746 while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
749 if (!GetFileFromWrapper(filename, source)) continue;
751 // Normally we should connect data from the copied file to the
752 // corresponding output container, but it is not obvious how to do this
753 // automatically if several objects in file...
754 TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
755 if (!f) f = TFile::Open(filename, "READ");
757 Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
762 // Cd to the directory pointed by the container
763 TString folder = cont->GetFolderName();
764 if (!folder.IsNull()) f->cd(folder);
765 // Try to fetch first an object having the container name.
766 obj = gDirectory->Get(cont->GetName());
768 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",
769 cont->GetType()->GetName(), cont->GetName(), filename, cont->GetFolderName(), cont->GetName());
772 wrap = new AliAnalysisDataWrapper(obj);
773 wrap->SetDeleteData(kFALSE);
775 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
777 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
782 printf(" Importing data for container %s\n", cont->GetName());
783 if (strlen(filename)) printf(" -> file %s\n", filename);
786 cont->ImportData(wrap);
788 if (cdir) cdir->cd();
789 if (fDebug > 1) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
792 //______________________________________________________________________________
793 void AliAnalysisManager::UnpackOutput(TList *source)
795 // Called by AliAnalysisSelector::Terminate only on the client.
796 if (fDebug > 1) printf("->AliAnalysisManager::UnpackOutput()\n");
798 Error("UnpackOutput", "No target. Exiting.");
801 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
803 if (fMode == kProofAnalysis) ImportWrappers(source);
805 TIter next(fOutputs);
806 AliAnalysisDataContainer *output;
807 while ((output=(AliAnalysisDataContainer*)next())) {
808 if (!output->GetData()) continue;
809 // Check if there are client tasks that run post event loop
810 if (output->HasConsumers()) {
811 // Disable event loop semaphore
812 output->SetPostEventLoop(kTRUE);
813 TObjArray *list = output->GetConsumers();
814 Int_t ncons = list->GetEntriesFast();
815 for (Int_t i=0; i<ncons; i++) {
816 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
817 task->CheckNotify(kTRUE);
818 // If task is active, execute it
819 if (task->IsPostEventLoop() && task->IsActive()) {
820 if (fDebug > 1) printf("== Executing post event loop task %s\n", task->GetName());
826 if (fDebug > 1) printf("<-AliAnalysisManager::UnpackOutput()\n");
829 //______________________________________________________________________________
830 void AliAnalysisManager::Terminate()
832 // The Terminate() function is the last function to be called during
833 // a query. It always runs on the client, it can be used to present
834 // the results graphically.
835 if (fDebug > 1) printf("->AliAnalysisManager::Terminate()\n");
836 TDirectory *cdir = gDirectory;
838 AliAnalysisTask *task;
839 AliAnalysisDataContainer *output;
842 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
843 // Call Terminate() for tasks
845 while (!IsSkipTerminate() && (task=(AliAnalysisTask*)next())) {
846 // Save all the canvases produced by the Terminate
847 TString pictname = Form("%s_%s", task->GetName(), task->ClassName());
851 AliSysInfo::AddStamp(Form("%s_TERMINATE",task->ClassName()),0, itask, 2);
853 if (TObject::TestBit(kSaveCanvases)) {
854 if (!gROOT->IsBatch()) {
855 if (fDebug>1) printf("Waiting 5 sec for %s::Terminate() to finish drawing ...\n", task->ClassName());
857 while (timer.CpuTime()<5) {
859 gSystem->ProcessEvents();
862 Int_t iend = gROOT->GetListOfCanvases()->GetEntries();
863 if (iend==0) continue;
865 for (Int_t ipict=0; ipict<iend; ipict++) {
866 canvas = (TCanvas*)gROOT->GetListOfCanvases()->At(ipict);
867 if (!canvas) continue;
868 canvas->SaveAs(Form("%s_%02d.gif", pictname.Data(),ipict));
870 gROOT->GetListOfCanvases()->Delete();
874 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
875 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
876 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
878 TObjArray *allOutputs = new TObjArray();
880 for (icont=0; icont<fOutputs->GetEntriesFast(); icont++) allOutputs->Add(fOutputs->At(icont));
881 if (!IsSkipTerminate())
882 for (icont=0; icont<fParamCont->GetEntriesFast(); icont++) allOutputs->Add(fParamCont->At(icont));
883 TIter next1(allOutputs);
884 TString handlerFile = "";
885 TString extraOutputs = "";
886 if (fOutputEventHandler) {
887 handlerFile = fOutputEventHandler->GetOutputFileName();
888 extraOutputs = fOutputEventHandler->GetExtraOutputs();
892 while ((output=(AliAnalysisDataContainer*)next1())) {
893 // Special outputs or grid files have the files already closed and written.
895 if (fMode == kGridAnalysis && icont<=fOutputs->GetEntriesFast()) continue;
896 if (fMode == kProofAnalysis) {
897 if (output->IsSpecialOutput() || output->IsRegisterDataset()) continue;
899 const char *filename = output->GetFileName();
900 TString openoption = "RECREATE";
901 if (!(strcmp(filename, "default"))) continue;
902 if (!strlen(filename)) continue;
903 if (!output->GetData()) continue;
904 TDirectory *opwd = gDirectory;
905 TFile *file = output->GetFile();
906 if (!file) file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
908 //if (handlerFile == filename && !gSystem->AccessPathName(filename)) openoption = "UPDATE";
909 Bool_t firsttime = kTRUE;
910 if (filestmp.FindObject(filename) || extraOutputs.Contains(filename)) {
913 filestmp.Add(new TNamed(filename,""));
915 if (!gSystem->AccessPathName(filename) && !firsttime) openoption = "UPDATE";
916 if (fDebug>1) printf("Opening file: %s option=%s\n",filename, openoption.Data());
917 file = new TFile(filename, openoption);
919 if (fDebug>1) printf("File <%s> already opened with option: <%s> \n", filename, file->GetOption());
920 openoption = file->GetOption();
921 if (openoption == "READ") {
922 if (fDebug>1) printf("...reopening in UPDATE mode\n");
923 file->ReOpen("UPDATE");
926 if (file->IsZombie()) {
927 Error("Terminate", "Cannot open output file %s", filename);
930 output->SetFile(file);
932 // Check for a folder request
933 TString dir = output->GetFolderName();
935 if (!file->GetDirectory(dir)) file->mkdir(dir);
938 if (fDebug > 1) printf("...writing container %s to file %s:%s\n", output->GetName(), file->GetName(), output->GetFolderName());
939 if (output->GetData()->InheritsFrom(TCollection::Class())) {
940 // If data is a collection, we set the name of the collection
941 // as the one of the container and we save as a single key.
942 TCollection *coll = (TCollection*)output->GetData();
943 coll->SetName(output->GetName());
944 coll->Write(output->GetName(), TObject::kSingleKey);
946 if (output->GetData()->InheritsFrom(TTree::Class())) {
947 TTree *tree = (TTree*)output->GetData();
948 tree->SetDirectory(gDirectory);
951 output->GetData()->Write();
954 if (opwd) opwd->cd();
958 while ((output=(AliAnalysisDataContainer*)next1())) {
959 // Close all files at output
960 TDirectory *opwd = gDirectory;
961 if (output->GetFile()) {
962 // Clear file list to release object ownership to user.
963 // output->GetFile()->Clear();
964 output->GetFile()->Close();
965 output->SetFile(NULL);
966 // Copy merged outputs in alien if requested
967 if (fSpecialOutputLocation.Length() &&
968 fSpecialOutputLocation.BeginsWith("alien://")) {
969 Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data());
970 TFile::Cp(output->GetFile()->GetName(),
971 Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
974 if (opwd) opwd->cd();
977 //Write statistics information on the client
978 if (fStatistics) WriteStatisticsMsg(fNcalls);
980 TDirectory *crtdir = gDirectory;
981 TFile f("syswatch.root", "RECREATE");
985 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
986 tree->SetName("syswatch");
987 tree->SetMarkerStyle(kCircle);
988 tree->SetMarkerColor(kBlue);
989 tree->SetMarkerSize(0.5);
990 if (!gROOT->IsBatch()) {
991 tree->SetAlias("event", "id0");
992 tree->SetAlias("task", "id1");
993 tree->SetAlias("stage", "id2");
994 // Already defined aliases
995 // tree->SetAlias("deltaT","stampSec-stampOldSec");
996 // tree->SetAlias("T","stampSec-first");
997 // tree->SetAlias("deltaVM","(pI.fMemVirtual-pIOld.fMemVirtual)");
998 // tree->SetAlias("VM","pI.fMemVirtual");
999 TCanvas *canvas = new TCanvas("SysInfo","SysInfo",10,10,1200,800);
1000 Int_t npads = 1 /*COO plot for all tasks*/ +
1001 fTopTasks->GetEntries() /*Exec plot per task*/ +
1002 1 /*Terminate plot for all tasks*/ +
1005 Int_t iopt = (Int_t)TMath::Sqrt((Double_t)npads);
1006 if (npads<iopt*(iopt+1))
1007 canvas->Divide(iopt, iopt+1, 0.01, 0.01);
1009 canvas->Divide(iopt+1, iopt+1, 0.01, 0.01);
1011 // draw the plot of deltaVM for Exec for each task
1012 for (itask=0; itask<fTopTasks->GetEntriesFast(); itask++) {
1013 task = (AliAnalysisTask*)fTopTasks->At(itask);
1015 cut = Form("task==%d && stage==1", itask);
1016 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1017 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1019 hist->SetTitle(Form("%s: Exec dVM[MB]/event", task->GetName()));
1020 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1023 // Draw the plot of deltaVM for CreateOutputObjects for all tasks
1025 tree->SetMarkerStyle(kFullTriangleUp);
1026 tree->SetMarkerColor(kRed);
1027 tree->SetMarkerSize(0.8);
1028 cut = "task>=0 && task<1000 && stage==0";
1029 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1030 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1032 hist->SetTitle("Memory in CreateOutputObjects()");
1033 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1034 hist->GetXaxis()->SetTitle("task");
1036 // draw the plot of deltaVM for Terminate for all tasks
1038 tree->SetMarkerStyle(kOpenSquare);
1039 tree->SetMarkerColor(kMagenta);
1040 cut = "task>=0 && task<1000 && stage==2";
1041 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1042 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1044 hist->SetTitle("Memory in Terminate()");
1045 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1046 hist->GetXaxis()->SetTitle("task");
1050 tree->SetMarkerStyle(kFullCircle);
1051 tree->SetMarkerColor(kGreen);
1052 cut = Form("task==%d && stage==1",fTopTasks->GetEntriesFast()-1);
1053 tree->Draw("VM:event",cut,"", 1234567890, 0);
1054 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1056 hist->SetTitle("Virtual memory");
1057 hist->GetYaxis()->SetTitle("VM [MB]");
1061 tree->SetMarkerStyle(kCircle);
1062 tree->SetMarkerColor(kBlue);
1063 tree->SetMarkerSize(0.5);
1068 if (crtdir) crtdir->cd();
1070 // Validate the output files
1071 if (ValidateOutputFiles() && fIsRemote && fMode!=kProofAnalysis) {
1073 out.open("outputs_valid", ios::out);
1077 if (fDebug > 1) printf("<-AliAnalysisManager::Terminate()\n");
1079 //______________________________________________________________________________
1080 void AliAnalysisManager::ProfileTask(Int_t itop, const char *option) const
1082 // Profiles the task having the itop index in the list of top (first level) tasks.
1083 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->At(itop);
1085 Error("ProfileTask", "There are only %d top tasks in the manager", fTopTasks->GetEntries());
1088 ProfileTask(task->GetName(), option);
1091 //______________________________________________________________________________
1092 void AliAnalysisManager::ProfileTask(const char *name, const char */*option*/) const
1094 // Profile a managed task after the execution of the analysis in case NSysInfo
1096 if (gSystem->AccessPathName("syswatch.root")) {
1097 Error("ProfileTask", "No file syswatch.root found in the current directory");
1100 if (gROOT->IsBatch()) return;
1101 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->FindObject(name);
1103 Error("ProfileTask", "No top task named %s known by the manager.", name);
1106 Int_t itop = fTopTasks->IndexOf(task);
1107 Int_t itask = fTasks->IndexOf(task);
1108 // Create canvas with 2 pads: first draw COO + Terminate, second Exec
1109 TDirectory *cdir = gDirectory;
1110 TFile f("syswatch.root");
1111 TTree *tree = (TTree*)f.Get("syswatch");
1113 Error("ProfileTask", "No tree named <syswatch> found in file syswatch.root");
1116 if (fDebug > 1) printf("=== Profiling task %s (class %s)\n", name, task->ClassName());
1117 TCanvas *canvas = new TCanvas(Form("profile_%d",itop),Form("Profile of task %s (class %s)",name,task->ClassName()),10,10,800,600);
1118 canvas->Divide(2, 2, 0.01, 0.01);
1122 // VM profile for COO and Terminate methods
1124 cut = Form("task==%d && (stage==0 || stage==2)",itask);
1125 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1126 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1128 hist->SetTitle("Alocated VM[MB] for COO and Terminate");
1129 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1130 hist->GetXaxis()->SetTitle("method");
1132 // CPU profile per event
1134 cut = Form("task==%d && stage==1",itop);
1135 tree->Draw("deltaT:event",cut,"", 1234567890, 0);
1136 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1138 hist->SetTitle("Execution time per event");
1139 hist->GetYaxis()->SetTitle("CPU/event [s]");
1141 // VM profile for Exec
1143 cut = Form("task==%d && stage==1",itop);
1144 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1145 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1147 hist->SetTitle("Alocated VM[MB] per event");
1148 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1153 if (cdir) cdir->cd();
1156 //______________________________________________________________________________
1157 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1159 // Adds a user task to the global list of tasks.
1161 Error("AddTask", "Cannot add task %s since InitAnalysis was already called", task->GetName());
1165 if (fTasks->FindObject(task)) {
1166 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1169 task->SetActive(kFALSE);
1173 //______________________________________________________________________________
1174 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1176 // Retreive task by name.
1177 if (!fTasks) return NULL;
1178 return (AliAnalysisTask*)fTasks->FindObject(name);
1181 //______________________________________________________________________________
1182 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
1183 TClass *datatype, EAliAnalysisContType type, const char *filename)
1185 // Create a data container of a certain type. Types can be:
1186 // kExchangeContainer = 0, used to exchange data between tasks
1187 // kInputContainer = 1, used to store input data
1188 // kOutputContainer = 2, used for writing result to a file
1189 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1190 // the output object to a folder inside the output file
1191 if (fContainers->FindObject(name)) {
1192 Error("CreateContainer","A container named %s already defined !",name);
1195 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1196 fContainers->Add(cont);
1198 case kInputContainer:
1201 case kOutputContainer:
1202 fOutputs->Add(cont);
1203 if (filename && strlen(filename)) {
1204 cont->SetFileName(filename);
1205 cont->SetDataOwned(kFALSE); // data owned by the file
1208 case kParamContainer:
1209 fParamCont->Add(cont);
1210 if (filename && strlen(filename)) {
1211 cont->SetFileName(filename);
1212 cont->SetDataOwned(kFALSE); // data owned by the file
1215 case kExchangeContainer:
1221 //______________________________________________________________________________
1222 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1223 AliAnalysisDataContainer *cont)
1225 // Connect input of an existing task to a data container.
1227 Error("ConnectInput", "Task pointer is NULL");
1230 if (!fTasks->FindObject(task)) {
1232 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1234 Bool_t connected = task->ConnectInput(islot, cont);
1238 //______________________________________________________________________________
1239 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1240 AliAnalysisDataContainer *cont)
1242 // Connect output of an existing task to a data container.
1244 Error("ConnectOutput", "Task pointer is NULL");
1247 if (!fTasks->FindObject(task)) {
1249 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1251 Bool_t connected = task->ConnectOutput(islot, cont);
1255 //______________________________________________________________________________
1256 void AliAnalysisManager::CleanContainers()
1258 // Clean data from all containers that have already finished all client tasks.
1259 TIter next(fContainers);
1260 AliAnalysisDataContainer *cont;
1261 while ((cont=(AliAnalysisDataContainer *)next())) {
1262 if (cont->IsOwnedData() &&
1263 cont->IsDataReady() &&
1264 cont->ClientsExecuted()) cont->DeleteData();
1268 //______________________________________________________________________________
1269 Bool_t AliAnalysisManager::InitAnalysis()
1271 // Initialization of analysis chain of tasks. Should be called after all tasks
1272 // and data containers are properly connected
1273 // Reset flag and remove valid_outputs file if exists
1274 if (fInitOK) return kTRUE;
1275 if (!gSystem->AccessPathName("outputs_valid"))
1276 gSystem->Unlink("outputs_valid");
1277 // Check for top tasks (depending only on input data containers)
1278 if (!fTasks->First()) {
1279 Error("InitAnalysis", "Analysis has no tasks !");
1283 AliAnalysisTask *task;
1284 AliAnalysisDataContainer *cont;
1287 Bool_t iszombie = kFALSE;
1288 Bool_t istop = kTRUE;
1290 while ((task=(AliAnalysisTask*)next())) {
1293 Int_t ninputs = task->GetNinputs();
1294 for (i=0; i<ninputs; i++) {
1295 cont = task->GetInputSlot(i)->GetContainer();
1299 fZombies->Add(task);
1303 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
1304 i, task->GetName());
1306 if (iszombie) continue;
1307 // Check if cont is an input container
1308 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1309 // Connect to parent task
1313 fTopTasks->Add(task);
1317 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1320 // Check now if there are orphan tasks
1321 for (i=0; i<ntop; i++) {
1322 task = (AliAnalysisTask*)fTopTasks->At(i);
1327 while ((task=(AliAnalysisTask*)next())) {
1328 if (!task->IsUsed()) {
1330 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1333 // Check the task hierarchy (no parent task should depend on data provided
1334 // by a daughter task)
1335 for (i=0; i<ntop; i++) {
1336 task = (AliAnalysisTask*)fTopTasks->At(i);
1337 if (task->CheckCircularDeps()) {
1338 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1343 // Check that all containers feeding post-event loop tasks are in the outputs list
1344 TIter nextcont(fContainers); // loop over all containers
1345 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1346 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1347 if (cont->HasConsumers()) {
1348 // Check if one of the consumers is post event loop
1349 TIter nextconsumer(cont->GetConsumers());
1350 while ((task=(AliAnalysisTask*)nextconsumer())) {
1351 if (task->IsPostEventLoop()) {
1352 fOutputs->Add(cont);
1359 // Check if all special output containers have a file name provided
1360 TIter nextout(fOutputs);
1361 while ((cont=(AliAnalysisDataContainer*)nextout())) {
1362 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1363 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1367 // Initialize requested branch list if needed
1368 if (!fAutoBranchHandling) {
1370 while ((task=(AliAnalysisTask*)next())) {
1371 if (!task->HasBranches()) {
1372 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\"",
1373 task->GetName(), task->ClassName());
1376 if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1377 Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1380 TString taskbranches;
1381 task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1382 if (taskbranches.IsNull()) {
1383 Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1384 task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1387 AddBranches(taskbranches);
1394 //______________________________________________________________________________
1395 void AliAnalysisManager::AddBranches(const char *branches)
1397 // Add branches to the existing fRequestedBranches.
1398 TString br(branches);
1399 TObjArray *arr = br.Tokenize(",");
1402 while ((obj=next())) {
1403 if (!fRequestedBranches.Contains(obj->GetName())) {
1404 if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1405 fRequestedBranches += obj->GetName();
1411 //______________________________________________________________________________
1412 void AliAnalysisManager::CheckBranches(Bool_t load)
1414 // The method checks the input branches to be loaded during the analysis.
1415 if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;
1416 TObjArray *arr = fRequestedBranches.Tokenize(",");
1419 while ((obj=next())) {
1420 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1422 br = fTree->GetBranch(obj->GetName());
1424 Error("CheckBranches", "Could not find branch %s",obj->GetName());
1429 if (load && br->GetReadEntry()!=GetCurrentEntry()) br->GetEntry(GetCurrentEntry());
1434 //______________________________________________________________________________
1435 Bool_t AliAnalysisManager::CheckTasks() const
1437 // Check consistency of tasks.
1438 Int_t ntasks = fTasks->GetEntries();
1440 Error("CheckTasks", "No tasks connected to the manager. This may be due to forgetting to compile the task or to load their library.");
1443 // Get the pointer to AliAnalysisTaskSE::Class()
1444 TClass *badptr = (TClass*)gROOT->ProcessLine("AliAnalysisTaskSE::Class()");
1445 // Loop all tasks to check if their corresponding library was loaded
1448 while ((obj=next())) {
1449 if (obj->IsA() == badptr) {
1450 Error("CheckTasks", "##################\n \
1451 Class for task %s NOT loaded. You probably forgot to load the library for this task (or compile it dynamically).\n###########################\n",obj->GetName());
1458 //______________________________________________________________________________
1459 void AliAnalysisManager::PrintStatus(Option_t *option) const
1461 // Print task hierarchy.
1463 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1466 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1468 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1469 TIter next(fTopTasks);
1470 AliAnalysisTask *task;
1471 while ((task=(AliAnalysisTask*)next()))
1472 task->PrintTask(option);
1474 if (!fAutoBranchHandling && !fRequestedBranches.IsNull())
1475 printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1477 TString sopt(option);
1480 if (sopt.Contains("ALL"))
1482 if ( fOutputEventHandler )
1484 cout << TString('_',78) << endl;
1485 cout << "OutputEventHandler:" << endl;
1486 fOutputEventHandler->Print(" ");
1491 //______________________________________________________________________________
1492 void AliAnalysisManager::ResetAnalysis()
1494 // Reset all execution flags and clean containers.
1498 //______________________________________________________________________________
1499 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1501 // Start analysis having a grid handler.
1502 if (!fGridHandler) {
1503 Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1504 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1508 return StartAnalysis(type, tree, nentries, firstentry);
1511 //______________________________________________________________________________
1512 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1514 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1515 // MIX. Process nentries starting from firstentry
1517 // Backup current directory and make sure gDirectory points to gROOT
1518 TDirectory *cdir = gDirectory;
1521 Error("StartAnalysis","Analysis manager was not initialized !");
1525 if (!CheckTasks()) Fatal("StartAnalysis", "Not all needed libraries were loaded");
1527 printf("StartAnalysis %s\n",GetName());
1528 AliLog::SetGlobalLogLevel(AliLog::kInfo);
1530 fMaxEntries = nentries;
1532 TString anaType = type;
1534 fMode = kLocalAnalysis;
1535 Bool_t runlocalinit = kTRUE;
1536 if (anaType.Contains("file")) {
1537 runlocalinit = kFALSE;
1540 if (anaType.Contains("proof")) fMode = kProofAnalysis;
1541 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1542 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
1544 if (fMode == kGridAnalysis) {
1546 if (!anaType.Contains("terminate")) {
1547 if (!fGridHandler) {
1548 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1549 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1553 // Write analysis manager in the analysis file
1554 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1555 // run local task configuration
1556 TIter nextTask(fTasks);
1557 AliAnalysisTask *task;
1558 while ((task=(AliAnalysisTask*)nextTask())) {
1562 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1563 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1568 // Terminate grid analysis
1569 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1570 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1571 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1572 if (!fGridHandler->MergeOutputs()) {
1573 // Return if outputs could not be merged or if it alien handler
1574 // was configured for offline mode or local testing.
1579 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1580 ImportWrappers(NULL);
1586 SetEventLoop(kFALSE);
1587 // Enable event loop mode if a tree was provided
1588 if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1591 TString ttype = "TTree";
1592 if (tree && tree->IsA() == TChain::Class()) {
1593 chain = (TChain*)tree;
1594 if (!chain || !chain->GetListOfFiles()->First()) {
1595 Error("StartAnalysis", "Cannot process null or empty chain...");
1602 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1603 if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1604 // Initialize locally all tasks (happens for all modes)
1606 AliAnalysisTask *task;
1608 while ((task=(AliAnalysisTask*)next())) {
1612 if (getsysInfo) AliSysInfo::AddStamp("LocalInit_all", 0);
1616 case kLocalAnalysis:
1617 if (!tree && !fGridHandler) {
1618 TIter nextT(fTasks);
1619 // Call CreateOutputObjects for all tasks
1621 Bool_t dirStatus = TH1::AddDirectoryStatus();
1622 while ((task=(AliAnalysisTask*)nextT())) {
1623 TH1::AddDirectory(kFALSE);
1624 task->CreateOutputObjects();
1625 if (!task->CheckPostData()) {
1626 Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
1627 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
1628 ########### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ###########", task->GetName(), task->ClassName());
1630 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1634 TH1::AddDirectory(dirStatus);
1635 if (IsExternalLoop()) {
1636 Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1637 \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1644 fSelector = new AliAnalysisSelector(this);
1645 // Check if a plugin handler is used
1647 // Get the chain from the plugin
1648 TString dataType = "esdTree";
1649 if (fInputEventHandler) {
1650 dataType = fInputEventHandler->GetDataType();
1654 chain = fGridHandler->GetChainForTestMode(dataType);
1656 Error("StartAnalysis", "No chain for test mode. Aborting.");
1659 cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1660 retv = chain->Process(fSelector, "", nentries, firstentry);
1663 // Run tree-based analysis via AliAnalysisSelector
1664 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1665 retv = tree->Process(fSelector, "", nentries, firstentry);
1667 case kProofAnalysis:
1669 // Check if the plugin is used
1671 return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1673 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1674 Error("StartAnalysis", "No PROOF!!! Exiting.");
1678 line = Form("gProof->AddInput((TObject*)%p);", this);
1679 gROOT->ProcessLine(line);
1682 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1683 retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1685 Error("StartAnalysis", "No chain!!! Exiting.");
1692 if (!anaType.Contains("terminate")) {
1693 if (!fGridHandler) {
1694 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1695 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1699 // Write analysis manager in the analysis file
1700 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1701 // Start the analysis via the handler
1702 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1703 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1708 // Terminate grid analysis
1709 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1710 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1711 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1712 if (!fGridHandler->MergeOutputs()) {
1713 // Return if outputs could not be merged or if it alien handler
1714 // was configured for offline mode or local testing.
1719 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1720 ImportWrappers(NULL);
1724 case kMixingAnalysis:
1725 // Run event mixing analysis
1727 Error("StartAnalysis", "Cannot run event mixing without event pool");
1731 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1732 fSelector = new AliAnalysisSelector(this);
1733 while ((chain=fEventPool->GetNextChain())) {
1735 // Call NotifyBinChange for all tasks
1736 while ((task=(AliAnalysisTask*)next()))
1737 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1738 retv = chain->Process(fSelector);
1740 Error("StartAnalysis", "Mixing analysis failed");
1745 PackOutput(fSelector->GetOutputList());
1752 //______________________________________________________________________________
1753 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1755 // Start analysis for this manager on a given dataset. Analysis task can be:
1756 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1758 Error("StartAnalysis","Analysis manager was not initialized !");
1762 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1763 TString anaType = type;
1765 if (!anaType.Contains("proof")) {
1766 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1769 fMode = kProofAnalysis;
1771 SetEventLoop(kTRUE);
1772 // Set the dataset flag
1773 TObject::SetBit(kUseDataSet);
1776 // Start proof analysis using the grid handler
1777 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1778 Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
1781 // Check if the plugin is in test mode
1782 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
1783 dataset = "test_collection";
1785 dataset = fGridHandler->GetProofDataSet();
1789 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1790 Error("StartAnalysis", "No PROOF!!! Exiting.");
1794 // Initialize locally all tasks
1796 AliAnalysisTask *task;
1797 while ((task=(AliAnalysisTask*)next())) {
1801 line = Form("gProof->AddInput((TObject*)%p);", this);
1802 gROOT->ProcessLine(line);
1804 line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1805 dataset, nentries, firstentry);
1806 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1807 retv = (Long_t)gROOT->ProcessLine(line);
1811 //______________________________________________________________________________
1812 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
1814 // Opens according the option the file specified by cont->GetFileName() and changes
1815 // current directory to cont->GetFolderName(). If the file was already opened, it
1816 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
1817 // be optionally ignored.
1818 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1819 TString filename = cont->GetFileName();
1821 if (filename.IsNull()) {
1822 ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
1825 if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
1827 f = mgr->OpenProofFile(cont,option);
1829 // Check first if the file is already opened
1830 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1832 // Check if option "UPDATE" was preserved
1833 TString opt(option);
1835 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1836 ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1838 f = TFile::Open(filename, option);
1841 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1845 // Check for a folder request
1846 TString dir = cont->GetFolderName();
1847 if (!dir.IsNull()) {
1848 if (!f->GetDirectory(dir)) f->mkdir(dir);
1853 ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
1854 cont->SetFile(NULL);
1858 //______________________________________________________________________________
1859 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
1861 // Opens a special output file used in PROOF.
1863 TString filename = cont->GetFileName();
1864 if (cont == fCommonOutput) {
1865 if (fOutputEventHandler) {
1866 if (strlen(extaod)) filename = extaod;
1867 filename = fOutputEventHandler->GetOutputFileName();
1869 else Fatal("OpenProofFile","No output container. Exiting.");
1872 if (fMode!=kProofAnalysis || !fSelector) {
1873 Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
1876 if (fSpecialOutputLocation.Length()) {
1877 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1879 // Check if option "UPDATE" was preserved
1880 TString opt(option);
1882 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1883 ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1885 f = new TFile(filename, option);
1887 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1891 // Check for a folder request
1892 TString dir = cont->GetFolderName();
1894 if (!f->GetDirectory(dir)) f->mkdir(dir);
1899 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1900 cont->SetFile(NULL);
1903 // Check if there is already a proof output file in the output list
1904 TObject *pof = fSelector->GetOutputList()->FindObject(filename);
1906 // Get the actual file
1907 line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
1908 filename = (const char*)gROOT->ProcessLine(line);
1910 printf("File: %s already booked via TProofOutputFile\n", filename.Data());
1912 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1914 Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
1917 // Check if option "UPDATE" was preserved
1918 TString opt(option);
1920 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1921 Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
1923 if (cont->IsRegisterDataset()) {
1924 TString dsetName = filename;
1925 dsetName.ReplaceAll(".root", cont->GetTitle());
1926 dsetName.ReplaceAll(":","_");
1927 if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
1928 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
1930 if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
1931 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
1933 if (fDebug > 1) printf("=== %s\n", line.Data());
1934 gROOT->ProcessLine(line);
1935 line = Form("pf->OpenFile(\"%s\");", option);
1936 gROOT->ProcessLine(line);
1939 gROOT->ProcessLine("pf->Print()");
1940 printf(" == proof file name: %s", f->GetName());
1942 // Add to proof output list
1943 line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
1944 if (fDebug > 1) printf("=== %s\n", line.Data());
1945 gROOT->ProcessLine(line);
1947 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1951 // Check for a folder request
1952 TString dir = cont->GetFolderName();
1953 if (!dir.IsNull()) {
1954 if (!f->GetDirectory(dir)) f->mkdir(dir);
1959 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1960 cont->SetFile(NULL);
1964 //______________________________________________________________________________
1965 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1967 // Execute analysis.
1968 static Long64_t nentries = 0;
1969 static TTree *lastTree = 0;
1970 static TStopwatch *timer = new TStopwatch();
1971 // Only the first call to Process will trigger a true Notify. Other Notify
1972 // coming before is ignored.
1973 if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) {
1974 TObject::SetBit(AliAnalysisManager::kTrueNotify);
1977 if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
1979 if (fTree && (fTree != lastTree)) {
1980 nentries += fTree->GetEntries();
1983 if (!fNcalls) timer->Start();
1984 if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, TMath::Min(fMaxEntries,nentries), timer, kFALSE);
1987 TDirectory *cdir = gDirectory;
1988 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1989 if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
1991 Error("ExecAnalysis", "Analysis manager was not initialized !");
1996 AliAnalysisTask *task;
1997 // Check if the top tree is active.
1999 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2000 AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
2002 // De-activate all tasks
2003 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
2004 AliAnalysisDataContainer *cont = fCommonInput;
2005 if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
2007 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
2011 cont->SetData(fTree); // This will notify all consumers
2012 Long64_t entry = fTree->GetTree()->GetReadEntry();
2014 // Call BeginEvent() for optional input/output and MC services
2015 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
2016 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
2017 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
2019 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2020 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2022 // Execute the tasks
2023 // TIter next1(cont->GetConsumers());
2024 TIter next1(fTopTasks);
2026 while ((task=(AliAnalysisTask*)next1())) {
2028 cout << " Executing task " << task->GetName() << endl;
2030 task->ExecuteTask(option);
2032 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2033 AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
2037 // Call FinishEvent() for optional output and MC services
2038 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2039 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2040 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2041 // Gather system information if requested
2042 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2043 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
2047 // The event loop is not controlled by TSelector
2049 // Call BeginEvent() for optional input/output and MC services
2050 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
2051 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
2052 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
2054 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2055 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2056 TIter next2(fTopTasks);
2057 while ((task=(AliAnalysisTask*)next2())) {
2058 task->SetActive(kTRUE);
2060 cout << " Executing task " << task->GetName() << endl;
2062 task->ExecuteTask(option);
2066 // Call FinishEvent() for optional output and MC services
2067 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2068 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2069 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2070 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2071 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2075 //______________________________________________________________________________
2076 Bool_t AliAnalysisManager::IsPipe(std::ostream &out)
2078 // Check if the stdout is connected to a pipe (C.Holm)
2079 Bool_t ispipe = kFALSE;
2080 out.seekp(0, std::ios_base::cur);
2083 if (errno == ESPIPE) ispipe = kTRUE;
2088 //______________________________________________________________________________
2089 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2091 // Set the input event handler and create a container for it.
2092 fInputEventHandler = handler;
2093 if (!fCommonInput) fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2096 //______________________________________________________________________________
2097 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2099 // Set the input event handler and create a container for it.
2100 fOutputEventHandler = handler;
2101 if (!fCommonOutput) fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2102 fCommonOutput->SetSpecialOutput();
2105 //______________________________________________________________________________
2106 void AliAnalysisManager::SetDebugLevel(UInt_t level)
2108 // Set verbosity of the analysis manager. If the progress bar is used, the call is ignored
2109 if (TObject::TestBit(kUseProgressBar)) {
2110 Info("SetDebugLevel","Ignored. Disable the progress bar first.");
2116 //______________________________________________________________________________
2117 void AliAnalysisManager::SetUseProgressBar(Bool_t flag, Int_t freq)
2119 // Enable a text mode progress bar. Resets debug level to 0.
2120 Info("SetUseProgressBar", "Progress bar enabled, updated every %d events.\n ### NOTE: Debug level reset to 0 ###", freq);
2121 TObject::SetBit(kUseProgressBar,flag);
2122 fPBUpdateFreq = freq;
2126 //______________________________________________________________________________
2127 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2129 // This method is used externally to register output files which are not
2130 // connected to any output container, so that the manager can properly register,
2131 // retrieve or merge them when running in distributed mode. The file names are
2132 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2133 if (fExtraFiles.Contains(fname)) return;
2134 if (fExtraFiles.Length()) fExtraFiles += " ";
2135 fExtraFiles += fname;
2138 //______________________________________________________________________________
2139 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2141 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2145 TObject *pof = source->FindObject(filename);
2146 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2147 Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2150 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2151 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2152 TString clientUrl(chUrl);
2153 TString fullPath_str(fullPath);
2154 if (clientUrl.Contains("localhost")){
2155 TObjArray* array = fullPath_str.Tokenize ( "//" );
2156 TObjString *strobj = ( TObjString *)array->At(1);
2157 TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2158 TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2159 fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2160 fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2161 if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2165 else if (clientUrl.Contains("__lite__")) {
2166 // Special case for ProofLite environement - get file info and copy.
2167 gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2168 fullPath_str = Form("%s/%s", tmp, fullPath);
2171 Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2172 Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename);
2174 Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2178 //______________________________________________________________________________
2179 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2181 // Fill analysis type in the provided string.
2183 case kLocalAnalysis:
2186 case kProofAnalysis:
2192 case kMixingAnalysis:
2197 //______________________________________________________________________________
2198 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2200 // Validate all output files.
2201 TIter next(fOutputs);
2202 AliAnalysisDataContainer *output;
2203 TDirectory *cdir = gDirectory;
2204 TString openedFiles;
2205 while ((output=(AliAnalysisDataContainer*)next())) {
2206 if (output->IsRegisterDataset()) continue;
2207 TString filename = output->GetFileName();
2208 if (filename == "default") {
2209 if (!fOutputEventHandler) continue;
2210 filename = fOutputEventHandler->GetOutputFileName();
2211 // Main AOD may not be there
2212 if (gSystem->AccessPathName(filename)) continue;
2214 // Check if the file is closed
2215 if (openedFiles.Contains(filename)) continue;;
2216 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2218 Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2219 // Clear file list to release object ownership to user.
2223 file = TFile::Open(filename);
2224 if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2225 Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2230 openedFiles += filename;
2237 //______________________________________________________________________________
2238 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2240 // Implements a nice text mode progress bar.
2241 static Long64_t icount = 0;
2242 static TString oname;
2243 static TString nname;
2244 static Long64_t ocurrent = 0;
2245 static Long64_t osize = 0;
2246 static Int_t oseconds = 0;
2247 static TStopwatch *owatch = 0;
2248 static Bool_t oneoftwo = kFALSE;
2249 static Int_t nrefresh = 0;
2250 static Int_t nchecks = 0;
2251 static char lastChar = 0;
2252 const char symbol[4] = {'-','\\','|','/'};
2254 if (!lastChar) lastChar = (IsPipe(std::cerr))?'\r':'\n';
2260 ocurrent = TMath::Abs(current);
2261 osize = TMath::Abs(size);
2262 if (ocurrent > osize) ocurrent=osize;
2267 if ((current % fPBUpdateFreq) != 0) return;
2269 char progress[11] = " ";
2270 Int_t ichar = icount%4;
2275 if (owatch && !last) {
2277 time = owatch->RealTime();
2278 seconds = int(time) % 60;
2279 minutes = (int(time) / 60) % 60;
2280 hours = (int(time) / 60 / 60);
2282 if (oseconds==seconds) {
2286 oneoftwo = !oneoftwo;
2290 if (refresh && oneoftwo) {
2292 if (nchecks <= 0) nchecks = nrefresh+1;
2293 Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2294 oname = Form(" == %d%% ==", pctdone);
2296 Double_t percent = 100.0*ocurrent/osize;
2297 Int_t nchar = Int_t(percent/10);
2298 if (nchar>10) nchar=10;
2300 for (i=0; i<nchar; i++) progress[i] = '=';
2301 progress[nchar] = symbol[ichar];
2302 for (i=nchar+1; i<10; i++) progress[i] = ' ';
2303 progress[10] = '\0';
2306 if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2307 else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2308 else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2310 Int_t full = Int_t(ocurrent > 0 ?
2311 time * (float(osize)/ocurrent) + .5 :
2313 Int_t remain = Int_t(full - time);
2314 Int_t rsec = remain % 60;
2315 Int_t rmin = (remain / 60) % 60;
2316 Int_t rhour = (remain / 60 / 60);
2317 fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d ETA %.2d:%.2d:%.2d%c",
2318 percent, hours, minutes, seconds, rhour, rmin, rsec, lastChar);
2320 else fprintf(stderr, "[%6.2f %%]%c", percent, lastChar);
2321 if (refresh && oneoftwo) oname = nname;
2322 if (owatch) owatch->Continue();
2331 fprintf(stderr, "\n");
2335 //______________________________________________________________________________
2336 void AliAnalysisManager::DoLoadBranch(const char *name)
2338 // Get tree and load branch if needed.
2339 static Long64_t crtEntry = -100;
2341 if (fAutoBranchHandling || !fTree)
2344 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2346 br = fTree->GetBranch(name);
2348 Error("DoLoadBranch", "Could not find branch %s",name);
2353 if (br->GetReadEntry()==fCurrentEntry) return;
2354 Int_t ret = br->GetEntry(GetCurrentEntry());
2356 Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2357 if (crtEntry != fCurrentEntry) {
2358 CountEvent(1,0,1,0);
2359 crtEntry = fCurrentEntry;
2362 if (crtEntry != fCurrentEntry) {
2363 CountEvent(1,1,0,0);
2364 crtEntry = fCurrentEntry;
2369 //______________________________________________________________________________
2370 void AliAnalysisManager::AddStatisticsTask(UInt_t offlineMask)
2372 // Add the statistics task to the manager.
2374 Info("AddStatisticsTask", "Already added");
2377 TString line = Form("AliAnalysisTaskStat::AddToManager(%u);", offlineMask);
2378 gROOT->ProcessLine(line);
2381 //______________________________________________________________________________
2382 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2384 // Bookkeep current event;
2385 if (!fStatistics) return;
2386 fStatistics->AddInput(ninput);
2387 fStatistics->AddProcessed(nprocessed);
2388 fStatistics->AddFailed(nfailed);
2389 fStatistics->AddAccepted(naccepted);
2392 //______________________________________________________________________________
2393 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2395 // Add a line in the statistics message. If available, the statistics message is written
2396 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2398 if (!strlen(line)) return;
2399 if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2400 fStatisticsMsg += line;
2403 //______________________________________________________________________________
2404 void AliAnalysisManager::WriteStatisticsMsg(Int_t)
2406 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2407 static Bool_t done = kFALSE;
2410 if (!fStatistics) return;
2412 AddStatisticsMsg(Form("Number of input events: %lld",fStatistics->GetNinput()));
2413 AddStatisticsMsg(Form("Number of processed events: %lld",fStatistics->GetNprocessed()));
2414 AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2415 AddStatisticsMsg(Form("Number of accepted events for mask %s: %lld", AliAnalysisStatistics::GetMaskAsString(fStatistics->GetOfflineMask()), fStatistics->GetNaccepted()));
2416 out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2417 fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2418 fStatistics->GetNaccepted()), ios::out);
2419 out << fStatisticsMsg << endl;
2423 //______________________________________________________________________________
2424 const char* AliAnalysisManager::GetOADBPath()
2426 // returns the path of the OADB
2427 // this static function just depends on environment variables
2429 static TString oadbPath;
2431 if (gSystem->Getenv("OADB_PATH"))
2432 oadbPath = gSystem->Getenv("OADB_PATH");
2433 else if (gSystem->Getenv("ALICE_ROOT"))
2434 oadbPath.Form("%s/OADB", gSystem->Getenv("ALICE_ROOT"));
2436 ::Fatal("AliAnalysisManager::GetOADBPath", "Cannot figure out AODB path. Define ALICE_ROOT or OADB_PATH!");