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"
30 #include <Riostream.h>
36 #include <TMethodCall.h>
41 #include <TStopwatch.h>
43 #include "AliAnalysisSelector.h"
44 #include "AliAnalysisGrid.h"
45 #include "AliAnalysisTask.h"
46 #include "AliAnalysisDataContainer.h"
47 #include "AliAnalysisDataSlot.h"
48 #include "AliVEventHandler.h"
49 #include "AliVEventPool.h"
50 #include "AliSysInfo.h"
51 #include "AliAnalysisStatistics.h"
53 ClassImp(AliAnalysisManager)
55 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
56 TString AliAnalysisManager::fgCommonFileName = "";
58 //______________________________________________________________________________
59 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
62 fInputEventHandler(NULL),
63 fOutputEventHandler(NULL),
64 fMCtruthEventHandler(NULL),
68 fMode(kLocalAnalysis),
72 fSpecialOutputLocation(""),
85 fAutoBranchHandling(kTRUE),
93 // Default constructor.
94 fgAnalysisManager = this;
95 fgCommonFileName = "AnalysisResults.root";
96 fTasks = new TObjArray();
97 fTopTasks = new TObjArray();
98 fZombies = new TObjArray();
99 fContainers = new TObjArray();
100 fInputs = new TObjArray();
101 fOutputs = new TObjArray();
102 fParamCont = new TObjArray();
104 TObject::SetObjectStat(kFALSE);
107 //______________________________________________________________________________
108 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
111 fInputEventHandler(NULL),
112 fOutputEventHandler(NULL),
113 fMCtruthEventHandler(NULL),
118 fInitOK(other.fInitOK),
119 fIsRemote(other.fIsRemote),
120 fDebug(other.fDebug),
121 fSpecialOutputLocation(""),
134 fAutoBranchHandling(other.fAutoBranchHandling),
137 fNcalls(other.fNcalls),
138 fStatisticsMsg(other.fStatisticsMsg),
139 fRequestedBranches(other.fRequestedBranches),
140 fStatistics(other.fStatistics)
143 fTasks = new TObjArray(*other.fTasks);
144 fTopTasks = new TObjArray(*other.fTopTasks);
145 fZombies = new TObjArray(*other.fZombies);
146 fContainers = new TObjArray(*other.fContainers);
147 fInputs = new TObjArray(*other.fInputs);
148 fOutputs = new TObjArray(*other.fOutputs);
149 fParamCont = new TObjArray(*other.fParamCont);
150 fgCommonFileName = "AnalysisResults.root";
151 fgAnalysisManager = this;
152 TObject::SetObjectStat(kFALSE);
155 //______________________________________________________________________________
156 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
159 if (&other != this) {
160 TNamed::operator=(other);
161 fInputEventHandler = other.fInputEventHandler;
162 fOutputEventHandler = other.fOutputEventHandler;
163 fMCtruthEventHandler = other.fMCtruthEventHandler;
164 fEventPool = other.fEventPool;
167 fNSysInfo = other.fNSysInfo;
169 fInitOK = other.fInitOK;
170 fIsRemote = other.fIsRemote;
171 fDebug = other.fDebug;
172 fTasks = new TObjArray(*other.fTasks);
173 fTopTasks = new TObjArray(*other.fTopTasks);
174 fZombies = new TObjArray(*other.fZombies);
175 fContainers = new TObjArray(*other.fContainers);
176 fInputs = new TObjArray(*other.fInputs);
177 fOutputs = new TObjArray(*other.fOutputs);
178 fParamCont = new TObjArray(*other.fParamCont);
180 fCommonOutput = NULL;
183 fExtraFiles = other.fExtraFiles;
184 fgCommonFileName = "AnalysisResults.root";
185 fgAnalysisManager = this;
186 fAutoBranchHandling = other.fAutoBranchHandling;
187 fTable.Clear("nodelete");
188 fRunFromPath = other.fRunFromPath;
189 fNcalls = other. fNcalls;
190 fStatisticsMsg = other.fStatisticsMsg;
191 fRequestedBranches = other.fRequestedBranches;
192 fStatistics = other.fStatistics;
197 //______________________________________________________________________________
198 AliAnalysisManager::~AliAnalysisManager()
201 if (fTasks) {fTasks->Delete(); delete fTasks;}
202 if (fTopTasks) delete fTopTasks;
203 if (fZombies) delete fZombies;
204 if (fContainers) {fContainers->Delete(); delete fContainers;}
205 if (fInputs) delete fInputs;
206 if (fOutputs) delete fOutputs;
207 if (fParamCont) delete fParamCont;
208 if (fGridHandler) delete fGridHandler;
209 if (fInputEventHandler) delete fInputEventHandler;
210 if (fOutputEventHandler) delete fOutputEventHandler;
211 if (fMCtruthEventHandler) delete fMCtruthEventHandler;
212 if (fEventPool) delete fEventPool;
213 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
214 TObject::SetObjectStat(kTRUE);
217 //______________________________________________________________________________
218 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
220 // Read one entry of the tree or a whole branch.
221 fCurrentEntry = entry;
222 if (!fAutoBranchHandling)
224 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : -1;
227 //______________________________________________________________________________
228 Int_t AliAnalysisManager::GetRunFromAlienPath(const char *path)
230 // Attempt to extract run number from input data path. Works only for paths to
231 // alice data in alien.
232 // sim: /alice/sim/<production>/run_no/...
233 // data: /alice/data/year/period/000run_no/... (ESD or AOD)
237 Int_t index = s.Index("/alice/sim");
239 for (Int_t i=0; i<3; i++) {
240 index = s.Index("/", index+1);
241 if (index<0) return 0;
246 index = s.Index("/alice/data");
248 for (Int_t i=0; i<4; i++) {
249 index = s.Index("/", index+1);
250 if (index<0) return 0;
258 //______________________________________________________________________________
259 Bool_t AliAnalysisManager::Init(TTree *tree)
261 // The Init() function is called when the selector needs to initialize
262 // a new tree or chain. Typically here the branch addresses of the tree
263 // will be set. It is normaly not necessary to make changes to the
264 // generated code, but the routine can be extended by the user if needed.
265 // Init() will be called many times when running with PROOF.
266 Bool_t init = kFALSE;
267 if (!tree) return kFALSE; // Should not happen - protected in selector caller
269 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
271 // Call InitTree of EventHandler
272 if (fOutputEventHandler) {
273 if (fMode == kProofAnalysis) {
274 init = fOutputEventHandler->Init(0x0, "proof");
276 init = fOutputEventHandler->Init(0x0, "local");
279 Error("Init", "Output event handler failed to initialize");
284 if (fInputEventHandler) {
285 if (fMode == kProofAnalysis) {
286 init = fInputEventHandler->Init(tree, "proof");
288 init = fInputEventHandler->Init(tree, "local");
291 Error("Init", "Input event handler failed to initialize tree");
295 // If no input event handler we need to get the tree once
297 if(!tree->GetTree()) {
298 Long64_t readEntry = tree->LoadTree(0);
299 if (readEntry == -2) {
300 Error("Init", "Input tree has no entry. Exiting");
306 if (fMCtruthEventHandler) {
307 if (fMode == kProofAnalysis) {
308 init = fMCtruthEventHandler->Init(0x0, "proof");
310 init = fMCtruthEventHandler->Init(0x0, "local");
313 Error("Init", "MC event handler failed to initialize");
318 if (!fInitOK) InitAnalysis();
319 if (!fInitOK) return kFALSE;
322 AliAnalysisDataContainer *top = fCommonInput;
323 if (!top) top = (AliAnalysisDataContainer*)fInputs->At(0);
325 Error("Init","No top input container !");
329 CheckBranches(kFALSE);
331 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
336 //______________________________________________________________________________
337 void AliAnalysisManager::SlaveBegin(TTree *tree)
339 // The SlaveBegin() function is called after the Begin() function.
340 // When running with PROOF SlaveBegin() is called on each slave server.
341 // The tree argument is deprecated (on PROOF 0 is passed).
342 if (fDebug > 1) printf("->AliAnalysisManager::SlaveBegin()\n");
343 static Bool_t isCalled = kFALSE;
344 Bool_t init = kFALSE;
345 Bool_t initOK = kTRUE;
347 TDirectory *curdir = gDirectory;
348 // Call SlaveBegin only once in case of mixing
349 if (isCalled && fMode==kMixingAnalysis) return;
351 // Call Init of EventHandler
352 if (fOutputEventHandler) {
353 if (fMode == kProofAnalysis) {
354 // Merging AOD's in PROOF via TProofOutputFile
355 if (fDebug > 1) printf(" Initializing AOD output file %s...\n", fOutputEventHandler->GetOutputFileName());
356 init = fOutputEventHandler->Init("proof");
357 if (!init) msg = "Failed to initialize output handler on worker";
359 init = fOutputEventHandler->Init("local");
360 if (!init) msg = "Failed to initialize output handler";
363 if (!fSelector) Error("SlaveBegin", "Selector not set");
364 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
367 if (fInputEventHandler) {
368 fInputEventHandler->SetInputTree(tree);
369 if (fMode == kProofAnalysis) {
370 init = fInputEventHandler->Init("proof");
371 if (!init) msg = "Failed to initialize input handler on worker";
373 init = fInputEventHandler->Init("local");
374 if (!init) msg = "Failed to initialize input handler";
377 if (!fSelector) Error("SlaveBegin", "Selector not set");
378 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
381 if (fMCtruthEventHandler) {
382 if (fMode == kProofAnalysis) {
383 init = fMCtruthEventHandler->Init("proof");
384 if (!init) msg = "Failed to initialize MC handler on worker";
386 init = fMCtruthEventHandler->Init("local");
387 if (!init) msg = "Failed to initialize MC handler";
390 if (!fSelector) Error("SlaveBegin", "Selector not set");
391 else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
393 if (curdir) curdir->cd();
397 AliAnalysisTask *task;
398 // Call CreateOutputObjects for all tasks
399 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
400 Bool_t dirStatus = TH1::AddDirectoryStatus();
402 while ((task=(AliAnalysisTask*)next())) {
404 // Start with memory as current dir and make sure by default histograms do not get attached to files.
405 TH1::AddDirectory(kFALSE);
406 task->CreateOutputObjects();
407 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
410 TH1::AddDirectory(dirStatus);
411 if (curdir) curdir->cd();
412 if (fDebug > 1) printf("<-AliAnalysisManager::SlaveBegin()\n");
415 //______________________________________________________________________________
416 Bool_t AliAnalysisManager::Notify()
418 // The Notify() function is called when a new file is opened. This
419 // can be either for a new TTree in a TChain or when when a new TTree
420 // is started when using PROOF. It is normaly not necessary to make changes
421 // to the generated code, but the routine can be extended by the
422 // user if needed. The return value is currently not used.
423 if (!fTree) return kFALSE;
425 fTable.Clear("nodelete"); // clearing the hash table may not be needed -> C.L.
426 if (fMode == kProofAnalysis) fIsRemote = kTRUE;
428 TFile *curfile = fTree->GetCurrentFile();
430 Error("Notify","No current file");
434 if (fDebug > 1) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
435 Int_t run = AliAnalysisManager::GetRunFromAlienPath(curfile->GetName());
436 if (run) SetRunFromPath(run);
437 if (fDebug > 1) printf(" ### run found from path: %d\n", run);
439 AliAnalysisTask *task;
441 // Call Notify of the event handlers
442 if (fInputEventHandler) {
443 fInputEventHandler->Notify(curfile->GetName());
446 if (fOutputEventHandler) {
447 fOutputEventHandler->Notify(curfile->GetName());
450 if (fMCtruthEventHandler) {
451 fMCtruthEventHandler->Notify(curfile->GetName());
454 // Call Notify for all tasks
455 while ((task=(AliAnalysisTask*)next()))
458 if (fDebug > 1) printf("<-AliAnalysisManager::Notify()\n");
462 //______________________________________________________________________________
463 Bool_t AliAnalysisManager::Process(Long64_t entry)
465 // The Process() function is called for each entry in the tree (or possibly
466 // keyed object in the case of PROOF) to be processed. The entry argument
467 // specifies which entry in the currently loaded tree is to be processed.
468 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
469 // to read either all or the required parts of the data. When processing
470 // keyed objects with PROOF, the object is already loaded and is available
471 // via the fObject pointer.
473 // This function should contain the "body" of the analysis. It can contain
474 // simple or elaborate selection criteria, run algorithms on the data
475 // of the event and typically fill histograms.
477 // WARNING when a selector is used with a TChain, you must use
478 // the pointer to the current TTree to call GetEntry(entry).
479 // The entry is always the local entry number in the current tree.
480 // Assuming that fChain is the pointer to the TChain being processed,
481 // use fChain->GetTree()->GetEntry(entry).
482 if (fDebug > 1) printf("->AliAnalysisManager::Process(%lld)\n", entry);
484 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
485 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
486 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
490 if (fInputEventHandler) fInputEventHandler ->GetEntry();
493 if (fDebug > 1) printf("<-AliAnalysisManager::Process()\n");
497 //______________________________________________________________________________
498 void AliAnalysisManager::PackOutput(TList *target)
500 // Pack all output data containers in the output list. Called at SlaveTerminate
501 // stage in PROOF case for each slave.
502 if (fDebug > 1) printf("->AliAnalysisManager::PackOutput()\n");
504 Error("PackOutput", "No target. Exiting.");
507 TDirectory *cdir = gDirectory;
509 if (fInputEventHandler) fInputEventHandler ->Terminate();
510 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
511 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
514 // Call FinishTaskOutput() for each event loop task (not called for
515 // post-event loop tasks - use Terminate() fo those)
516 TIter nexttask(fTasks);
517 AliAnalysisTask *task;
518 while ((task=(AliAnalysisTask*)nexttask())) {
519 if (!task->IsPostEventLoop()) {
520 if (fDebug > 1) printf("->FinishTaskOutput: task %s\n", task->GetName());
521 task->FinishTaskOutput();
523 if (fDebug > 1) printf("<-FinishTaskOutput: task %s\n", task->GetName());
526 // Write statistics message on the workers.
527 WriteStatisticsMsg(fNcalls);
529 if (fMode == kProofAnalysis) {
530 TIter next(fOutputs);
531 AliAnalysisDataContainer *output;
532 Bool_t isManagedByHandler = kFALSE;
535 while ((output=(AliAnalysisDataContainer*)next())) {
536 // Do not consider outputs of post event loop tasks
537 isManagedByHandler = kFALSE;
538 if (output->GetProducer() && output->GetProducer()->IsPostEventLoop()) continue;
539 const char *filename = output->GetFileName();
540 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
541 isManagedByHandler = kTRUE;
542 printf("#### Handler output. Extra: %s\n", fExtraFiles.Data());
543 filename = fOutputEventHandler->GetOutputFileName();
545 // Check if data was posted to this container. If not, issue an error.
546 if (!output->GetData() && !isManagedByHandler) {
547 Error("PackOutput", "No data for output container %s. Forgot to PostData ?", output->GetName());
550 if (!output->IsSpecialOutput()) {
552 if (strlen(filename) && !isManagedByHandler) {
553 // Backup current folder
554 TDirectory *opwd = gDirectory;
555 // File resident outputs.
556 // Check first if the file exists.
557 TString openoption = "RECREATE";
558 Bool_t firsttime = kTRUE;
559 if (filestmp.FindObject(output->GetFileName())) {
562 filestmp.Add(new TNamed(output->GetFileName(),""));
564 if (!gSystem->AccessPathName(output->GetFileName()) && !firsttime) openoption = "UPDATE";
565 // TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
566 // Save data to file, then close.
567 if (output->GetData()->InheritsFrom(TCollection::Class())) {
568 // If data is a collection, we set the name of the collection
569 // as the one of the container and we save as a single key.
570 TCollection *coll = (TCollection*)output->GetData();
571 coll->SetName(output->GetName());
572 // coll->Write(output->GetName(), TObject::kSingleKey);
574 if (output->GetData()->InheritsFrom(TTree::Class())) {
575 TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
576 // Save data to file, then close.
577 TTree *tree = (TTree*)output->GetData();
578 // Check if tree is in memory
579 if (tree->GetDirectory()==gROOT) tree->SetDirectory(gDirectory);
583 // output->GetData()->Write();
586 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
588 // printf(" file %s listing content:\n", filename);
591 // Clear file list to release object ownership to user.
594 output->SetFile(NULL);
595 // Restore current directory
596 if (opwd) opwd->cd();
598 // Memory-resident outputs
599 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
601 AliAnalysisDataWrapper *wrap = 0;
602 if (isManagedByHandler) {
603 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
604 wrap->SetName(output->GetName());
606 else wrap =output->ExportData();
607 // Output wrappers must NOT delete data after merging - the user owns them
608 wrap->SetDeleteData(kFALSE);
611 // Special outputs. The file must be opened and connected to the container.
612 TDirectory *opwd = gDirectory;
613 TFile *file = output->GetFile();
615 AliAnalysisTask *producer = output->GetProducer();
617 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
618 output->GetFileName(), output->GetName(), producer->ClassName());
621 TString outFilename = file->GetName();
622 if (fDebug > 1) printf("PackOutput %s: special output\n", output->GetName());
623 if (isManagedByHandler) {
624 // Terminate IO for files managed by the output handler
625 // file->Write() moved to AOD handler (A.G. 11.01.10)
626 // if (file) file->Write();
627 if (file && fDebug > 2) {
628 printf(" handled file %s listing content:\n", file->GetName());
631 fOutputEventHandler->TerminateIO();
634 // Release object ownership to users after writing data to file
635 if (output->GetData()->InheritsFrom(TCollection::Class())) {
636 // If data is a collection, we set the name of the collection
637 // as the one of the container and we save as a single key.
638 TCollection *coll = (TCollection*)output->GetData();
639 coll->SetName(output->GetName());
640 coll->Write(output->GetName(), TObject::kSingleKey);
642 if (output->GetData()->InheritsFrom(TTree::Class())) {
643 TTree *tree = (TTree*)output->GetData();
644 tree->SetDirectory(file);
647 output->GetData()->Write();
651 printf(" file %s listing content:\n", output->GetFileName());
654 // Clear file list to release object ownership to user.
657 output->SetFile(NULL);
659 // Restore current directory
660 if (opwd) opwd->cd();
661 // Check if a special output location was provided or the output files have to be merged
662 if (strlen(fSpecialOutputLocation.Data())) {
663 TString remote = fSpecialOutputLocation;
665 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
666 if (remote.BeginsWith("alien:")) {
667 gROOT->ProcessLine("TGrid::Connect(\"alien:\", gProofServ->GetUser());");
668 remote += outFilename;
669 remote.ReplaceAll(".root", Form("_%d.root", gid));
671 remote += Form("%s_%d_", gSystem->HostName(), gid);
672 remote += outFilename;
675 Info("PackOutput", "Output file for container %s to be copied \n at: %s. No merging.",
676 output->GetName(), remote.Data());
677 TFile::Cp ( outFilename.Data(), remote.Data() );
678 // Copy extra outputs
679 if (fExtraFiles.Length() && isManagedByHandler) {
680 TObjArray *arr = fExtraFiles.Tokenize(" ");
682 TIter nextfilename(arr);
683 while ((os=(TObjString*)nextfilename())) {
684 outFilename = os->GetString();
685 remote = fSpecialOutputLocation;
687 if (remote.BeginsWith("alien://")) {
688 remote += outFilename;
689 remote.ReplaceAll(".root", Form("_%d.root", gid));
691 remote += Form("%s_%d_", gSystem->HostName(), gid);
692 remote += outFilename;
695 Info("PackOutput", "Extra AOD file %s to be copied \n at: %s. No merging.",
696 outFilename.Data(), remote.Data());
697 TFile::Cp ( outFilename.Data(), remote.Data() );
702 // No special location specified-> use TProofOutputFile as merging utility
703 // The file at this output slot must be opened in CreateOutputObjects
704 if (fDebug > 1) printf(" File for container %s to be merged via file merger...\n", output->GetName());
710 if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
713 //______________________________________________________________________________
714 void AliAnalysisManager::ImportWrappers(TList *source)
716 // Import data in output containers from wrappers coming in source.
717 if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
718 TIter next(fOutputs);
719 AliAnalysisDataContainer *cont;
720 AliAnalysisDataWrapper *wrap;
722 Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
723 TDirectory *cdir = gDirectory;
724 while ((cont=(AliAnalysisDataContainer*)next())) {
726 if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
727 if (cont->IsRegisterDataset()) continue;
728 const char *filename = cont->GetFileName();
729 Bool_t isManagedByHandler = kFALSE;
730 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
731 isManagedByHandler = kTRUE;
732 filename = fOutputEventHandler->GetOutputFileName();
734 if (cont->IsSpecialOutput() || inGrid) {
735 if (strlen(fSpecialOutputLocation.Data())) continue;
736 // Copy merged file from PROOF scratch space.
737 // In case of grid the files are already in the current directory.
739 if (isManagedByHandler && fExtraFiles.Length()) {
740 // Copy extra registered dAOD files.
741 TObjArray *arr = fExtraFiles.Tokenize(" ");
743 TIter nextfilename(arr);
744 while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
747 if (!GetFileFromWrapper(filename, source)) continue;
749 // Normally we should connect data from the copied file to the
750 // corresponding output container, but it is not obvious how to do this
751 // automatically if several objects in file...
752 TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
753 if (!f) f = TFile::Open(filename, "READ");
755 Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
759 // Cd to the directory pointed by the container
760 TString folder = cont->GetFolderName();
761 if (!folder.IsNull()) f->cd(folder);
762 // Try to fetch first an object having the container name.
763 obj = gDirectory->Get(cont->GetName());
765 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",
766 cont->GetType()->GetName(), cont->GetName(), filename, cont->GetFolderName(), cont->GetName());
769 wrap = new AliAnalysisDataWrapper(obj);
770 wrap->SetDeleteData(kFALSE);
772 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
774 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
779 printf(" Importing data for container %s\n", cont->GetName());
780 if (strlen(filename)) printf(" -> file %s\n", filename);
783 cont->ImportData(wrap);
785 if (cdir) cdir->cd();
786 if (fDebug > 1) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
789 //______________________________________________________________________________
790 void AliAnalysisManager::UnpackOutput(TList *source)
792 // Called by AliAnalysisSelector::Terminate only on the client.
793 if (fDebug > 1) printf("->AliAnalysisManager::UnpackOutput()\n");
795 Error("UnpackOutput", "No target. Exiting.");
798 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
800 if (fMode == kProofAnalysis) ImportWrappers(source);
802 TIter next(fOutputs);
803 AliAnalysisDataContainer *output;
804 while ((output=(AliAnalysisDataContainer*)next())) {
805 if (!output->GetData()) continue;
806 // Check if there are client tasks that run post event loop
807 if (output->HasConsumers()) {
808 // Disable event loop semaphore
809 output->SetPostEventLoop(kTRUE);
810 TObjArray *list = output->GetConsumers();
811 Int_t ncons = list->GetEntriesFast();
812 for (Int_t i=0; i<ncons; i++) {
813 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
814 task->CheckNotify(kTRUE);
815 // If task is active, execute it
816 if (task->IsPostEventLoop() && task->IsActive()) {
817 if (fDebug > 1) printf("== Executing post event loop task %s\n", task->GetName());
823 if (fDebug > 1) printf("<-AliAnalysisManager::UnpackOutput()\n");
826 //______________________________________________________________________________
827 void AliAnalysisManager::Terminate()
829 // The Terminate() function is the last function to be called during
830 // a query. It always runs on the client, it can be used to present
831 // the results graphically.
832 if (fDebug > 1) printf("->AliAnalysisManager::Terminate()\n");
833 TDirectory *cdir = gDirectory;
835 AliAnalysisTask *task;
836 AliAnalysisDataContainer *output;
839 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
840 // Call Terminate() for tasks
842 while (!IsSkipTerminate() && (task=(AliAnalysisTask*)next())) {
843 // Save all the canvases produced by the Terminate
844 TString pictname = Form("%s_%s", task->GetName(), task->ClassName());
848 AliSysInfo::AddStamp(Form("%s_TERMINATE",task->ClassName()),0, itask, 2);
850 if (TObject::TestBit(kSaveCanvases)) {
851 if (!gROOT->IsBatch()) {
852 if (fDebug>1) printf("Waiting 5 sec for %s::Terminate() to finish drawing ...\n", task->ClassName());
854 while (timer.CpuTime()<5) {
856 gSystem->ProcessEvents();
859 Int_t iend = gROOT->GetListOfCanvases()->GetEntries();
860 if (iend==0) continue;
862 for (Int_t ipict=0; ipict<iend; ipict++) {
863 canvas = (TCanvas*)gROOT->GetListOfCanvases()->At(ipict);
864 if (!canvas) continue;
865 canvas->SaveAs(Form("%s_%02d.gif", pictname.Data(),ipict));
867 gROOT->GetListOfCanvases()->Delete();
871 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
872 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
873 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
875 TObjArray *allOutputs = new TObjArray();
877 for (icont=0; icont<fOutputs->GetEntriesFast(); icont++) allOutputs->Add(fOutputs->At(icont));
878 if (!IsSkipTerminate())
879 for (icont=0; icont<fParamCont->GetEntriesFast(); icont++) allOutputs->Add(fParamCont->At(icont));
880 TIter next1(allOutputs);
881 TString handlerFile = "";
882 TString extraOutputs = "";
883 if (fOutputEventHandler) {
884 handlerFile = fOutputEventHandler->GetOutputFileName();
885 extraOutputs = fOutputEventHandler->GetExtraOutputs();
889 while ((output=(AliAnalysisDataContainer*)next1())) {
890 // Special outputs or grid files have the files already closed and written.
892 if (fMode == kGridAnalysis && icont<=fOutputs->GetEntriesFast()) continue;
893 if (fMode == kProofAnalysis) {
894 if (output->IsSpecialOutput() || output->IsRegisterDataset()) continue;
896 const char *filename = output->GetFileName();
897 TString openoption = "RECREATE";
898 if (!(strcmp(filename, "default"))) continue;
899 if (!strlen(filename)) continue;
900 if (!output->GetData()) continue;
901 TDirectory *opwd = gDirectory;
902 TFile *file = output->GetFile();
903 if (!file) file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
905 //if (handlerFile == filename && !gSystem->AccessPathName(filename)) openoption = "UPDATE";
906 Bool_t firsttime = kTRUE;
907 if (filestmp.FindObject(filename) || extraOutputs.Contains(filename)) {
910 filestmp.Add(new TNamed(filename,""));
912 if (!gSystem->AccessPathName(filename) && !firsttime) openoption = "UPDATE";
913 if (fDebug>1) printf("Opening file: %s option=%s\n",filename, openoption.Data());
914 file = new TFile(filename, openoption);
916 if (fDebug>1) printf("File <%s> already opened with option: <%s> \n", filename, file->GetOption());
917 openoption = file->GetOption();
918 if (openoption == "READ") {
919 if (fDebug>1) printf("...reopening in UPDATE mode\n");
920 file->ReOpen("UPDATE");
923 if (file->IsZombie()) {
924 Error("Terminate", "Cannot open output file %s", filename);
927 output->SetFile(file);
929 // Check for a folder request
930 TString dir = output->GetFolderName();
932 if (!file->GetDirectory(dir)) file->mkdir(dir);
935 if (fDebug > 1) printf("...writing container %s to file %s:%s\n", output->GetName(), file->GetName(), output->GetFolderName());
936 if (output->GetData()->InheritsFrom(TCollection::Class())) {
937 // If data is a collection, we set the name of the collection
938 // as the one of the container and we save as a single key.
939 TCollection *coll = (TCollection*)output->GetData();
940 coll->SetName(output->GetName());
941 coll->Write(output->GetName(), TObject::kSingleKey);
943 if (output->GetData()->InheritsFrom(TTree::Class())) {
944 TTree *tree = (TTree*)output->GetData();
945 tree->SetDirectory(gDirectory);
948 output->GetData()->Write();
951 if (opwd) opwd->cd();
955 while ((output=(AliAnalysisDataContainer*)next1())) {
956 // Close all files at output
957 TDirectory *opwd = gDirectory;
958 if (output->GetFile()) {
959 // Clear file list to release object ownership to user.
960 // output->GetFile()->Clear();
961 output->GetFile()->Close();
962 output->SetFile(NULL);
963 // Copy merged outputs in alien if requested
964 if (fSpecialOutputLocation.Length() &&
965 fSpecialOutputLocation.BeginsWith("alien://")) {
966 Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data());
967 TFile::Cp(output->GetFile()->GetName(),
968 Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
971 if (opwd) opwd->cd();
974 //Write statistics information on the client
975 WriteStatisticsMsg(fNcalls);
977 TDirectory *crtdir = gDirectory;
978 TFile f("syswatch.root", "RECREATE");
982 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
983 tree->SetName("syswatch");
984 tree->SetMarkerStyle(kCircle);
985 tree->SetMarkerColor(kBlue);
986 tree->SetMarkerSize(0.5);
987 if (!gROOT->IsBatch()) {
988 tree->SetAlias("event", "id0");
989 tree->SetAlias("task", "id1");
990 tree->SetAlias("stage", "id2");
991 // Already defined aliases
992 // tree->SetAlias("deltaT","stampSec-stampOldSec");
993 // tree->SetAlias("T","stampSec-first");
994 // tree->SetAlias("deltaVM","(pI.fMemVirtual-pIOld.fMemVirtual)");
995 // tree->SetAlias("VM","pI.fMemVirtual");
996 TCanvas *canvas = new TCanvas("SysInfo","SysInfo",10,10,1200,800);
997 Int_t npads = 1 /*COO plot for all tasks*/ +
998 fTopTasks->GetEntries() /*Exec plot per task*/ +
999 1 /*Terminate plot for all tasks*/ +
1002 Int_t iopt = (Int_t)TMath::Sqrt((Double_t)npads);
1003 if (npads<iopt*(iopt+1))
1004 canvas->Divide(iopt, iopt+1, 0.01, 0.01);
1006 canvas->Divide(iopt+1, iopt+1, 0.01, 0.01);
1008 // draw the plot of deltaVM for Exec for each task
1009 for (itask=0; itask<fTopTasks->GetEntriesFast(); itask++) {
1010 task = (AliAnalysisTask*)fTopTasks->At(itask);
1012 cut = Form("task==%d && stage==1", itask);
1013 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1014 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1016 hist->SetTitle(Form("%s: Exec dVM[kB]/event", task->GetName()));
1017 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1020 // Draw the plot of deltaVM for CreateOutputObjects for all tasks
1022 tree->SetMarkerStyle(kFullTriangleUp);
1023 tree->SetMarkerColor(kRed);
1024 tree->SetMarkerSize(0.8);
1025 cut = "task>=0 && task<1000 && stage==0";
1026 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1027 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1029 hist->SetTitle("Memory in CreateOutputObjects()");
1030 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1031 hist->GetXaxis()->SetTitle("task");
1033 // draw the plot of deltaVM for Terminate for all tasks
1035 tree->SetMarkerStyle(kOpenSquare);
1036 tree->SetMarkerColor(kMagenta);
1037 cut = "task>=0 && task<1000 && stage==2";
1038 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1039 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1041 hist->SetTitle("Memory in Terminate()");
1042 hist->GetYaxis()->SetTitle("deltaVM [MB]");
1043 hist->GetXaxis()->SetTitle("task");
1047 tree->SetMarkerStyle(kFullCircle);
1048 tree->SetMarkerColor(kGreen);
1049 cut = Form("task==%d && stage==1",fTopTasks->GetEntriesFast()-1);
1050 tree->Draw("VM:event",cut,"", 1234567890, 0);
1051 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1053 hist->SetTitle("Virtual memory");
1054 hist->GetYaxis()->SetTitle("VM [MB]");
1058 tree->SetMarkerStyle(kCircle);
1059 tree->SetMarkerColor(kBlue);
1060 tree->SetMarkerSize(0.5);
1065 if (crtdir) crtdir->cd();
1067 // Validate the output files
1068 if (ValidateOutputFiles()) {
1070 out.open("outputs_valid", ios::out);
1074 if (fDebug > 1) printf("<-AliAnalysisManager::Terminate()\n");
1076 //______________________________________________________________________________
1077 void AliAnalysisManager::ProfileTask(Int_t itop, const char *option) const
1079 // Profiles the task having the itop index in the list of top (first level) tasks.
1080 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->At(itop);
1082 Error("ProfileTask", "There are only %d top tasks in the manager", fTopTasks->GetEntries());
1085 ProfileTask(task->GetName(), option);
1088 //______________________________________________________________________________
1089 void AliAnalysisManager::ProfileTask(const char *name, const char */*option*/) const
1091 // Profile a managed task after the execution of the analysis in case NSysInfo
1093 if (gSystem->AccessPathName("syswatch.root")) {
1094 Error("ProfileTask", "No file syswatch.root found in the current directory");
1097 if (gROOT->IsBatch()) return;
1098 AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->FindObject(name);
1100 Error("ProfileTask", "No top task named %s known by the manager.", name);
1103 Int_t itop = fTopTasks->IndexOf(task);
1104 Int_t itask = fTasks->IndexOf(task);
1105 // Create canvas with 2 pads: first draw COO + Terminate, second Exec
1106 TDirectory *cdir = gDirectory;
1107 TFile f("syswatch.root");
1108 TTree *tree = (TTree*)f.Get("syswatch");
1110 Error("ProfileTask", "No tree named <syswatch> found in file syswatch.root");
1113 if (fDebug > 1) printf("=== Profiling task %s (class %s)\n", name, task->ClassName());
1114 TCanvas *canvas = new TCanvas(Form("profile_%d",itop),Form("Profile of task %s (class %s)",name,task->ClassName()),10,10,800,600);
1115 canvas->Divide(2, 2, 0.01, 0.01);
1119 // VM profile for COO and Terminate methods
1121 cut = Form("task==%d && (stage==0 || stage==2)",itask);
1122 tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1123 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1125 hist->SetTitle("Alocated VM[kB] for COO and Terminate");
1126 hist->GetYaxis()->SetTitle("deltaVM [kB]");
1127 hist->GetXaxis()->SetTitle("method");
1129 // CPU profile per event
1131 cut = Form("task==%d && stage==1",itop);
1132 tree->Draw("deltaT:event",cut,"", 1234567890, 0);
1133 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1135 hist->SetTitle("Execution time per event");
1136 hist->GetYaxis()->SetTitle("CPU/event [s]");
1138 // VM profile for Exec
1140 cut = Form("task==%d && stage==1",itop);
1141 tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1142 hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1144 hist->SetTitle("Alocated VM[kB] per event");
1145 hist->GetYaxis()->SetTitle("deltaVM [kB]");
1150 if (cdir) cdir->cd();
1153 //______________________________________________________________________________
1154 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1156 // Adds a user task to the global list of tasks.
1157 if (fTasks->FindObject(task)) {
1158 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1161 task->SetActive(kFALSE);
1165 //______________________________________________________________________________
1166 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1168 // Retreive task by name.
1169 if (!fTasks) return NULL;
1170 return (AliAnalysisTask*)fTasks->FindObject(name);
1173 //______________________________________________________________________________
1174 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
1175 TClass *datatype, EAliAnalysisContType type, const char *filename)
1177 // Create a data container of a certain type. Types can be:
1178 // kExchangeContainer = 0, used to exchange data between tasks
1179 // kInputContainer = 1, used to store input data
1180 // kOutputContainer = 2, used for writing result to a file
1181 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1182 // the output object to a folder inside the output file
1183 if (fContainers->FindObject(name)) {
1184 Error("CreateContainer","A container named %s already defined !",name);
1187 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1188 fContainers->Add(cont);
1190 case kInputContainer:
1193 case kOutputContainer:
1194 fOutputs->Add(cont);
1195 if (filename && strlen(filename)) {
1196 cont->SetFileName(filename);
1197 cont->SetDataOwned(kFALSE); // data owned by the file
1200 case kParamContainer:
1201 fParamCont->Add(cont);
1202 if (filename && strlen(filename)) {
1203 cont->SetFileName(filename);
1204 cont->SetDataOwned(kFALSE); // data owned by the file
1207 case kExchangeContainer:
1213 //______________________________________________________________________________
1214 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1215 AliAnalysisDataContainer *cont)
1217 // Connect input of an existing task to a data container.
1219 Error("ConnectInput", "Task pointer is NULL");
1222 if (!fTasks->FindObject(task)) {
1224 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1226 Bool_t connected = task->ConnectInput(islot, cont);
1230 //______________________________________________________________________________
1231 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1232 AliAnalysisDataContainer *cont)
1234 // Connect output of an existing task to a data container.
1236 Error("ConnectOutput", "Task pointer is NULL");
1239 if (!fTasks->FindObject(task)) {
1241 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1243 Bool_t connected = task->ConnectOutput(islot, cont);
1247 //______________________________________________________________________________
1248 void AliAnalysisManager::CleanContainers()
1250 // Clean data from all containers that have already finished all client tasks.
1251 TIter next(fContainers);
1252 AliAnalysisDataContainer *cont;
1253 while ((cont=(AliAnalysisDataContainer *)next())) {
1254 if (cont->IsOwnedData() &&
1255 cont->IsDataReady() &&
1256 cont->ClientsExecuted()) cont->DeleteData();
1260 //______________________________________________________________________________
1261 Bool_t AliAnalysisManager::InitAnalysis()
1263 // Initialization of analysis chain of tasks. Should be called after all tasks
1264 // and data containers are properly connected
1265 // Reset flag and remove valid_outputs file if exists
1267 if (!gSystem->AccessPathName("outputs_valid"))
1268 gSystem->Unlink("outputs_valid");
1269 // Check for top tasks (depending only on input data containers)
1270 if (!fTasks->First()) {
1271 Error("InitAnalysis", "Analysis has no tasks !");
1275 AliAnalysisTask *task;
1276 AliAnalysisDataContainer *cont;
1279 Bool_t iszombie = kFALSE;
1280 Bool_t istop = kTRUE;
1282 while ((task=(AliAnalysisTask*)next())) {
1285 Int_t ninputs = task->GetNinputs();
1286 for (i=0; i<ninputs; i++) {
1287 cont = task->GetInputSlot(i)->GetContainer();
1291 fZombies->Add(task);
1295 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
1296 i, task->GetName());
1298 if (iszombie) continue;
1299 // Check if cont is an input container
1300 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1301 // Connect to parent task
1305 fTopTasks->Add(task);
1309 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1312 // Check now if there are orphan tasks
1313 for (i=0; i<ntop; i++) {
1314 task = (AliAnalysisTask*)fTopTasks->At(i);
1319 while ((task=(AliAnalysisTask*)next())) {
1320 if (!task->IsUsed()) {
1322 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1325 // Check the task hierarchy (no parent task should depend on data provided
1326 // by a daughter task)
1327 for (i=0; i<ntop; i++) {
1328 task = (AliAnalysisTask*)fTopTasks->At(i);
1329 if (task->CheckCircularDeps()) {
1330 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1335 // Check that all containers feeding post-event loop tasks are in the outputs list
1336 TIter nextcont(fContainers); // loop over all containers
1337 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1338 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1339 if (cont->HasConsumers()) {
1340 // Check if one of the consumers is post event loop
1341 TIter nextconsumer(cont->GetConsumers());
1342 while ((task=(AliAnalysisTask*)nextconsumer())) {
1343 if (task->IsPostEventLoop()) {
1344 fOutputs->Add(cont);
1351 // Check if all special output containers have a file name provided
1352 TIter nextout(fOutputs);
1353 while ((cont=(AliAnalysisDataContainer*)nextout())) {
1354 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1355 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1359 // Initialize requested branch list if needed
1360 if (!fAutoBranchHandling) {
1362 while ((task=(AliAnalysisTask*)next())) {
1363 if (!task->HasBranches()) {
1364 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\"",
1365 task->GetName(), task->ClassName());
1368 if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1369 Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1372 TString taskbranches;
1373 task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1374 if (taskbranches.IsNull()) {
1375 Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1376 task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1379 AddBranches(taskbranches);
1386 //______________________________________________________________________________
1387 void AliAnalysisManager::AddBranches(const char *branches)
1389 // Add branches to the existing fRequestedBranches.
1390 TString br(branches);
1391 TObjArray *arr = br.Tokenize(",");
1394 while ((obj=next())) {
1395 if (!fRequestedBranches.Contains(obj->GetName())) {
1396 if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1397 fRequestedBranches += obj->GetName();
1400 if (arr) delete arr;
1403 //______________________________________________________________________________
1404 void AliAnalysisManager::CheckBranches(Bool_t load)
1406 // The method checks the input branches to be loaded during the analysis.
1407 if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;
1408 TObjArray *arr = fRequestedBranches.Tokenize(",");
1411 while ((obj=next())) {
1412 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1414 br = fTree->GetBranch(obj->GetName());
1416 Error("CheckBranches", "Could not find branch %s",obj->GetName());
1421 if (load && br->GetReadEntry()!=GetCurrentEntry()) br->GetEntry(GetCurrentEntry());
1425 //______________________________________________________________________________
1426 void AliAnalysisManager::PrintStatus(Option_t *option) const
1428 // Print task hierarchy.
1430 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1433 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1435 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1436 TIter next(fTopTasks);
1437 AliAnalysisTask *task;
1438 while ((task=(AliAnalysisTask*)next()))
1439 task->PrintTask(option);
1440 if (!fAutoBranchHandling && !fRequestedBranches.IsNull())
1441 printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1444 //______________________________________________________________________________
1445 void AliAnalysisManager::ResetAnalysis()
1447 // Reset all execution flags and clean containers.
1451 //______________________________________________________________________________
1452 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1454 // Start analysis having a grid handler.
1455 if (!fGridHandler) {
1456 Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1457 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1461 return StartAnalysis(type, tree, nentries, firstentry);
1464 //______________________________________________________________________________
1465 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1467 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1468 // MIX. Process nentries starting from firstentry
1470 // Backup current directory and make sure gDirectory points to gROOT
1471 TDirectory *cdir = gDirectory;
1474 Error("StartAnalysis","Analysis manager was not initialized !");
1478 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1480 TString anaType = type;
1482 fMode = kLocalAnalysis;
1483 Bool_t runlocalinit = kTRUE;
1484 if (anaType.Contains("file")) {
1485 runlocalinit = kFALSE;
1488 if (anaType.Contains("proof")) fMode = kProofAnalysis;
1489 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1490 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
1492 if (fMode == kGridAnalysis) {
1494 if (!anaType.Contains("terminate")) {
1495 if (!fGridHandler) {
1496 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1497 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1501 // Write analysis manager in the analysis file
1502 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1503 // run local task configuration
1504 TIter nextTask(fTasks);
1505 AliAnalysisTask *task;
1506 while ((task=(AliAnalysisTask*)nextTask())) {
1510 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1511 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1516 // Terminate grid analysis
1517 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1518 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1519 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1520 if (!fGridHandler->MergeOutputs()) {
1521 // Return if outputs could not be merged or if it alien handler
1522 // was configured for offline mode or local testing.
1527 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1528 ImportWrappers(NULL);
1534 SetEventLoop(kFALSE);
1535 // Enable event loop mode if a tree was provided
1536 if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1539 TString ttype = "TTree";
1540 if (tree && tree->IsA() == TChain::Class()) {
1541 chain = (TChain*)tree;
1542 if (!chain || !chain->GetListOfFiles()->First()) {
1543 Error("StartAnalysis", "Cannot process null or empty chain...");
1550 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1551 if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1552 // Initialize locally all tasks (happens for all modes)
1554 AliAnalysisTask *task;
1556 while ((task=(AliAnalysisTask*)next())) {
1560 if (getsysInfo) AliSysInfo::AddStamp("LocalInit_all", 0);
1564 case kLocalAnalysis:
1565 if (!tree && !fGridHandler) {
1566 TIter nextT(fTasks);
1567 // Call CreateOutputObjects for all tasks
1569 Bool_t dirStatus = TH1::AddDirectoryStatus();
1570 while ((task=(AliAnalysisTask*)nextT())) {
1571 TH1::AddDirectory(kFALSE);
1572 task->CreateOutputObjects();
1573 if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1577 TH1::AddDirectory(dirStatus);
1578 if (IsExternalLoop()) {
1579 Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1580 \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1587 fSelector = new AliAnalysisSelector(this);
1588 // Check if a plugin handler is used
1590 // Get the chain from the plugin
1591 TString dataType = "esdTree";
1592 if (fInputEventHandler) {
1593 dataType = fInputEventHandler->GetDataType();
1597 chain = fGridHandler->GetChainForTestMode(dataType);
1599 Error("StartAnalysis", "No chain for test mode. Aborting.");
1602 cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1603 retv = chain->Process(fSelector, "", nentries, firstentry);
1606 // Run tree-based analysis via AliAnalysisSelector
1607 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1608 retv = tree->Process(fSelector, "", nentries, firstentry);
1610 case kProofAnalysis:
1612 // Check if the plugin is used
1614 return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1616 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1617 Error("StartAnalysis", "No PROOF!!! Exiting.");
1621 line = Form("gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1622 gROOT->ProcessLine(line);
1625 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1626 retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1628 Error("StartAnalysis", "No chain!!! Exiting.");
1635 if (!anaType.Contains("terminate")) {
1636 if (!fGridHandler) {
1637 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1638 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1642 // Write analysis manager in the analysis file
1643 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1644 // Start the analysis via the handler
1645 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1646 Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1651 // Terminate grid analysis
1652 if (fSelector && fSelector->GetStatus() == -1) {cdir->cd(); return -1;}
1653 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {cdir->cd(); return 0;}
1654 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1655 if (!fGridHandler->MergeOutputs()) {
1656 // Return if outputs could not be merged or if it alien handler
1657 // was configured for offline mode or local testing.
1662 cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1663 ImportWrappers(NULL);
1667 case kMixingAnalysis:
1668 // Run event mixing analysis
1670 Error("StartAnalysis", "Cannot run event mixing without event pool");
1674 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1675 fSelector = new AliAnalysisSelector(this);
1676 while ((chain=fEventPool->GetNextChain())) {
1678 // Call NotifyBinChange for all tasks
1679 while ((task=(AliAnalysisTask*)next()))
1680 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1681 retv = chain->Process(fSelector);
1683 Error("StartAnalysis", "Mixing analysis failed");
1688 PackOutput(fSelector->GetOutputList());
1695 //______________________________________________________________________________
1696 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1698 // Start analysis for this manager on a given dataset. Analysis task can be:
1699 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1701 Error("StartAnalysis","Analysis manager was not initialized !");
1705 if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1706 TString anaType = type;
1708 if (!anaType.Contains("proof")) {
1709 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1712 fMode = kProofAnalysis;
1714 SetEventLoop(kTRUE);
1715 // Set the dataset flag
1716 TObject::SetBit(kUseDataSet);
1720 // Start proof analysis using the grid handler
1721 if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1722 Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
1725 // Check if the plugin is in test mode
1726 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
1727 dataset = "test_collection";
1729 dataset = fGridHandler->GetProofDataSet();
1733 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1734 Error("StartAnalysis", "No PROOF!!! Exiting.");
1738 // Initialize locally all tasks
1740 AliAnalysisTask *task;
1741 while ((task=(AliAnalysisTask*)next())) {
1745 line = Form("gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1746 gROOT->ProcessLine(line);
1749 // chain->SetProof();
1750 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON TEST CHAIN " << chain->GetName() << endl;
1751 retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1753 line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1754 dataset, nentries, firstentry);
1755 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1756 retv = (Long_t)gROOT->ProcessLine(line);
1761 //______________________________________________________________________________
1762 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
1764 // Opens according the option the file specified by cont->GetFileName() and changes
1765 // current directory to cont->GetFolderName(). If the file was already opened, it
1766 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
1767 // be optionally ignored.
1768 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1769 TString filename = cont->GetFileName();
1771 if (filename.IsNull()) {
1772 ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
1775 if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
1777 f = mgr->OpenProofFile(cont,option);
1779 // Check first if the file is already opened
1780 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1782 // Check if option "UPDATE" was preserved
1783 TString opt(option);
1785 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1786 ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1788 f = TFile::Open(filename, option);
1791 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1795 // Check for a folder request
1796 TString dir = cont->GetFolderName();
1797 if (!dir.IsNull()) {
1798 if (!f->GetDirectory(dir)) f->mkdir(dir);
1803 ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
1804 cont->SetFile(NULL);
1808 //______________________________________________________________________________
1809 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
1811 // Opens a special output file used in PROOF.
1813 TString filename = cont->GetFileName();
1814 if (cont == fCommonOutput) {
1815 if (fOutputEventHandler) {
1816 if (strlen(extaod)) filename = extaod;
1817 filename = fOutputEventHandler->GetOutputFileName();
1819 else Fatal("OpenProofFile","No output container. Exiting.");
1822 if (fMode!=kProofAnalysis || !fSelector) {
1823 Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
1826 if (fSpecialOutputLocation.Length()) {
1827 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1829 // Check if option "UPDATE" was preserved
1830 TString opt(option);
1832 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1833 ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1835 f = new TFile(filename, option);
1837 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1841 // Check for a folder request
1842 TString dir = cont->GetFolderName();
1844 if (!f->GetDirectory(dir)) f->mkdir(dir);
1849 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1850 cont->SetFile(NULL);
1853 // Check if there is already a proof output file in the output list
1854 TObject *pof = fSelector->GetOutputList()->FindObject(filename);
1856 // Get the actual file
1857 line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
1858 filename = (const char*)gROOT->ProcessLine(line);
1860 printf("File: %s already booked via TProofOutputFile\n", filename.Data());
1862 f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1864 Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
1867 // Check if option "UPDATE" was preserved
1868 TString opt(option);
1870 if ((opt=="UPDATE") && (opt!=f->GetOption()))
1871 Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
1873 if (cont->IsRegisterDataset()) {
1874 TString dsetName = filename;
1875 dsetName.ReplaceAll(".root", cont->GetTitle());
1876 dsetName.ReplaceAll(":","_");
1877 if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
1878 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
1880 if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
1881 line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
1883 if (fDebug > 1) printf("=== %s\n", line.Data());
1884 gROOT->ProcessLine(line);
1885 line = Form("pf->OpenFile(\"%s\");", option);
1886 gROOT->ProcessLine(line);
1889 gROOT->ProcessLine("pf->Print()");
1890 printf(" == proof file name: %s", f->GetName());
1892 // Add to proof output list
1893 line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
1894 if (fDebug > 1) printf("=== %s\n", line.Data());
1895 gROOT->ProcessLine(line);
1897 if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1901 // Check for a folder request
1902 TString dir = cont->GetFolderName();
1903 if (!dir.IsNull()) {
1904 if (!f->GetDirectory(dir)) f->mkdir(dir);
1909 Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
1910 cont->SetFile(NULL);
1914 //______________________________________________________________________________
1915 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1917 // Execute analysis.
1918 static Long64_t nentries = 0;
1919 static TTree *lastTree = 0;
1920 static TStopwatch *timer = new TStopwatch();
1921 if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
1923 if (fTree && (fTree != lastTree)) {
1924 nentries += fTree->GetEntries();
1927 if (!fNcalls) timer->Start();
1928 if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, nentries, timer, kFALSE);
1931 TDirectory *cdir = gDirectory;
1932 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1933 if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
1935 Error("ExecAnalysis", "Analysis manager was not initialized !");
1940 AliAnalysisTask *task;
1941 // Check if the top tree is active.
1943 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1944 AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
1946 // De-activate all tasks
1947 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1948 AliAnalysisDataContainer *cont = fCommonInput;
1949 if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
1951 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1955 cont->SetData(fTree); // This will notify all consumers
1956 Long64_t entry = fTree->GetTree()->GetReadEntry();
1958 // Call BeginEvent() for optional input/output and MC services
1959 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1960 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1961 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1963 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1964 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
1966 // Execute the tasks
1967 // TIter next1(cont->GetConsumers());
1968 TIter next1(fTopTasks);
1970 while ((task=(AliAnalysisTask*)next1())) {
1972 cout << " Executing task " << task->GetName() << endl;
1974 task->ExecuteTask(option);
1976 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1977 AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
1981 // Call FinishEvent() for optional output and MC services
1982 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1983 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1984 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1985 // Gather system information if requested
1986 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1987 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
1991 // The event loop is not controlled by TSelector
1993 // Call BeginEvent() for optional input/output and MC services
1994 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
1995 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
1996 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1998 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
1999 AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2000 TIter next2(fTopTasks);
2001 while ((task=(AliAnalysisTask*)next2())) {
2002 task->SetActive(kTRUE);
2004 cout << " Executing task " << task->GetName() << endl;
2006 task->ExecuteTask(option);
2010 // Call FinishEvent() for optional output and MC services
2011 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
2012 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
2013 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2014 if (getsysInfo && ((fNcalls%fNSysInfo)==0))
2015 AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2019 //______________________________________________________________________________
2020 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2022 // Set the input event handler and create a container for it.
2023 fInputEventHandler = handler;
2024 fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2025 // Warning("SetInputEventHandler", " An automatic input container for the input chain was created.\nPlease use: mgr->GetCommonInputContainer() to access it.");
2028 //______________________________________________________________________________
2029 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2031 // Set the input event handler and create a container for it.
2032 fOutputEventHandler = handler;
2033 fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2034 fCommonOutput->SetSpecialOutput();
2035 // Warning("SetOutputEventHandler", " An automatic output container for the output tree was created.\nPlease use: mgr->GetCommonOutputContainer() to access it.");
2038 //______________________________________________________________________________
2039 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2041 // This method is used externally to register output files which are not
2042 // connected to any output container, so that the manager can properly register,
2043 // retrieve or merge them when running in distributed mode. The file names are
2044 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2045 if (fExtraFiles.Contains(fname)) return;
2046 if (fExtraFiles.Length()) fExtraFiles += " ";
2047 fExtraFiles += fname;
2050 //______________________________________________________________________________
2051 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2053 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2057 TObject *pof = source->FindObject(filename);
2058 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2059 Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2062 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2063 gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2064 TString clientUrl(chUrl);
2065 TString fullPath_str(fullPath);
2066 if (clientUrl.Contains("localhost")){
2067 TObjArray* array = fullPath_str.Tokenize ( "//" );
2068 TObjString *strobj = ( TObjString *)array->At(1);
2069 TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2070 TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2071 fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2072 fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2073 if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2077 else if (clientUrl.Contains("__lite__")) {
2078 // Special case for ProofLite environement - get file info and copy.
2079 gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2080 fullPath_str = Form("%s/%s", tmp, fullPath);
2083 Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2084 Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename);
2086 Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2090 //______________________________________________________________________________
2091 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2093 // Fill analysis type in the provided string.
2095 case kLocalAnalysis:
2098 case kProofAnalysis:
2104 case kMixingAnalysis:
2109 //______________________________________________________________________________
2110 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2112 // Validate all output files.
2113 TIter next(fOutputs);
2114 AliAnalysisDataContainer *output;
2115 TDirectory *cdir = gDirectory;
2116 TString openedFiles;
2117 while ((output=(AliAnalysisDataContainer*)next())) {
2118 if (output->IsRegisterDataset()) continue;
2119 TString filename = output->GetFileName();
2120 if (filename == "default") {
2121 if (!fOutputEventHandler) continue;
2122 filename = fOutputEventHandler->GetOutputFileName();
2123 // Main AOD may not be there
2124 if (gSystem->AccessPathName(filename)) continue;
2126 // Check if the file is closed
2127 if (openedFiles.Contains(filename)) continue;;
2128 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2130 Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2131 // Clear file list to release object ownership to user.
2135 file = TFile::Open(filename);
2136 if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2137 Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2142 openedFiles += filename;
2149 //______________________________________________________________________________
2150 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2152 // Implements a nice text mode progress bar.
2153 static Long64_t icount = 0;
2154 static TString oname;
2155 static TString nname;
2156 static Long64_t ocurrent = 0;
2157 static Long64_t osize = 0;
2158 static Int_t oseconds = 0;
2159 static TStopwatch *owatch = 0;
2160 static Bool_t oneoftwo = kFALSE;
2161 static Int_t nrefresh = 0;
2162 static Int_t nchecks = 0;
2163 const char symbol[4] = {'=','\\','|','/'};
2164 char progress[11] = " ";
2165 Int_t ichar = icount%4;
2172 ocurrent = TMath::Abs(current);
2173 osize = TMath::Abs(size);
2174 if (ocurrent > osize) ocurrent=osize;
2184 if (owatch && !last) {
2186 time = owatch->RealTime();
2187 hours = (Int_t)(time/3600.);
2189 minutes = (Int_t)(time/60.);
2191 seconds = (Int_t)time;
2193 if (oseconds==seconds) {
2197 oneoftwo = !oneoftwo;
2201 if (refresh && oneoftwo) {
2203 if (nchecks <= 0) nchecks = nrefresh+1;
2204 Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2205 oname = Form(" == %d%% ==", pctdone);
2207 Double_t percent = 100.0*ocurrent/osize;
2208 Int_t nchar = Int_t(percent/10);
2209 if (nchar>10) nchar=10;
2211 for (i=0; i<nchar; i++) progress[i] = '=';
2212 progress[nchar] = symbol[ichar];
2213 for (i=nchar+1; i<10; i++) progress[i] = ' ';
2214 progress[10] = '\0';
2217 if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2218 else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2219 else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2220 if (time>0.) fprintf(stderr, "[%6.2f %%] TIME %.2d:%.2d:%.2d \r", percent, hours, minutes, seconds);
2221 else fprintf(stderr, "[%6.2f %%]\r", percent);
2222 if (refresh && oneoftwo) oname = nname;
2223 if (owatch) owatch->Continue();
2232 fprintf(stderr, "\n");
2236 //______________________________________________________________________________
2237 void AliAnalysisManager::DoLoadBranch(const char *name)
2239 // Get tree and load branch if needed.
2240 static Long64_t crtEntry = -100;
2242 if (fAutoBranchHandling || !fTree)
2245 TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2247 br = fTree->GetBranch(name);
2249 Error("DoLoadBranch", "Could not find branch %s",name);
2254 if (br->GetReadEntry()==fCurrentEntry) return;
2255 Int_t ret = br->GetEntry(GetCurrentEntry());
2257 Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2258 if (crtEntry != fCurrentEntry) {
2259 CountEvent(1,0,1,0);
2260 crtEntry = fCurrentEntry;
2263 if (crtEntry != fCurrentEntry) {
2264 CountEvent(1,1,0,0);
2265 crtEntry = fCurrentEntry;
2270 //______________________________________________________________________________
2271 void AliAnalysisManager::AddStatisticsTask()
2273 // Add the statistics task to the manager.
2275 Info("AddStatisticsTask", "Already added");
2278 fStatistics = (AliAnalysisStatistics*)gROOT->ProcessLine("AliAnalysisTaskStat::AddToManager()->GetStatistics();");
2281 //______________________________________________________________________________
2282 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2284 // Bookkeep current event;
2285 if (!fStatistics) return;
2286 fStatistics->AddInput(ninput);
2287 fStatistics->AddProcessed(nprocessed);
2288 fStatistics->AddFailed(nfailed);
2289 fStatistics->AddAccepted(naccepted);
2292 //______________________________________________________________________________
2293 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2295 // Add a line in the statistics message. If available, the statistics message is written
2296 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2298 if (!strlen(line)) return;
2299 if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2300 fStatisticsMsg += line;
2303 //______________________________________________________________________________
2304 void AliAnalysisManager::WriteStatisticsMsg(Int_t nevents)
2306 // Write the statistics message in a file named <nevents.stat>.
2307 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2310 AddStatisticsMsg(Form("Number of input events: %lld",fStatistics->GetNinput()));
2311 AddStatisticsMsg(Form("Number of processed events: %lld",fStatistics->GetNprocessed()));
2312 AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2313 AddStatisticsMsg(Form("Number of accepted events: %lld",fStatistics->GetNaccepted()));
2314 out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2315 fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2316 fStatistics->GetNaccepted()), ios::out);
2317 out << fStatisticsMsg << endl;
2319 if (!nevents) return;
2320 out.open(Form("%09d.stat", nevents), ios::out);
2321 if (!fStatisticsMsg.IsNull()) out << fStatisticsMsg << endl;