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 // AliAnalysysManager - 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 <Riostream.h>
32 #include <TMethodCall.h>
38 #include "AliAnalysisSelector.h"
39 #include "AliAnalysisTask.h"
40 #include "AliAnalysisDataContainer.h"
41 #include "AliAnalysisDataSlot.h"
42 #include "AliVEventHandler.h"
43 #include "AliVEventPool.h"
44 #include "AliSysInfo.h"
45 #include "AliAnalysisManager.h"
47 ClassImp(AliAnalysisManager)
49 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
51 //______________________________________________________________________________
52 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
55 fInputEventHandler(NULL),
56 fOutputEventHandler(NULL),
57 fMCtruthEventHandler(NULL),
61 fMode(kLocalAnalysis),
64 fSpecialOutputLocation(""),
73 // Default constructor.
74 fgAnalysisManager = this;
75 fTasks = new TObjArray();
76 fTopTasks = new TObjArray();
77 fZombies = new TObjArray();
78 fContainers = new TObjArray();
79 fInputs = new TObjArray();
80 fOutputs = new TObjArray();
84 //______________________________________________________________________________
85 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
88 fInputEventHandler(NULL),
89 fOutputEventHandler(NULL),
90 fMCtruthEventHandler(NULL),
95 fInitOK(other.fInitOK),
97 fSpecialOutputLocation(""),
107 fTasks = new TObjArray(*other.fTasks);
108 fTopTasks = new TObjArray(*other.fTopTasks);
109 fZombies = new TObjArray(*other.fZombies);
110 fContainers = new TObjArray(*other.fContainers);
111 fInputs = new TObjArray(*other.fInputs);
112 fOutputs = new TObjArray(*other.fOutputs);
113 fgAnalysisManager = this;
116 //______________________________________________________________________________
117 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
120 if (&other != this) {
121 TNamed::operator=(other);
122 fInputEventHandler = other.fInputEventHandler;
123 fOutputEventHandler = other.fOutputEventHandler;
124 fMCtruthEventHandler = other.fMCtruthEventHandler;
125 fEventPool = other.fEventPool;
128 fNSysInfo = other.fNSysInfo;
130 fInitOK = other.fInitOK;
131 fDebug = other.fDebug;
132 fTasks = new TObjArray(*other.fTasks);
133 fTopTasks = new TObjArray(*other.fTopTasks);
134 fZombies = new TObjArray(*other.fZombies);
135 fContainers = new TObjArray(*other.fContainers);
136 fInputs = new TObjArray(*other.fInputs);
137 fOutputs = new TObjArray(*other.fOutputs);
139 fgAnalysisManager = this;
144 //______________________________________________________________________________
145 AliAnalysisManager::~AliAnalysisManager()
148 if (fTasks) {fTasks->Delete(); delete fTasks;}
149 if (fTopTasks) delete fTopTasks;
150 if (fZombies) delete fZombies;
151 if (fContainers) {fContainers->Delete(); delete fContainers;}
152 if (fInputs) delete fInputs;
153 if (fOutputs) delete fOutputs;
154 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
157 //______________________________________________________________________________
158 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
160 // Read one entry of the tree or a whole branch.
161 if (fDebug > 0) printf("== AliAnalysisManager::GetEntry(%lld)\n", entry);
162 fCurrentEntry = entry;
163 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
166 //______________________________________________________________________________
167 void AliAnalysisManager::Init(TTree *tree)
169 // The Init() function is called when the selector needs to initialize
170 // a new tree or chain. Typically here the branch addresses of the tree
171 // will be set. It is normaly not necessary to make changes to the
172 // generated code, but the routine can be extended by the user if needed.
173 // Init() will be called many times when running with PROOF.
176 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
179 // Call InitTree of EventHandler
180 if (fOutputEventHandler) {
181 if (fMode == kProofAnalysis) {
182 fOutputEventHandler->Init(0x0, "proof");
184 fOutputEventHandler->Init(0x0, "local");
188 if (fInputEventHandler) {
189 if (fMode == kProofAnalysis) {
190 fInputEventHandler->Init(tree, "proof");
192 fInputEventHandler->Init(tree, "local");
195 // If no input event handler we need to get the tree once
197 if(!tree->GetTree()) tree->LoadTree(0);
201 if (fMCtruthEventHandler) {
202 if (fMode == kProofAnalysis) {
203 fMCtruthEventHandler->Init(0x0, "proof");
205 fMCtruthEventHandler->Init(0x0, "local");
209 if (!fInitOK) InitAnalysis();
210 if (!fInitOK) return;
212 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
214 Error("Init","No top input container !");
219 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
223 //______________________________________________________________________________
224 void AliAnalysisManager::SlaveBegin(TTree *tree)
226 // The SlaveBegin() function is called after the Begin() function.
227 // When running with PROOF SlaveBegin() is called on each slave server.
228 // The tree argument is deprecated (on PROOF 0 is passed).
229 if (fDebug > 0) printf("->AliAnalysisManager::SlaveBegin()\n");
231 // Call Init of EventHandler
232 if (fOutputEventHandler) {
233 if (fMode == kProofAnalysis) {
234 fOutputEventHandler->Init("proof");
236 fOutputEventHandler->Init("local");
240 if (fInputEventHandler) {
241 fInputEventHandler->SetInputTree(tree);
242 if (fMode == kProofAnalysis) {
243 fInputEventHandler->Init("proof");
245 fInputEventHandler->Init("local");
249 if (fMCtruthEventHandler) {
250 if (fMode == kProofAnalysis) {
251 fMCtruthEventHandler->Init("proof");
253 fMCtruthEventHandler->Init("local");
258 AliAnalysisTask *task;
259 // Call CreateOutputObjects for all tasks
260 while ((task=(AliAnalysisTask*)next())) {
261 TDirectory *curdir = gDirectory;
262 task->CreateOutputObjects();
263 if (curdir) curdir->cd();
266 if (fDebug > 0) printf("<-AliAnalysisManager::SlaveBegin()\n");
269 //______________________________________________________________________________
270 Bool_t AliAnalysisManager::Notify()
272 // The Notify() function is called when a new file is opened. This
273 // can be either for a new TTree in a TChain or when when a new TTree
274 // is started when using PROOF. It is normaly not necessary to make changes
275 // to the generated code, but the routine can be extended by the
276 // user if needed. The return value is currently not used.
278 Error("Notify","No current tree.");
281 TFile *curfile = fTree->GetCurrentFile();
283 Error("Notify","No current file");
287 if (fDebug > 0) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
289 AliAnalysisTask *task;
290 // Call Notify for all tasks
291 while ((task=(AliAnalysisTask*)next()))
294 // Call Notify of the event handlers
295 if (fInputEventHandler) {
296 fInputEventHandler->Notify(curfile->GetName());
299 if (fOutputEventHandler) {
300 fOutputEventHandler->Notify(curfile->GetName());
303 if (fMCtruthEventHandler) {
304 fMCtruthEventHandler->Notify(curfile->GetName());
306 if (fDebug > 0) printf("<-AliAnalysisManager::Notify()\n");
310 //______________________________________________________________________________
311 Bool_t AliAnalysisManager::Process(Long64_t entry)
313 // The Process() function is called for each entry in the tree (or possibly
314 // keyed object in the case of PROOF) to be processed. The entry argument
315 // specifies which entry in the currently loaded tree is to be processed.
316 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
317 // to read either all or the required parts of the data. When processing
318 // keyed objects with PROOF, the object is already loaded and is available
319 // via the fObject pointer.
321 // This function should contain the "body" of the analysis. It can contain
322 // simple or elaborate selection criteria, run algorithms on the data
323 // of the event and typically fill histograms.
325 // WARNING when a selector is used with a TChain, you must use
326 // the pointer to the current TTree to call GetEntry(entry).
327 // The entry is always the local entry number in the current tree.
328 // Assuming that fChain is the pointer to the TChain being processed,
329 // use fChain->GetTree()->GetEntry(entry).
330 if (fDebug > 0) printf("->AliAnalysisManager::Process(%lld)\n", entry);
332 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
333 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
334 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
338 if (fDebug > 0) printf("<-AliAnalysisManager::Process()\n");
342 //______________________________________________________________________________
343 void AliAnalysisManager::PackOutput(TList *target)
345 // Pack all output data containers in the output list. Called at SlaveTerminate
346 // stage in PROOF case for each slave.
347 if (fDebug > 0) printf("->AliAnalysisManager::PackOutput()\n");
349 Error("PackOutput", "No target. Aborting.");
352 if (fInputEventHandler) fInputEventHandler ->Terminate();
353 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
354 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
356 // Call FinishTaskOutput() for each event loop task (not called for
357 // post-event loop tasks - use Terminate() fo those)
358 TIter nexttask(fTasks);
359 AliAnalysisTask *task;
360 while ((task=(AliAnalysisTask*)nexttask())) {
361 if (!task->IsPostEventLoop()) {
362 if (fDebug > 0) printf("->FinishTaskOutput: task %s\n", task->GetName());
363 task->FinishTaskOutput();
364 if (fDebug > 0) printf("<-FinishTaskOutput: task %s\n", task->GetName());
368 if (fMode == kProofAnalysis) {
369 TIter next(fOutputs);
370 AliAnalysisDataContainer *output;
371 while ((output=(AliAnalysisDataContainer*)next())) {
372 // Do not consider outputs of post event loop tasks
373 if (output->GetProducer()->IsPostEventLoop()) continue;
374 // Check if data was posted to this container. If not, issue an error.
375 if (!output->GetData() ) {
376 Error("PackOutput", "No data for output container %s. Forgot to PostData ?\n", output->GetName());
379 if (!output->IsSpecialOutput()) {
381 const char *filename = output->GetFileName();
382 if (!(strcmp(filename, "default"))) {
383 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
385 if (strlen(filename)) {
386 // File resident outputs
387 TFile *file = output->GetFile();
388 // Backup current folder
389 TDirectory *opwd = gDirectory;
390 // Create file if not existing and register to container.
391 if (file) file->cd();
392 else file = new TFile(filename, "RECREATE");
393 if (file->IsZombie()) {
394 Fatal("PackOutput", "Could not recreate file %s\n", filename);
397 output->SetFile(file);
398 // Clear file list to release object ownership to user.
400 // Save data to file, then close.
401 if (output->GetData()->InheritsFrom(TCollection::Class())) {
402 // If data is a collection, we set the name of the collection
403 // as the one of the container and we save as a single key.
404 TCollection *coll = (TCollection*)output->GetData();
405 coll->SetName(output->GetName());
406 coll->Write(output->GetName(), TObject::kSingleKey);
408 output->GetData()->Write();
410 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
412 printf(" file %s listing content:\n", filename);
416 // Restore current directory
417 if (opwd) opwd->cd();
419 // Memory-resident outputs
420 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", output->GetName());
422 AliAnalysisDataWrapper *wrap = output->ExportData();
423 // Output wrappers must delete data after merging (AG 13/11/07)
424 wrap->SetDeleteData(kTRUE);
428 if (output->IsSpecialOutput()) {
429 TDirectory *opwd = gDirectory;
430 TFile *file = output->GetFile();
432 AliAnalysisTask *producer = output->GetProducer();
434 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
435 output->GetFileName(), output->GetName(), producer->ClassName());
439 // Release object ownership to users after writing data to file
440 if (output->GetData()->InheritsFrom(TCollection::Class())) {
441 // If data is a collection, we set the name of the collection
442 // as the one of the container and we save as a single key.
443 TCollection *coll = (TCollection*)output->GetData();
444 coll->SetName(output->GetName());
445 coll->Write(output->GetName(), TObject::kSingleKey);
447 output->GetData()->Write();
451 printf(" file %s listing content:\n", output->GetFileName());
455 // Restore current directory
456 if (opwd) opwd->cd();
457 // Check if a special output location was provided or the output files have to be merged
458 if (strlen(fSpecialOutputLocation.Data())) {
459 TString remote = fSpecialOutputLocation;
461 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
462 remote += Form("%s_%d_", gSystem->HostName(), gid);
463 remote += output->GetFileName();
464 TFile::Cp(output->GetFileName(), remote.Data());
466 // No special location specified-> use TProofOutputFile as merging utility
467 // The file at this output slot must be opened in CreateOutputObjects
468 if (fDebug > 1) printf(" File %s to be merged...\n", output->GetFileName());
473 if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
476 //______________________________________________________________________________
477 void AliAnalysisManager::ImportWrappers(TList *source)
479 // Import data in output containers from wrappers coming in source.
480 if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
481 TIter next(fOutputs);
482 AliAnalysisDataContainer *cont;
483 AliAnalysisDataWrapper *wrap;
485 while ((cont=(AliAnalysisDataContainer*)next())) {
486 if (cont->GetProducer()->IsPostEventLoop()) continue;
487 if (cont->IsSpecialOutput()) {
488 if (strlen(fSpecialOutputLocation.Data())) continue;
489 // Copy merged file from PROOF scratch space
491 printf(" Copying file %s from PROOF scratch space\n", cont->GetFileName());
492 Bool_t gotit = TFile::Cp(Form("root://lxb6045.cern.ch:11094//pool/scratch/%s",cont->GetFileName()),
493 cont->GetFileName());
495 Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
497 // Normally we should connect data from the copied file to the
498 // corresponding output container, but it is not obvious how to do this
499 // automatically if several objects in file...
502 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
504 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
509 printf(" Importing data for container %s", cont->GetName());
510 if (strlen(cont->GetFileName())) printf(" -> file %s\n", cont->GetFileName());
513 cont->ImportData(wrap);
515 if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
518 //______________________________________________________________________________
519 void AliAnalysisManager::UnpackOutput(TList *source)
521 // Called by AliAnalysisSelector::Terminate only on the client.
522 if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
524 Error("UnpackOutput", "No target. Aborting.");
527 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
529 if (fMode == kProofAnalysis) ImportWrappers(source);
531 TIter next(fOutputs);
532 AliAnalysisDataContainer *output;
533 while ((output=(AliAnalysisDataContainer*)next())) {
534 if (!output->GetData()) continue;
535 // Check if there are client tasks that run post event loop
536 if (output->HasConsumers()) {
537 // Disable event loop semaphore
538 output->SetPostEventLoop(kTRUE);
539 TObjArray *list = output->GetConsumers();
540 Int_t ncons = list->GetEntriesFast();
541 for (Int_t i=0; i<ncons; i++) {
542 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
543 task->CheckNotify(kTRUE);
544 // If task is active, execute it
545 if (task->IsPostEventLoop() && task->IsActive()) {
546 if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
552 if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
555 //______________________________________________________________________________
556 void AliAnalysisManager::Terminate()
558 // The Terminate() function is the last function to be called during
559 // a query. It always runs on the client, it can be used to present
560 // the results graphically.
561 if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
562 AliAnalysisTask *task;
564 // Call Terminate() for tasks
565 while ((task=(AliAnalysisTask*)next())) task->Terminate();
567 TIter next1(fOutputs);
568 AliAnalysisDataContainer *output;
569 while ((output=(AliAnalysisDataContainer*)next1())) {
570 // Close all files at output
571 // Special outputs have the files already closed and written.
572 if (output->IsSpecialOutput()) continue;
573 const char *filename = output->GetFileName();
574 if (!(strcmp(filename, "default"))) {
575 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
576 TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
578 if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
583 if (!strlen(filename)) continue;
584 if (!output->GetData()) continue;
585 TFile *file = output->GetFile();
586 TDirectory *opwd = gDirectory;
590 file = new TFile(filename, "RECREATE");
591 if (file->IsZombie()) continue;
592 output->SetFile(file);
594 if (fDebug > 1) printf(" writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
595 if (output->GetData()->InheritsFrom(TCollection::Class())) {
596 // If data is a collection, we set the name of the collection
597 // as the one of the container and we save as a single key.
598 TCollection *coll = (TCollection*)output->GetData();
599 coll->SetName(output->GetName());
600 coll->Write(output->GetName(), TObject::kSingleKey);
602 output->GetData()->Write();
605 if (opwd) opwd->cd();
608 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
609 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
610 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
612 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
614 TDirectory *cdir = gDirectory;
615 TFile f("syswatch.root", "RECREATE");
617 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
618 tree->SetMarkerStyle(kCircle);
619 tree->SetMarkerColor(kBlue);
620 tree->SetMarkerSize(0.5);
621 if (!gROOT->IsBatch()) {
622 tree->SetAlias("event", "id0");
623 tree->SetAlias("memUSED", "pI.fMemVirtual");
624 tree->SetAlias("userCPU", "pI.fCpuUser");
625 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
626 c->Divide(2,1,0.01,0.01);
628 tree->Draw("memUSED:event","","", 1234567890, 0);
630 tree->Draw("userCPU:event","","", 1234567890, 0);
636 if (cdir) cdir->cd();
638 if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
641 //______________________________________________________________________________
642 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
644 // Adds a user task to the global list of tasks.
645 if (fTasks->FindObject(task)) {
646 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
649 task->SetActive(kFALSE);
653 //______________________________________________________________________________
654 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
656 // Retreive task by name.
657 if (!fTasks) return NULL;
658 return (AliAnalysisTask*)fTasks->FindObject(name);
661 //______________________________________________________________________________
662 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
663 TClass *datatype, EAliAnalysisContType type, const char *filename)
665 // Create a data container of a certain type. Types can be:
666 // kExchangeContainer = 0, used to exchange date between tasks
667 // kInputContainer = 1, used to store input data
668 // kOutputContainer = 2, used for posting results
669 if (fContainers->FindObject(name)) {
670 Error("CreateContainer","A container named %s already defined !\n",name);
673 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
674 fContainers->Add(cont);
676 case kInputContainer:
679 case kOutputContainer:
681 if (filename && strlen(filename)) {
682 cont->SetFileName(filename);
683 cont->SetDataOwned(kFALSE); // data owned by the file
686 case kExchangeContainer:
692 //______________________________________________________________________________
693 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
694 AliAnalysisDataContainer *cont)
696 // Connect input of an existing task to a data container.
697 if (!fTasks->FindObject(task)) {
699 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
701 Bool_t connected = task->ConnectInput(islot, cont);
705 //______________________________________________________________________________
706 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
707 AliAnalysisDataContainer *cont)
709 // Connect output of an existing task to a data container.
710 if (!fTasks->FindObject(task)) {
712 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
714 Bool_t connected = task->ConnectOutput(islot, cont);
718 //______________________________________________________________________________
719 void AliAnalysisManager::CleanContainers()
721 // Clean data from all containers that have already finished all client tasks.
722 TIter next(fContainers);
723 AliAnalysisDataContainer *cont;
724 while ((cont=(AliAnalysisDataContainer *)next())) {
725 if (cont->IsOwnedData() &&
726 cont->IsDataReady() &&
727 cont->ClientsExecuted()) cont->DeleteData();
731 //______________________________________________________________________________
732 Bool_t AliAnalysisManager::InitAnalysis()
734 // Initialization of analysis chain of tasks. Should be called after all tasks
735 // and data containers are properly connected
736 // Check for input/output containers
738 // Check for top tasks (depending only on input data containers)
739 if (!fTasks->First()) {
740 Error("InitAnalysis", "Analysis has no tasks !");
744 AliAnalysisTask *task;
745 AliAnalysisDataContainer *cont;
748 Bool_t iszombie = kFALSE;
749 Bool_t istop = kTRUE;
751 while ((task=(AliAnalysisTask*)next())) {
754 Int_t ninputs = task->GetNinputs();
755 for (i=0; i<ninputs; i++) {
756 cont = task->GetInputSlot(i)->GetContainer();
764 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
767 if (iszombie) continue;
768 // Check if cont is an input container
769 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
770 // Connect to parent task
774 fTopTasks->Add(task);
778 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
781 // Check now if there are orphan tasks
782 for (i=0; i<ntop; i++) {
783 task = (AliAnalysisTask*)fTopTasks->At(i);
788 while ((task=(AliAnalysisTask*)next())) {
789 if (!task->IsUsed()) {
791 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
794 // Check the task hierarchy (no parent task should depend on data provided
795 // by a daughter task)
796 for (i=0; i<ntop; i++) {
797 task = (AliAnalysisTask*)fTopTasks->At(i);
798 if (task->CheckCircularDeps()) {
799 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
804 // Check that all containers feeding post-event loop tasks are in the outputs list
805 TIter nextcont(fContainers); // loop over all containers
806 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
807 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
808 if (cont->HasConsumers()) {
809 // Check if one of the consumers is post event loop
810 TIter nextconsumer(cont->GetConsumers());
811 while ((task=(AliAnalysisTask*)nextconsumer())) {
812 if (task->IsPostEventLoop()) {
820 // Check if all special output containers have a file name provided
821 TIter nextout(fOutputs);
822 while ((cont=(AliAnalysisDataContainer*)nextout())) {
823 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
824 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
832 //______________________________________________________________________________
833 void AliAnalysisManager::PrintStatus(Option_t *option) const
835 // Print task hierarchy.
837 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
840 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
842 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
843 TIter next(fTopTasks);
844 AliAnalysisTask *task;
845 while ((task=(AliAnalysisTask*)next()))
846 task->PrintTask(option);
849 //______________________________________________________________________________
850 void AliAnalysisManager::ResetAnalysis()
852 // Reset all execution flags and clean containers.
856 //______________________________________________________________________________
857 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
859 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
860 // Process nentries starting from firstentry
862 Error("StartAnalysis","Analysis manager was not initialized !");
865 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
866 TString anaType = type;
868 fMode = kLocalAnalysis;
870 if (anaType.Contains("proof")) fMode = kProofAnalysis;
871 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
873 if (fMode == kGridAnalysis) {
874 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
875 fMode = kLocalAnalysis;
878 SetEventLoop(kFALSE);
879 // Enable event loop mode if a tree was provided
880 if (tree) SetEventLoop(kTRUE);
883 TString ttype = "TTree";
884 if (tree->IsA() == TChain::Class()) {
885 chain = (TChain*)tree;
886 if (!chain || !chain->GetListOfFiles()->First()) {
887 Error("StartAnalysis", "Cannot process null or empty chain...");
893 // Initialize locally all tasks
895 AliAnalysisTask *task;
896 while ((task=(AliAnalysisTask*)next())) {
904 // Call CreateOutputObjects for all tasks
905 while ((task=(AliAnalysisTask*)nextT())) {
906 TDirectory *curdir = gDirectory;
907 task->CreateOutputObjects();
908 if (curdir) curdir->cd();
914 // Run tree-based analysis via AliAnalysisSelector
915 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
916 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
917 gROOT->ProcessLine(line);
918 sprintf(line, "((%s*)0x%lx)->Process(selector, \"\",%lld, %lld);",ttype.Data(),(ULong_t)tree, nentries, firstentry);
919 gROOT->ProcessLine(line);
922 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
923 printf("StartAnalysis: no PROOF!!!\n");
926 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
927 gROOT->ProcessLine(line);
930 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
931 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
933 printf("StartAnalysis: no chain\n");
938 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
942 //______________________________________________________________________________
943 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
945 // Start analysis for this manager on a given dataset. Analysis task can be:
946 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
948 Error("StartAnalysis","Analysis manager was not initialized !");
951 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
952 TString anaType = type;
954 if (!anaType.Contains("proof")) {
955 Error("Cannot process datasets in %s mode. Try PROOF.", type);
958 fMode = kProofAnalysis;
961 // Set the dataset flag
962 TObject::SetBit(kUseDataSet);
965 // Initialize locally all tasks
967 AliAnalysisTask *task;
968 while ((task=(AliAnalysisTask*)next())) {
972 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
973 printf("StartAnalysis: no PROOF!!!\n");
976 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
977 gROOT->ProcessLine(line);
978 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
979 if (!gROOT->ProcessLine(line)) {
980 Error("StartAnalysis", "Dataset %s not found", dataset);
983 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
984 dataset, nentries, firstentry);
985 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
986 gROOT->ProcessLine(line);
989 //______________________________________________________________________________
990 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
992 // Opens a special output file used in PROOF.
994 if (fMode!=kProofAnalysis || !fSelector) {
995 Error("OpenProofFile","Cannot open PROOF file %s",filename);
998 sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
999 if (fDebug > 1) printf("=== %s\n", line);
1000 gROOT->ProcessLine(line);
1001 sprintf(line, "pf->OpenFile(\"%s\");", option);
1002 gROOT->ProcessLine(line);
1004 gROOT->ProcessLine("pf->Print()");
1005 printf(" == proof file name: %s\n", gFile->GetName());
1007 sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1008 if (fDebug > 1) printf("=== %s\n", line);
1009 gROOT->ProcessLine(line);
1013 //______________________________________________________________________________
1014 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1016 // Execute analysis.
1017 static Long64_t ncalls = 0;
1018 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1019 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1022 Error("ExecAnalysis", "Analysis manager was not initialized !");
1025 AliAnalysisTask *task;
1026 // Check if the top tree is active.
1029 // De-activate all tasks
1030 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1031 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1033 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1036 cont->SetData(fTree); // This will notify all consumers
1037 Long64_t entry = fTree->GetTree()->GetReadEntry();
1040 // Call BeginEvent() for optional input/output and MC services
1041 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1042 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1043 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1045 // Execute the tasks
1046 // TIter next1(cont->GetConsumers());
1047 TIter next1(fTopTasks);
1048 while ((task=(AliAnalysisTask*)next1())) {
1050 cout << " Executing task " << task->GetName() << endl;
1053 task->ExecuteTask(option);
1056 // Call FinishEvent() for optional output and MC services
1057 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1058 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1059 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1060 // Gather system information if requested
1061 if (getsysInfo && ((ncalls%fNSysInfo)==0))
1062 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1065 // The event loop is not controlled by TSelector
1067 // Call BeginEvent() for optional input/output and MC services
1068 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
1069 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
1070 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1071 TIter next2(fTopTasks);
1072 while ((task=(AliAnalysisTask*)next2())) {
1073 task->SetActive(kTRUE);
1075 cout << " Executing task " << task->GetName() << endl;
1077 task->ExecuteTask(option);
1080 // Call FinishEvent() for optional output and MC services
1081 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1082 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1083 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1086 //______________________________________________________________________________
1087 void AliAnalysisManager::FinishAnalysis()