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>
37 #include "AliAnalysisTask.h"
38 #include "AliAnalysisDataContainer.h"
39 #include "AliAnalysisDataSlot.h"
40 #include "AliVEventHandler.h"
41 #include "AliAnalysisManager.h"
43 ClassImp(AliAnalysisManager)
45 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
47 //______________________________________________________________________________
48 AliAnalysisManager::AliAnalysisManager()
51 fInputEventHandler(NULL),
52 fOutputEventHandler(NULL),
53 fMCtruthEventHandler(NULL),
55 fMode(kLocalAnalysis),
66 fgAnalysisManager = this;
70 //______________________________________________________________________________
71 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
74 fInputEventHandler(NULL),
75 fOutputEventHandler(NULL),
76 fMCtruthEventHandler(NULL),
78 fMode(kLocalAnalysis),
88 // Default constructor.
89 fgAnalysisManager = this;
90 fTasks = new TObjArray();
91 fTopTasks = new TObjArray();
92 fZombies = new TObjArray();
93 fContainers = new TObjArray();
94 fInputs = new TObjArray();
95 fOutputs = new TObjArray();
99 //______________________________________________________________________________
100 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
103 fInputEventHandler(NULL),
104 fOutputEventHandler(NULL),
105 fMCtruthEventHandler(NULL),
108 fInitOK(other.fInitOK),
109 fDebug(other.fDebug),
118 fTasks = new TObjArray(*other.fTasks);
119 fTopTasks = new TObjArray(*other.fTopTasks);
120 fZombies = new TObjArray(*other.fZombies);
121 fContainers = new TObjArray(*other.fContainers);
122 fInputs = new TObjArray(*other.fInputs);
123 fOutputs = new TObjArray(*other.fOutputs);
124 fgAnalysisManager = this;
127 //______________________________________________________________________________
128 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
131 if (&other != this) {
132 TNamed::operator=(other);
133 fInputEventHandler = other.fInputEventHandler;
134 fOutputEventHandler = other.fOutputEventHandler;
135 fMCtruthEventHandler = other.fMCtruthEventHandler;
139 fInitOK = other.fInitOK;
140 fDebug = other.fDebug;
141 fTasks = new TObjArray(*other.fTasks);
142 fTopTasks = new TObjArray(*other.fTopTasks);
143 fZombies = new TObjArray(*other.fZombies);
144 fContainers = new TObjArray(*other.fContainers);
145 fInputs = new TObjArray(*other.fInputs);
146 fOutputs = new TObjArray(*other.fOutputs);
147 fgAnalysisManager = this;
152 //______________________________________________________________________________
153 AliAnalysisManager::~AliAnalysisManager()
156 if (fTasks) {fTasks->Delete(); delete fTasks;}
157 if (fTopTasks) delete fTopTasks;
158 if (fZombies) delete fZombies;
159 if (fContainers) {fContainers->Delete(); delete fContainers;}
160 if (fInputs) delete fInputs;
161 if (fOutputs) delete fOutputs;
162 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
165 //______________________________________________________________________________
166 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
168 // Read one entry of the tree or a whole branch.
170 cout << "== AliAnalysisManager::GetEntry()" << endl;
172 fCurrentEntry = entry;
173 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
176 //______________________________________________________________________________
177 void AliAnalysisManager::Init(TTree *tree)
179 // The Init() function is called when the selector needs to initialize
180 // a new tree or chain. Typically here the branch addresses of the tree
181 // will be set. It is normaly not necessary to make changes to the
182 // generated code, but the routine can be extended by the user if needed.
183 // Init() will be called many times when running with PROOF.
186 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
189 if (fInputEventHandler) {
190 fInputEventHandler->SetInputTree(tree);
191 fInputEventHandler->InitIO("");
194 if (!fInitOK) InitAnalysis();
195 if (!fInitOK) return;
197 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
199 cout<<"Error: No top input container !" <<endl;
204 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
208 //______________________________________________________________________________
209 void AliAnalysisManager::Begin(TTree *tree)
211 // The Begin() function is called at the start of the query.
212 // When running with PROOF Begin() is only called on the client.
213 // The tree argument is deprecated (on PROOF 0 is passed).
215 cout << "AliAnalysisManager::Begin()" << endl;
220 //______________________________________________________________________________
221 void AliAnalysisManager::SlaveBegin(TTree *tree)
223 // The SlaveBegin() function is called after the Begin() function.
224 // When running with PROOF SlaveBegin() is called on each slave server.
225 // The tree argument is deprecated (on PROOF 0 is passed).
227 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
229 // Call InitIO of EventHandler
230 if (fOutputEventHandler) {
231 if (fMode == kProofAnalysis) {
232 fOutputEventHandler->InitIO("proof");
234 fOutputEventHandler->InitIO("local");
237 if (fInputEventHandler && fMode == kLocalAnalysis) {
238 fInputEventHandler->SetInputTree(tree);
239 fInputEventHandler->InitIO("");
244 AliAnalysisTask *task;
245 // Call CreateOutputObjects for all tasks
246 while ((task=(AliAnalysisTask*)next())) {
247 TDirectory *curdir = gDirectory;
248 task->CreateOutputObjects();
249 if (curdir) curdir->cd();
251 if (fMode == kLocalAnalysis)
254 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
258 //______________________________________________________________________________
259 Bool_t AliAnalysisManager::Notify()
261 // The Notify() function is called when a new file is opened. This
262 // can be either for a new TTree in a TChain or when when a new TTree
263 // is started when using PROOF. It is normaly not necessary to make changes
264 // to the generated code, but the routine can be extended by the
265 // user if needed. The return value is currently not used.
267 TFile *curfile = fTree->GetCurrentFile();
268 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
270 AliAnalysisTask *task;
271 // Call Notify for all tasks
272 while ((task=(AliAnalysisTask*)next()))
275 // Call Notify of the MC truth handler
276 if (fMCtruthEventHandler) {
277 fMCtruthEventHandler->Notify(curfile->GetName());
284 //______________________________________________________________________________
285 Bool_t AliAnalysisManager::Process(Long64_t entry)
287 // The Process() function is called for each entry in the tree (or possibly
288 // keyed object in the case of PROOF) to be processed. The entry argument
289 // specifies which entry in the currently loaded tree is to be processed.
290 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
291 // to read either all or the required parts of the data. When processing
292 // keyed objects with PROOF, the object is already loaded and is available
293 // via the fObject pointer.
295 // This function should contain the "body" of the analysis. It can contain
296 // simple or elaborate selection criteria, run algorithms on the data
297 // of the event and typically fill histograms.
299 // WARNING when a selector is used with a TChain, you must use
300 // the pointer to the current TTree to call GetEntry(entry).
301 // The entry is always the local entry number in the current tree.
302 // Assuming that fChain is the pointer to the TChain being processed,
303 // use fChain->GetTree()->GetEntry(entry).
305 cout << "->AliAnalysisManager::Process()" << endl;
307 if (fInputEventHandler) fInputEventHandler ->BeginEvent();
308 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent();
309 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
314 cout << "<-AliAnalysisManager::Process()" << endl;
319 //______________________________________________________________________________
320 void AliAnalysisManager::PackOutput(TList *target)
322 // Pack all output data containers in the output list. Called at SlaveTerminate
323 // stage in PROOF case for each slave.
325 cout << "->AliAnalysisManager::PackOutput()" << endl;
328 Error("PackOutput", "No target. Aborting.");
332 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
333 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
335 if (fMode == kProofAnalysis) {
336 TIter next(fOutputs);
337 AliAnalysisDataContainer *output;
338 while ((output=(AliAnalysisDataContainer*)next())) {
339 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
340 if (output->GetData()) target->Add(output->ExportData());
344 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
348 //______________________________________________________________________________
349 void AliAnalysisManager::ImportWrappers(TList *source)
351 // Import data in output containers from wrappers coming in source.
353 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
355 TIter next(fOutputs);
356 AliAnalysisDataContainer *cont;
357 AliAnalysisDataWrapper *wrap;
359 while ((cont=(AliAnalysisDataContainer*)next())) {
360 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
361 if (!wrap && fDebug>1) {
362 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
366 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
367 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
368 cont->ImportData(wrap);
371 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
375 //______________________________________________________________________________
376 void AliAnalysisManager::UnpackOutput(TList *source)
378 // Called by AliAnalysisSelector::Terminate. Output containers should
379 // be in source in the same order as in fOutputs.
381 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
384 Error("UnpackOutput", "No target. Aborting.");
388 printf(" Source list contains %d containers\n", source->GetSize());
391 if (fMode == kProofAnalysis) ImportWrappers(source);
393 TIter next(fOutputs);
394 AliAnalysisDataContainer *output;
395 while ((output=(AliAnalysisDataContainer*)next())) {
396 if (!output->GetData()) continue;
397 // Check if there are client tasks that run post event loop
398 if (output->HasConsumers()) {
399 // Disable event loop semaphore
400 output->SetPostEventLoop(kTRUE);
401 TObjArray *list = output->GetConsumers();
402 Int_t ncons = list->GetEntriesFast();
403 for (Int_t i=0; i<ncons; i++) {
404 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
405 task->CheckNotify(kTRUE);
406 // If task is active, execute it
407 if (task->IsPostEventLoop() && task->IsActive()) {
409 cout << "== Executing post event loop task " << task->GetName() << endl;
415 // Check if the output need to be written to a file.
416 const char *filename = output->GetFileName();
417 if (!(strcmp(filename, "default"))) {
418 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
421 if (!filename || !strlen(filename)) continue;
422 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
423 if (file) file->cd();
424 else file = new TFile(filename, "RECREATE");
425 if (file->IsZombie()) continue;
426 // Reparent data to this file
428 if (output->GetData()->IsA())
429 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
430 if (callEnv.IsValid()) {
431 callEnv.SetParam((Long_t) file);
432 callEnv.Execute(output->GetData());
434 output->GetData()->Write();
437 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
441 //______________________________________________________________________________
442 void AliAnalysisManager::Terminate()
444 // The Terminate() function is the last function to be called during
445 // a query. It always runs on the client, it can be used to present
446 // the results graphically.
448 cout << "->AliAnalysisManager::Terminate()" << endl;
450 AliAnalysisTask *task;
452 // Call Terminate() for tasks
453 while ((task=(AliAnalysisTask*)next())) task->Terminate();
455 cout << "<-AliAnalysisManager::Terminate()" << endl;
458 if (fOutputEventHandler) fOutputEventHandler->TerminateIO();
461 //______________________________________________________________________________
462 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
464 // Adds a user task to the global list of tasks.
465 task->SetActive(kFALSE);
469 //______________________________________________________________________________
470 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
472 // Retreive task by name.
473 if (!fTasks) return NULL;
474 return (AliAnalysisTask*)fTasks->FindObject(name);
477 //______________________________________________________________________________
478 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
479 TClass *datatype, EAliAnalysisContType type, const char *filename)
481 // Create a data container of a certain type. Types can be:
482 // kExchangeContainer = 0, used to exchange date between tasks
483 // kInputContainer = 1, used to store input data
484 // kOutputContainer = 2, used for posting results
485 if (fContainers->FindObject(name)) {
486 Error("CreateContainer","A container named %s already defined !\n",name);
489 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
490 fContainers->Add(cont);
492 case kInputContainer:
495 case kOutputContainer:
497 if (filename && strlen(filename)) cont->SetFileName(filename);
499 case kExchangeContainer:
505 //______________________________________________________________________________
506 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
507 AliAnalysisDataContainer *cont)
509 // Connect input of an existing task to a data container.
510 if (!fTasks->FindObject(task)) {
512 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
514 Bool_t connected = task->ConnectInput(islot, cont);
518 //______________________________________________________________________________
519 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
520 AliAnalysisDataContainer *cont)
522 // Connect output of an existing task to a data container.
523 if (!fTasks->FindObject(task)) {
525 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
527 Bool_t connected = task->ConnectOutput(islot, cont);
531 //______________________________________________________________________________
532 void AliAnalysisManager::CleanContainers()
534 // Clean data from all containers that have already finished all client tasks.
535 TIter next(fContainers);
536 AliAnalysisDataContainer *cont;
537 while ((cont=(AliAnalysisDataContainer *)next())) {
538 if (cont->IsOwnedData() &&
539 cont->IsDataReady() &&
540 cont->ClientsExecuted()) cont->DeleteData();
544 //______________________________________________________________________________
545 Bool_t AliAnalysisManager::InitAnalysis()
547 // Initialization of analysis chain of tasks. Should be called after all tasks
548 // and data containers are properly connected
549 // Check for input/output containers
551 // Check for top tasks (depending only on input data containers)
552 if (!fTasks->First()) {
553 Error("InitAnalysis", "Analysis has no tasks !");
557 AliAnalysisTask *task;
558 AliAnalysisDataContainer *cont;
561 Bool_t iszombie = kFALSE;
562 Bool_t istop = kTRUE;
564 while ((task=(AliAnalysisTask*)next())) {
567 Int_t ninputs = task->GetNinputs();
568 for (i=0; i<ninputs; i++) {
569 cont = task->GetInputSlot(i)->GetContainer();
577 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
580 if (iszombie) continue;
581 // Check if cont is an input container
582 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
583 // Connect to parent task
587 fTopTasks->Add(task);
591 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
594 // Check now if there are orphan tasks
595 for (i=0; i<ntop; i++) {
596 task = (AliAnalysisTask*)fTopTasks->At(i);
601 while ((task=(AliAnalysisTask*)next())) {
602 if (!task->IsUsed()) {
604 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
607 // Check the task hierarchy (no parent task should depend on data provided
608 // by a daughter task)
609 for (i=0; i<ntop; i++) {
610 task = (AliAnalysisTask*)fTopTasks->At(i);
611 if (task->CheckCircularDeps()) {
612 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
617 // Check that all containers feeding post-event loop tasks are in the outputs list
618 TIter nextcont(fContainers); // loop over all containers
619 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
620 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
621 if (cont->HasConsumers()) {
622 // Check if one of the consumers is post event loop
623 TIter nextconsumer(cont->GetConsumers());
624 while ((task=(AliAnalysisTask*)nextconsumer())) {
625 if (task->IsPostEventLoop()) {
637 //______________________________________________________________________________
638 void AliAnalysisManager::PrintStatus(Option_t *option) const
640 // Print task hierarchy.
641 TIter next(fTopTasks);
642 AliAnalysisTask *task;
643 while ((task=(AliAnalysisTask*)next()))
644 task->PrintTask(option);
647 //______________________________________________________________________________
648 void AliAnalysisManager::ResetAnalysis()
650 // Reset all execution flags and clean containers.
654 //______________________________________________________________________________
655 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
657 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
659 Error("StartAnalysis","Analysis manager was not initialized !");
663 cout << "StartAnalysis: " << GetName() << endl;
665 TString anaType = type;
667 fMode = kLocalAnalysis;
669 if (anaType.Contains("proof")) fMode = kProofAnalysis;
670 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
672 if (fMode == kGridAnalysis) {
673 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
674 fMode = kLocalAnalysis;
677 SetEventLoop(kFALSE);
678 // Disable all branches if requested and set event loop mode
680 if (TestBit(kDisableBranches)) {
681 printf("Disabling all branches...\n");
682 // tree->SetBranchStatus("*",0); // not yet working
687 TChain *chain = dynamic_cast<TChain*>(tree);
689 // Initialize locally all tasks
691 AliAnalysisTask *task;
692 while ((task=(AliAnalysisTask*)next())) {
700 AliAnalysisTask *task;
701 // Call CreateOutputObjects for all tasks
702 while ((task=(AliAnalysisTask*)next())) {
703 TDirectory *curdir = gDirectory;
704 task->CreateOutputObjects();
705 if (curdir) curdir->cd();
711 // Run tree-based analysis via AliAnalysisSelector
712 // gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
713 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
714 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
715 gROOT->ProcessLine(line);
716 sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
717 gROOT->ProcessLine(line);
720 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
721 printf("StartAnalysis: no PROOF!!!\n");
724 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
725 gROOT->ProcessLine(line);
728 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
729 chain->Process("AliAnalysisSelector");
731 printf("StartAnalysis: no chain\n");
736 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
740 //______________________________________________________________________________
741 void AliAnalysisManager::ExecAnalysis(Option_t *option)
745 Error("ExecAnalysis", "Analysis manager was not initialized !");
748 AliAnalysisTask *task;
749 // Check if the top tree is active.
752 printf("AliAnalysisManager::ExecAnalysis\n");
755 // De-activate all tasks
756 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
757 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
759 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
762 cont->SetData(fTree); // This will notify all consumers
764 // Call BeginEvent() for optional input/output and MC services
765 if (fInputEventHandler) fInputEventHandler ->BeginEvent();
766 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent();
767 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
770 TIter next1(cont->GetConsumers());
771 while ((task=(AliAnalysisTask*)next1())) {
773 cout << " Executing task " << task->GetName() << endl;
776 task->ExecuteTask(option);
779 // Call FinishEvent() for optional output and MC services
780 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
781 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
785 // The event loop is not controlled by TSelector
787 // Call BeginEvent() for optional input/output and MC services
788 if (fInputEventHandler) fInputEventHandler ->BeginEvent();
789 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent();
790 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
791 TIter next2(fTopTasks);
792 while ((task=(AliAnalysisTask*)next2())) {
793 task->SetActive(kTRUE);
795 cout << " Executing task " << task->GetName() << endl;
797 task->ExecuteTask(option);
800 // Call FinishEvent() for optional output and MC services
801 if (fOutputEventHandler) fOutputEventHandler->FinishEvent();
802 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
805 //______________________________________________________________________________
806 void AliAnalysisManager::FinishAnalysis()