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 "AliVirtualEventHandler.h"
41 #include "AliAnalysisManager.h"
43 ClassImp(AliAnalysisManager)
45 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
47 //______________________________________________________________________________
48 AliAnalysisManager::AliAnalysisManager()
53 fMode(kLocalAnalysis),
63 fgAnalysisManager = this;
67 //______________________________________________________________________________
68 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
73 fMode(kLocalAnalysis),
83 // Default constructor.
84 fgAnalysisManager = this;
85 fTasks = new TObjArray();
86 fTopTasks = new TObjArray();
87 fZombies = new TObjArray();
88 fContainers = new TObjArray();
89 fInputs = new TObjArray();
90 fOutputs = new TObjArray();
94 //______________________________________________________________________________
95 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
101 fInitOK(other.fInitOK),
102 fDebug(other.fDebug),
111 fTasks = new TObjArray(*other.fTasks);
112 fTopTasks = new TObjArray(*other.fTopTasks);
113 fZombies = new TObjArray(*other.fZombies);
114 fContainers = new TObjArray(*other.fContainers);
115 fInputs = new TObjArray(*other.fInputs);
116 fOutputs = new TObjArray(*other.fOutputs);
117 fgAnalysisManager = this;
120 //______________________________________________________________________________
121 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
124 if (&other != this) {
125 TNamed::operator=(other);
127 fEventHandler = other.fEventHandler;
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);
138 fgAnalysisManager = this;
143 //______________________________________________________________________________
144 AliAnalysisManager::~AliAnalysisManager()
147 if (fTasks) {fTasks->Delete(); delete fTasks;}
148 if (fTopTasks) delete fTopTasks;
149 if (fZombies) delete fZombies;
150 if (fContainers) {fContainers->Delete(); delete fContainers;}
151 if (fInputs) delete fInputs;
152 if (fOutputs) delete fOutputs;
153 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
156 //______________________________________________________________________________
157 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
159 // Read one entry of the tree or a whole branch.
161 cout << "== AliAnalysisManager::GetEntry()" << endl;
163 fCurrentEntry = entry;
164 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
167 //______________________________________________________________________________
168 void AliAnalysisManager::Init(TTree *tree)
170 // The Init() function is called when the selector needs to initialize
171 // a new tree or chain. Typically here the branch addresses of the tree
172 // will be set. It is normaly not necessary to make changes to the
173 // generated code, but the routine can be extended by the user if needed.
174 // Init() will be called many times when running with PROOF.
177 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
179 if (!fInitOK) InitAnalysis();
180 if (!fInitOK) return;
182 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
184 cout<<"Error: No top input container !" <<endl;
189 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
193 //______________________________________________________________________________
194 void AliAnalysisManager::Begin(TTree *tree)
196 // The Begin() function is called at the start of the query.
197 // When running with PROOF Begin() is only called on the client.
198 // The tree argument is deprecated (on PROOF 0 is passed).
200 cout << "AliAnalysisManager::Begin()" << endl;
205 //______________________________________________________________________________
206 void AliAnalysisManager::SlaveBegin(TTree *tree)
208 // The SlaveBegin() function is called after the Begin() function.
209 // When running with PROOF SlaveBegin() is called on each slave server.
210 // The tree argument is deprecated (on PROOF 0 is passed).
212 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
214 // Call InitIO of EventHandler
216 if (fMode == kProofAnalysis) {
217 fEventHandler->InitIO("proof");
219 fEventHandler->InitIO("local");
226 AliAnalysisTask *task;
227 // Call CreateOutputObjects for all tasks
228 while ((task=(AliAnalysisTask*)next())) {
229 TDirectory *curdir = gDirectory;
230 task->CreateOutputObjects();
231 if (curdir) curdir->cd();
233 if (fMode == kLocalAnalysis) Init(tree);
235 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
239 //______________________________________________________________________________
240 Bool_t AliAnalysisManager::Notify()
242 // The Notify() function is called when a new file is opened. This
243 // can be either for a new TTree in a TChain or when when a new TTree
244 // is started when using PROOF. It is normaly not necessary to make changes
245 // to the generated code, but the routine can be extended by the
246 // user if needed. The return value is currently not used.
248 TFile *curfile = fTree->GetCurrentFile();
249 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
251 AliAnalysisTask *task;
252 // Call Notify for all tasks
253 while ((task=(AliAnalysisTask*)next()))
259 //______________________________________________________________________________
260 Bool_t AliAnalysisManager::Process(Long64_t entry)
262 // The Process() function is called for each entry in the tree (or possibly
263 // keyed object in the case of PROOF) to be processed. The entry argument
264 // specifies which entry in the currently loaded tree is to be processed.
265 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
266 // to read either all or the required parts of the data. When processing
267 // keyed objects with PROOF, the object is already loaded and is available
268 // via the fObject pointer.
270 // This function should contain the "body" of the analysis. It can contain
271 // simple or elaborate selection criteria, run algorithms on the data
272 // of the event and typically fill histograms.
274 // WARNING when a selector is used with a TChain, you must use
275 // the pointer to the current TTree to call GetEntry(entry).
276 // The entry is always the local entry number in the current tree.
277 // Assuming that fChain is the pointer to the TChain being processed,
278 // use fChain->GetTree()->GetEntry(entry).
280 cout << "->AliAnalysisManager::Process()" << endl;
285 cout << "<-AliAnalysisManager::Process()" << endl;
290 //______________________________________________________________________________
291 void AliAnalysisManager::PackOutput(TList *target)
293 // Pack all output data containers in the output list. Called at SlaveTerminate
294 // stage in PROOF case for each slave.
296 cout << "->AliAnalysisManager::PackOutput()" << endl;
299 Error("PackOutput", "No target. Aborting.");
303 if (fEventHandler) fEventHandler->Terminate();
305 if (fMode == kProofAnalysis) {
306 TIter next(fOutputs);
307 AliAnalysisDataContainer *output;
308 while ((output=(AliAnalysisDataContainer*)next())) {
309 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
310 if (output->GetData()) target->Add(output->ExportData());
314 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
318 //______________________________________________________________________________
319 void AliAnalysisManager::ImportWrappers(TList *source)
321 // Import data in output containers from wrappers coming in source.
323 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
325 TIter next(fOutputs);
326 AliAnalysisDataContainer *cont;
327 AliAnalysisDataWrapper *wrap;
329 while ((cont=(AliAnalysisDataContainer*)next())) {
330 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
331 if (!wrap && fDebug>1) {
332 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
336 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
337 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
338 cont->ImportData(wrap);
341 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
345 //______________________________________________________________________________
346 void AliAnalysisManager::UnpackOutput(TList *source)
348 // Called by AliAnalysisSelector::Terminate. Output containers should
349 // be in source in the same order as in fOutputs.
351 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
354 Error("UnpackOutput", "No target. Aborting.");
358 printf(" Source list contains %d containers\n", source->GetSize());
361 if (fMode == kProofAnalysis) ImportWrappers(source);
363 TIter next(fOutputs);
364 AliAnalysisDataContainer *output;
365 while ((output=(AliAnalysisDataContainer*)next())) {
366 if (!output->GetData()) continue;
367 // Check if there are client tasks that run post event loop
368 if (output->HasConsumers()) {
369 // Disable event loop semaphore
370 output->SetPostEventLoop(kTRUE);
371 TObjArray *list = output->GetConsumers();
372 Int_t ncons = list->GetEntriesFast();
373 for (Int_t i=0; i<ncons; i++) {
374 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
375 task->CheckNotify(kTRUE);
376 // If task is active, execute it
377 if (task->IsPostEventLoop() && task->IsActive()) {
379 cout << "== Executing post event loop task " << task->GetName() << endl;
385 // Check if the output need to be written to a file.
386 const char *filename = output->GetFileName();
387 if (!(strcmp(filename, "default"))) {
388 if (fEventHandler) filename = fEventHandler->GetOutputFileName();
391 if (!filename || !strlen(filename)) continue;
392 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
393 if (file) file->cd();
394 else file = new TFile(filename, "RECREATE");
395 if (file->IsZombie()) continue;
396 // Reparent data to this file
398 if (output->GetData()->IsA())
399 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
400 if (callEnv.IsValid()) {
401 callEnv.SetParam((Long_t) file);
402 callEnv.Execute(output->GetData());
404 output->GetData()->Write();
407 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
411 //______________________________________________________________________________
412 void AliAnalysisManager::Terminate()
414 // The Terminate() function is the last function to be called during
415 // a query. It always runs on the client, it can be used to present
416 // the results graphically.
418 cout << "->AliAnalysisManager::Terminate()" << endl;
420 AliAnalysisTask *task;
422 // Call Terminate() for tasks
423 while ((task=(AliAnalysisTask*)next())) task->Terminate();
425 cout << "<-AliAnalysisManager::Terminate()" << endl;
428 if (fEventHandler) fEventHandler->TerminateIO();
431 //______________________________________________________________________________
432 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
434 // Adds a user task to the global list of tasks.
435 task->SetActive(kFALSE);
439 //______________________________________________________________________________
440 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
442 // Retreive task by name.
443 if (!fTasks) return NULL;
444 return (AliAnalysisTask*)fTasks->FindObject(name);
447 //______________________________________________________________________________
448 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
449 TClass *datatype, EAliAnalysisContType type, const char *filename)
451 // Create a data container of a certain type. Types can be:
452 // kExchangeContainer = 0, used to exchange date between tasks
453 // kInputContainer = 1, used to store input data
454 // kOutputContainer = 2, used for posting results
455 if (fContainers->FindObject(name)) {
456 Error("CreateContainer","A container named %s already defined !\n",name);
459 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
460 fContainers->Add(cont);
462 case kInputContainer:
465 case kOutputContainer:
467 if (filename && strlen(filename)) cont->SetFileName(filename);
469 case kExchangeContainer:
475 //______________________________________________________________________________
476 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
477 AliAnalysisDataContainer *cont)
479 // Connect input of an existing task to a data container.
480 if (!fTasks->FindObject(task)) {
482 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
484 Bool_t connected = task->ConnectInput(islot, cont);
488 //______________________________________________________________________________
489 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
490 AliAnalysisDataContainer *cont)
492 // Connect output of an existing task to a data container.
493 if (!fTasks->FindObject(task)) {
495 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
497 Bool_t connected = task->ConnectOutput(islot, cont);
501 //______________________________________________________________________________
502 void AliAnalysisManager::CleanContainers()
504 // Clean data from all containers that have already finished all client tasks.
505 TIter next(fContainers);
506 AliAnalysisDataContainer *cont;
507 while ((cont=(AliAnalysisDataContainer *)next())) {
508 if (cont->IsOwnedData() &&
509 cont->IsDataReady() &&
510 cont->ClientsExecuted()) cont->DeleteData();
514 //______________________________________________________________________________
515 Bool_t AliAnalysisManager::InitAnalysis()
517 // Initialization of analysis chain of tasks. Should be called after all tasks
518 // and data containers are properly connected
519 // Check for input/output containers
521 // Check for top tasks (depending only on input data containers)
522 if (!fTasks->First()) {
523 Error("InitAnalysis", "Analysis has no tasks !");
527 AliAnalysisTask *task;
528 AliAnalysisDataContainer *cont;
531 Bool_t iszombie = kFALSE;
532 Bool_t istop = kTRUE;
534 while ((task=(AliAnalysisTask*)next())) {
537 Int_t ninputs = task->GetNinputs();
538 for (i=0; i<ninputs; i++) {
539 cont = task->GetInputSlot(i)->GetContainer();
547 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
550 if (iszombie) continue;
551 // Check if cont is an input container
552 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
553 // Connect to parent task
557 fTopTasks->Add(task);
561 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
564 // Check now if there are orphan tasks
565 for (i=0; i<ntop; i++) {
566 task = (AliAnalysisTask*)fTopTasks->At(i);
571 while ((task=(AliAnalysisTask*)next())) {
572 if (!task->IsUsed()) {
574 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
577 // Check the task hierarchy (no parent task should depend on data provided
578 // by a daughter task)
579 for (i=0; i<ntop; i++) {
580 task = (AliAnalysisTask*)fTopTasks->At(i);
581 if (task->CheckCircularDeps()) {
582 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
587 // Check that all containers feeding post-event loop tasks are in the outputs list
588 TIter nextcont(fContainers); // loop over all containers
589 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
590 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
591 if (cont->HasConsumers()) {
592 // Check if one of the consumers is post event loop
593 TIter nextconsumer(cont->GetConsumers());
594 while ((task=(AliAnalysisTask*)nextconsumer())) {
595 if (task->IsPostEventLoop()) {
607 //______________________________________________________________________________
608 void AliAnalysisManager::PrintStatus(Option_t *option) const
610 // Print task hierarchy.
611 TIter next(fTopTasks);
612 AliAnalysisTask *task;
613 while ((task=(AliAnalysisTask*)next()))
614 task->PrintTask(option);
617 //______________________________________________________________________________
618 void AliAnalysisManager::ResetAnalysis()
620 // Reset all execution flags and clean containers.
624 //______________________________________________________________________________
625 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
627 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
629 Error("StartAnalysis","Analysis manager was not initialized !");
633 cout << "StartAnalysis: " << GetName() << endl;
635 TString anaType = type;
637 fMode = kLocalAnalysis;
639 if (anaType.Contains("proof")) fMode = kProofAnalysis;
640 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
642 if (fMode == kGridAnalysis) {
643 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
644 fMode = kLocalAnalysis;
647 SetEventLoop(kFALSE);
648 // Disable all branches if requested and set event loop mode
650 if (TestBit(kDisableBranches)) {
651 printf("Disabling all branches...\n");
652 // tree->SetBranchStatus("*",0); // not yet working
657 TChain *chain = dynamic_cast<TChain*>(tree);
659 // Initialize locally all tasks
661 AliAnalysisTask *task;
662 while ((task=(AliAnalysisTask*)next())) {
670 AliAnalysisTask *task;
671 // Call CreateOutputObjects for all tasks
672 while ((task=(AliAnalysisTask*)next())) {
673 TDirectory *curdir = gDirectory;
674 task->CreateOutputObjects();
675 if (curdir) curdir->cd();
681 // Run tree-based analysis via AliAnalysisSelector
682 gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
683 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
684 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
685 gROOT->ProcessLine(line);
686 sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
687 gROOT->ProcessLine(line);
690 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
691 printf("StartAnalysis: no PROOF!!!\n");
694 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
695 gROOT->ProcessLine(line);
698 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
699 chain->Process(gSystem->ExpandPathName("$ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+"));
701 printf("StartAnalysis: no chain\n");
706 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
710 //______________________________________________________________________________
711 void AliAnalysisManager::ExecAnalysis(Option_t *option)
715 Error("ExecAnalysis", "Analysis manager was not initialized !");
718 AliAnalysisTask *task;
719 // Check if the top tree is active.
722 printf("AliAnalysisManager::ExecAnalysis\n");
725 // De-activate all tasks
726 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
727 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
729 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
732 cont->SetData(fTree); // This will notify all consumers
733 TIter next1(cont->GetConsumers());
734 while ((task=(AliAnalysisTask*)next1())) {
735 // task->SetActive(kTRUE);
737 cout << " Executing task " << task->GetName() << endl;
739 task->ExecuteTask(option);
741 if (fEventHandler) fEventHandler->Fill();
744 // The event loop is not controlled by TSelector
745 TIter next2(fTopTasks);
746 while ((task=(AliAnalysisTask*)next2())) {
747 task->SetActive(kTRUE);
749 cout << " Executing task " << task->GetName() << endl;
751 task->ExecuteTask(option);
753 if (fEventHandler) fEventHandler->Fill();
756 //______________________________________________________________________________
757 void AliAnalysisManager::FinishAnalysis()