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
215 if (fMode == kProofAnalysis) {
216 fEventHandler->InitIO("proof");
218 fEventHandler->InitIO("local");
223 AliAnalysisTask *task;
224 // Call CreateOutputObjects for all tasks
225 while ((task=(AliAnalysisTask*)next())) {
226 TDirectory *curdir = gDirectory;
227 task->CreateOutputObjects();
228 if (curdir) curdir->cd();
230 if (fMode == kLocalAnalysis) Init(tree);
232 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
236 //______________________________________________________________________________
237 Bool_t AliAnalysisManager::Notify()
239 // The Notify() function is called when a new file is opened. This
240 // can be either for a new TTree in a TChain or when when a new TTree
241 // is started when using PROOF. It is normaly not necessary to make changes
242 // to the generated code, but the routine can be extended by the
243 // user if needed. The return value is currently not used.
245 TFile *curfile = fTree->GetCurrentFile();
246 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
248 AliAnalysisTask *task;
249 // Call Notify for all tasks
250 while ((task=(AliAnalysisTask*)next()))
256 //______________________________________________________________________________
257 Bool_t AliAnalysisManager::Process(Long64_t entry)
259 // The Process() function is called for each entry in the tree (or possibly
260 // keyed object in the case of PROOF) to be processed. The entry argument
261 // specifies which entry in the currently loaded tree is to be processed.
262 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
263 // to read either all or the required parts of the data. When processing
264 // keyed objects with PROOF, the object is already loaded and is available
265 // via the fObject pointer.
267 // This function should contain the "body" of the analysis. It can contain
268 // simple or elaborate selection criteria, run algorithms on the data
269 // of the event and typically fill histograms.
271 // WARNING when a selector is used with a TChain, you must use
272 // the pointer to the current TTree to call GetEntry(entry).
273 // The entry is always the local entry number in the current tree.
274 // Assuming that fChain is the pointer to the TChain being processed,
275 // use fChain->GetTree()->GetEntry(entry).
277 cout << "->AliAnalysisManager::Process()" << endl;
282 cout << "<-AliAnalysisManager::Process()" << endl;
287 //______________________________________________________________________________
288 void AliAnalysisManager::PackOutput(TList *target)
290 // Pack all output data containers in the output list. Called at SlaveTerminate
291 // stage in PROOF case for each slave.
293 cout << "->AliAnalysisManager::PackOutput()" << endl;
296 Error("PackOutput", "No target. Aborting.");
300 fEventHandler->Terminate();
302 if (fMode == kProofAnalysis) {
303 TIter next(fOutputs);
304 AliAnalysisDataContainer *output;
305 while ((output=(AliAnalysisDataContainer*)next())) {
306 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
307 if (output->GetData()) target->Add(output->ExportData());
311 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
315 //______________________________________________________________________________
316 void AliAnalysisManager::ImportWrappers(TList *source)
318 // Import data in output containers from wrappers coming in source.
320 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
322 TIter next(fOutputs);
323 AliAnalysisDataContainer *cont;
324 AliAnalysisDataWrapper *wrap;
326 while ((cont=(AliAnalysisDataContainer*)next())) {
327 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
328 if (!wrap && fDebug>1) {
329 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
333 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
334 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
335 cont->ImportData(wrap);
338 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
342 //______________________________________________________________________________
343 void AliAnalysisManager::UnpackOutput(TList *source)
345 // Called by AliAnalysisSelector::Terminate. Output containers should
346 // be in source in the same order as in fOutputs.
348 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
351 Error("UnpackOutput", "No target. Aborting.");
355 printf(" Source list contains %d containers\n", source->GetSize());
358 if (fMode == kProofAnalysis) ImportWrappers(source);
360 TIter next(fOutputs);
361 AliAnalysisDataContainer *output;
362 while ((output=(AliAnalysisDataContainer*)next())) {
363 if (!output->GetData()) continue;
364 // Check if there are client tasks that run post event loop
365 if (output->HasConsumers()) {
366 // Disable event loop semaphore
367 output->SetPostEventLoop(kTRUE);
368 TObjArray *list = output->GetConsumers();
369 Int_t ncons = list->GetEntriesFast();
370 for (Int_t i=0; i<ncons; i++) {
371 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
372 task->CheckNotify(kTRUE);
373 // If task is active, execute it
374 if (task->IsPostEventLoop() && task->IsActive()) {
376 cout << "== Executing post event loop task " << task->GetName() << endl;
382 // Check if the output need to be written to a file.
383 const char *filename = output->GetFileName();
384 if (!(strcmp(filename, "default"))) filename = fEventHandler->GetOutputFileName();
386 if (!filename || !strlen(filename)) continue;
387 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
388 if (file) file->cd();
389 else file = new TFile(filename, "RECREATE");
390 if (file->IsZombie()) continue;
391 // Reparent data to this file
393 if (output->GetData()->IsA())
394 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
395 if (callEnv.IsValid()) {
396 callEnv.SetParam((Long_t) file);
397 callEnv.Execute(output->GetData());
399 output->GetData()->Write();
402 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
406 //______________________________________________________________________________
407 void AliAnalysisManager::Terminate()
409 // The Terminate() function is the last function to be called during
410 // a query. It always runs on the client, it can be used to present
411 // the results graphically.
413 cout << "->AliAnalysisManager::Terminate()" << endl;
415 AliAnalysisTask *task;
417 // Call Terminate() for tasks
418 while ((task=(AliAnalysisTask*)next())) task->Terminate();
420 cout << "<-AliAnalysisManager::Terminate()" << endl;
423 fEventHandler->TerminateIO();
426 //______________________________________________________________________________
427 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
429 // Adds a user task to the global list of tasks.
430 task->SetActive(kFALSE);
434 //______________________________________________________________________________
435 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
437 // Retreive task by name.
438 if (!fTasks) return NULL;
439 return (AliAnalysisTask*)fTasks->FindObject(name);
442 //______________________________________________________________________________
443 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
444 TClass *datatype, EAliAnalysisContType type, const char *filename)
446 // Create a data container of a certain type. Types can be:
447 // kExchangeContainer = 0, used to exchange date between tasks
448 // kInputContainer = 1, used to store input data
449 // kOutputContainer = 2, used for posting results
450 if (fContainers->FindObject(name)) {
451 Error("CreateContainer","A container named %s already defined !\n",name);
454 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
455 fContainers->Add(cont);
457 case kInputContainer:
460 case kOutputContainer:
462 if (filename && strlen(filename)) cont->SetFileName(filename);
464 case kExchangeContainer:
470 //______________________________________________________________________________
471 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
472 AliAnalysisDataContainer *cont)
474 // Connect input of an existing task to a data container.
475 if (!fTasks->FindObject(task)) {
477 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
479 Bool_t connected = task->ConnectInput(islot, cont);
483 //______________________________________________________________________________
484 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
485 AliAnalysisDataContainer *cont)
487 // Connect output of an existing task to a data container.
488 if (!fTasks->FindObject(task)) {
490 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
492 Bool_t connected = task->ConnectOutput(islot, cont);
496 //______________________________________________________________________________
497 void AliAnalysisManager::CleanContainers()
499 // Clean data from all containers that have already finished all client tasks.
500 TIter next(fContainers);
501 AliAnalysisDataContainer *cont;
502 while ((cont=(AliAnalysisDataContainer *)next())) {
503 if (cont->IsOwnedData() &&
504 cont->IsDataReady() &&
505 cont->ClientsExecuted()) cont->DeleteData();
509 //______________________________________________________________________________
510 Bool_t AliAnalysisManager::InitAnalysis()
512 // Initialization of analysis chain of tasks. Should be called after all tasks
513 // and data containers are properly connected
514 // Check for input/output containers
516 // Check for top tasks (depending only on input data containers)
517 if (!fTasks->First()) {
518 Error("InitAnalysis", "Analysis has no tasks !");
522 AliAnalysisTask *task;
523 AliAnalysisDataContainer *cont;
526 Bool_t iszombie = kFALSE;
527 Bool_t istop = kTRUE;
529 while ((task=(AliAnalysisTask*)next())) {
532 Int_t ninputs = task->GetNinputs();
533 for (i=0; i<ninputs; i++) {
534 cont = task->GetInputSlot(i)->GetContainer();
542 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
545 if (iszombie) continue;
546 // Check if cont is an input container
547 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
548 // Connect to parent task
552 fTopTasks->Add(task);
556 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
559 // Check now if there are orphan tasks
560 for (i=0; i<ntop; i++) {
561 task = (AliAnalysisTask*)fTopTasks->At(i);
566 while ((task=(AliAnalysisTask*)next())) {
567 if (!task->IsUsed()) {
569 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
572 // Check the task hierarchy (no parent task should depend on data provided
573 // by a daughter task)
574 for (i=0; i<ntop; i++) {
575 task = (AliAnalysisTask*)fTopTasks->At(i);
576 if (task->CheckCircularDeps()) {
577 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
582 // Check that all containers feeding post-event loop tasks are in the outputs list
583 TIter nextcont(fContainers); // loop over all containers
584 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
585 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
586 if (cont->HasConsumers()) {
587 // Check if one of the consumers is post event loop
588 TIter nextconsumer(cont->GetConsumers());
589 while ((task=(AliAnalysisTask*)nextconsumer())) {
590 if (task->IsPostEventLoop()) {
602 //______________________________________________________________________________
603 void AliAnalysisManager::PrintStatus(Option_t *option) const
605 // Print task hierarchy.
606 TIter next(fTopTasks);
607 AliAnalysisTask *task;
608 while ((task=(AliAnalysisTask*)next()))
609 task->PrintTask(option);
612 //______________________________________________________________________________
613 void AliAnalysisManager::ResetAnalysis()
615 // Reset all execution flags and clean containers.
619 //______________________________________________________________________________
620 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
622 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
624 Error("StartAnalysis","Analysis manager was not initialized !");
628 cout << "StartAnalysis: " << GetName() << endl;
630 TString anaType = type;
632 fMode = kLocalAnalysis;
634 if (anaType.Contains("proof")) fMode = kProofAnalysis;
635 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
637 if (fMode == kGridAnalysis) {
638 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
639 fMode = kLocalAnalysis;
642 SetEventLoop(kFALSE);
643 // Disable all branches if requested and set event loop mode
645 if (TestBit(kDisableBranches)) {
646 printf("Disabling all branches...\n");
647 // tree->SetBranchStatus("*",0); // not yet working
652 TChain *chain = dynamic_cast<TChain*>(tree);
654 // Initialize locally all tasks
656 AliAnalysisTask *task;
657 while ((task=(AliAnalysisTask*)next())) {
665 AliAnalysisTask *task;
666 // Call CreateOutputObjects for all tasks
667 while ((task=(AliAnalysisTask*)next())) {
668 TDirectory *curdir = gDirectory;
669 task->CreateOutputObjects();
670 if (curdir) curdir->cd();
676 // Run tree-based analysis via AliAnalysisSelector
677 gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
678 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
679 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
680 gROOT->ProcessLine(line);
681 sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
682 gROOT->ProcessLine(line);
685 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
686 printf("StartAnalysis: no PROOF!!!\n");
689 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
690 gROOT->ProcessLine(line);
693 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
694 chain->Process(gSystem->ExpandPathName("$ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+"));
696 printf("StartAnalysis: no chain\n");
701 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
705 //______________________________________________________________________________
706 void AliAnalysisManager::ExecAnalysis(Option_t *option)
710 Error("ExecAnalysis", "Analysis manager was not initialized !");
713 AliAnalysisTask *task;
714 // Check if the top tree is active.
717 printf("AliAnalysisManager::ExecAnalysis\n");
720 // De-activate all tasks
721 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
722 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
724 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
727 cont->SetData(fTree); // This will notify all consumers
728 TIter next1(cont->GetConsumers());
729 while ((task=(AliAnalysisTask*)next1())) {
730 // task->SetActive(kTRUE);
732 cout << " Executing task " << task->GetName() << endl;
734 task->ExecuteTask(option);
736 fEventHandler->Fill();
739 // The event loop is not controlled by TSelector
740 TIter next2(fTopTasks);
741 while ((task=(AliAnalysisTask*)next2())) {
742 task->SetActive(kTRUE);
744 cout << " Executing task " << task->GetName() << endl;
746 task->ExecuteTask(option);
748 fEventHandler->Fill();
751 //______________________________________________________________________________
752 void AliAnalysisManager::FinishAnalysis()