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 "AliAnalysisManager.h"
42 ClassImp(AliAnalysisManager)
44 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
46 //______________________________________________________________________________
47 AliAnalysisManager::AliAnalysisManager()
51 fMode(kLocalAnalysis),
61 fgAnalysisManager = this;
65 //______________________________________________________________________________
66 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
70 fMode(kLocalAnalysis),
80 // Default constructor.
81 fgAnalysisManager = this;
82 fTasks = new TObjArray();
83 fTopTasks = new TObjArray();
84 fZombies = new TObjArray();
85 fContainers = new TObjArray();
86 fInputs = new TObjArray();
87 fOutputs = new TObjArray();
91 //______________________________________________________________________________
92 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
97 fInitOK(other.fInitOK),
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);
125 fInitOK = other.fInitOK;
126 fDebug = other.fDebug;
127 fTasks = new TObjArray(*other.fTasks);
128 fTopTasks = new TObjArray(*other.fTopTasks);
129 fZombies = new TObjArray(*other.fZombies);
130 fContainers = new TObjArray(*other.fContainers);
131 fInputs = new TObjArray(*other.fInputs);
132 fOutputs = new TObjArray(*other.fOutputs);
133 fgAnalysisManager = this;
138 //______________________________________________________________________________
139 AliAnalysisManager::~AliAnalysisManager()
142 if (fTasks) {fTasks->Delete(); delete fTasks;}
143 if (fTopTasks) delete fTopTasks;
144 if (fZombies) delete fZombies;
145 if (fContainers) {fContainers->Delete(); delete fContainers;}
146 if (fInputs) delete fInputs;
147 if (fOutputs) delete fOutputs;
148 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
151 //______________________________________________________________________________
152 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
154 // Read one entry of the tree or a whole branch.
156 cout << "== AliAnalysisManager::GetEntry()" << endl;
158 fCurrentEntry = entry;
159 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
162 //______________________________________________________________________________
163 void AliAnalysisManager::Init(TTree *tree)
165 // The Init() function is called when the selector needs to initialize
166 // a new tree or chain. Typically here the branch addresses of the tree
167 // will be set. It is normaly not necessary to make changes to the
168 // generated code, but the routine can be extended by the user if needed.
169 // Init() will be called many times when running with PROOF.
172 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
174 if (!fInitOK) InitAnalysis();
175 if (!fInitOK) return;
177 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
179 cout<<"Error: No top input container !" <<endl;
184 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
188 //______________________________________________________________________________
189 void AliAnalysisManager::Begin(TTree *tree)
191 // The Begin() function is called at the start of the query.
192 // When running with PROOF Begin() is only called on the client.
193 // The tree argument is deprecated (on PROOF 0 is passed).
195 cout << "AliAnalysisManager::Begin()" << endl;
200 //______________________________________________________________________________
201 void AliAnalysisManager::SlaveBegin(TTree *tree)
203 // The SlaveBegin() function is called after the Begin() function.
204 // When running with PROOF SlaveBegin() is called on each slave server.
205 // The tree argument is deprecated (on PROOF 0 is passed).
207 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
211 AliAnalysisTask *task;
212 // Call CreateOutputObjects for all tasks
213 while ((task=(AliAnalysisTask*)next())) {
214 TDirectory *curdir = gDirectory;
215 task->CreateOutputObjects();
216 if (curdir) curdir->cd();
218 if (fMode == kLocalAnalysis) Init(tree);
220 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
224 //______________________________________________________________________________
225 Bool_t AliAnalysisManager::Notify()
227 // The Notify() function is called when a new file is opened. This
228 // can be either for a new TTree in a TChain or when when a new TTree
229 // is started when using PROOF. It is normaly not necessary to make changes
230 // to the generated code, but the routine can be extended by the
231 // user if needed. The return value is currently not used.
233 TFile *curfile = fTree->GetCurrentFile();
234 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
236 AliAnalysisTask *task;
237 // Call Notify for all tasks
238 while ((task=(AliAnalysisTask*)next()))
244 //______________________________________________________________________________
245 Bool_t AliAnalysisManager::Process(Long64_t entry)
247 // The Process() function is called for each entry in the tree (or possibly
248 // keyed object in the case of PROOF) to be processed. The entry argument
249 // specifies which entry in the currently loaded tree is to be processed.
250 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
251 // to read either all or the required parts of the data. When processing
252 // keyed objects with PROOF, the object is already loaded and is available
253 // via the fObject pointer.
255 // This function should contain the "body" of the analysis. It can contain
256 // simple or elaborate selection criteria, run algorithms on the data
257 // of the event and typically fill histograms.
259 // WARNING when a selector is used with a TChain, you must use
260 // the pointer to the current TTree to call GetEntry(entry).
261 // The entry is always the local entry number in the current tree.
262 // Assuming that fChain is the pointer to the TChain being processed,
263 // use fChain->GetTree()->GetEntry(entry).
265 cout << "->AliAnalysisManager::Process()" << endl;
270 cout << "<-AliAnalysisManager::Process()" << endl;
275 //______________________________________________________________________________
276 void AliAnalysisManager::PackOutput(TList *target)
278 // Pack all output data containers in the output list. Called at SlaveTerminate
279 // stage in PROOF case for each slave.
281 cout << "->AliAnalysisManager::PackOutput()" << endl;
284 Error("PackOutput", "No target. Aborting.");
288 if (fMode == kProofAnalysis) {
289 TIter next(fOutputs);
290 AliAnalysisDataContainer *output;
291 while ((output=(AliAnalysisDataContainer*)next())) {
292 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
293 if (output->GetData()) target->Add(output->ExportData());
297 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
301 //______________________________________________________________________________
302 void AliAnalysisManager::ImportWrappers(TList *source)
304 // Import data in output containers from wrappers coming in source.
306 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
308 TIter next(fOutputs);
309 AliAnalysisDataContainer *cont;
310 AliAnalysisDataWrapper *wrap;
312 while ((cont=(AliAnalysisDataContainer*)next())) {
313 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
314 if (!wrap && fDebug>1) {
315 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
319 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
320 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
321 cont->ImportData(wrap);
324 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
328 //______________________________________________________________________________
329 void AliAnalysisManager::UnpackOutput(TList *source)
331 // Called by AliAnalysisSelector::Terminate. Output containers should
332 // be in source in the same order as in fOutputs.
334 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
337 Error("UnpackOutput", "No target. Aborting.");
341 printf(" Source list contains %d containers\n", source->GetSize());
344 if (fMode == kProofAnalysis) ImportWrappers(source);
346 TIter next(fOutputs);
347 AliAnalysisDataContainer *output;
348 while ((output=(AliAnalysisDataContainer*)next())) {
349 if (!output->GetData()) continue;
350 // Check if there are client tasks that run post event loop
351 if (output->HasConsumers()) {
352 // Disable event loop semaphore
353 output->SetPostEventLoop(kTRUE);
354 TObjArray *list = output->GetConsumers();
355 Int_t ncons = list->GetEntriesFast();
356 for (Int_t i=0; i<ncons; i++) {
357 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
358 task->CheckNotify(kTRUE);
359 // If task is active, execute it
360 if (task->IsPostEventLoop() && task->IsActive()) {
362 cout << "== Executing post event loop task " << task->GetName() << endl;
368 // Check if the output need to be written to a file.
369 const char *filename = output->GetFileName();
370 if (!filename || !strlen(filename)) continue;
371 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
372 if (file) file->cd();
373 else file = new TFile(filename, "RECREATE");
374 if (file->IsZombie()) continue;
375 // Reparent data to this file
377 if (output->GetData()->IsA())
378 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
379 if (callEnv.IsValid()) {
380 callEnv.SetParam((Long_t) file);
381 callEnv.Execute(output->GetData());
383 output->GetData()->Write();
386 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
390 //______________________________________________________________________________
391 void AliAnalysisManager::Terminate()
393 // The Terminate() function is the last function to be called during
394 // a query. It always runs on the client, it can be used to present
395 // the results graphically.
397 cout << "->AliAnalysisManager::Terminate()" << endl;
399 AliAnalysisTask *task;
401 // Call Terminate() for tasks
402 while ((task=(AliAnalysisTask*)next())) task->Terminate();
404 cout << "<-AliAnalysisManager::Terminate()" << endl;
408 //______________________________________________________________________________
409 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
411 // Adds a user task to the global list of tasks.
412 task->SetActive(kFALSE);
416 //______________________________________________________________________________
417 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
419 // Retreive task by name.
420 if (!fTasks) return NULL;
421 return (AliAnalysisTask*)fTasks->FindObject(name);
424 //______________________________________________________________________________
425 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
426 TClass *datatype, EAliAnalysisContType type, const char *filename)
428 // Create a data container of a certain type. Types can be:
429 // kExchangeContainer = 0, used to exchange date between tasks
430 // kInputContainer = 1, used to store input data
431 // kOutputContainer = 2, used for posting results
432 if (fContainers->FindObject(name)) {
433 Error("CreateContainer","A container named %s already defined !\n",name);
436 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
437 fContainers->Add(cont);
439 case kInputContainer:
442 case kOutputContainer:
444 if (filename && strlen(filename)) cont->SetFileName(filename);
446 case kExchangeContainer:
452 //______________________________________________________________________________
453 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
454 AliAnalysisDataContainer *cont)
456 // Connect input of an existing task to a data container.
457 if (!fTasks->FindObject(task)) {
459 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
461 Bool_t connected = task->ConnectInput(islot, cont);
465 //______________________________________________________________________________
466 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
467 AliAnalysisDataContainer *cont)
469 // Connect output of an existing task to a data container.
470 if (!fTasks->FindObject(task)) {
472 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
474 Bool_t connected = task->ConnectOutput(islot, cont);
478 //______________________________________________________________________________
479 void AliAnalysisManager::CleanContainers()
481 // Clean data from all containers that have already finished all client tasks.
482 TIter next(fContainers);
483 AliAnalysisDataContainer *cont;
484 while ((cont=(AliAnalysisDataContainer *)next())) {
485 if (cont->IsOwnedData() &&
486 cont->IsDataReady() &&
487 cont->ClientsExecuted()) cont->DeleteData();
491 //______________________________________________________________________________
492 Bool_t AliAnalysisManager::InitAnalysis()
494 // Initialization of analysis chain of tasks. Should be called after all tasks
495 // and data containers are properly connected
496 // Check for input/output containers
498 // Check for top tasks (depending only on input data containers)
499 if (!fTasks->First()) {
500 Error("InitAnalysis", "Analysis has no tasks !");
504 AliAnalysisTask *task;
505 AliAnalysisDataContainer *cont;
508 Bool_t iszombie = kFALSE;
509 Bool_t istop = kTRUE;
511 while ((task=(AliAnalysisTask*)next())) {
514 Int_t ninputs = task->GetNinputs();
515 for (i=0; i<ninputs; i++) {
516 cont = task->GetInputSlot(i)->GetContainer();
524 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
527 if (iszombie) continue;
528 // Check if cont is an input container
529 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
530 // Connect to parent task
534 fTopTasks->Add(task);
538 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
541 // Check now if there are orphan tasks
542 for (i=0; i<ntop; i++) {
543 task = (AliAnalysisTask*)fTopTasks->At(i);
548 while ((task=(AliAnalysisTask*)next())) {
549 if (!task->IsUsed()) {
551 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
554 // Check the task hierarchy (no parent task should depend on data provided
555 // by a daughter task)
556 for (i=0; i<ntop; i++) {
557 task = (AliAnalysisTask*)fTopTasks->At(i);
558 if (task->CheckCircularDeps()) {
559 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
564 // Check that all containers feeding post-event loop tasks are in the outputs list
565 TIter nextcont(fContainers); // loop over all containers
566 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
567 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
568 if (cont->HasConsumers()) {
569 // Check if one of the consumers is post event loop
570 TIter nextconsumer(cont->GetConsumers());
571 while ((task=(AliAnalysisTask*)nextconsumer())) {
572 if (task->IsPostEventLoop()) {
584 //______________________________________________________________________________
585 void AliAnalysisManager::PrintStatus(Option_t *option) const
587 // Print task hierarchy.
588 TIter next(fTopTasks);
589 AliAnalysisTask *task;
590 while ((task=(AliAnalysisTask*)next()))
591 task->PrintTask(option);
594 //______________________________________________________________________________
595 void AliAnalysisManager::ResetAnalysis()
597 // Reset all execution flags and clean containers.
601 //______________________________________________________________________________
602 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
604 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
606 Error("StartAnalysis","Analysis manager was not initialized !");
610 cout << "StartAnalysis: " << GetName() << endl;
612 TString anaType = type;
614 fMode = kLocalAnalysis;
616 if (anaType.Contains("proof")) fMode = kProofAnalysis;
617 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
619 if (fMode == kGridAnalysis) {
620 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
621 fMode = kLocalAnalysis;
624 SetEventLoop(kFALSE);
625 // Disable all branches if requested and set event loop mode
627 if (TestBit(kDisableBranches)) {
628 printf("Disabling all branches...\n");
629 // tree->SetBranchStatus("*",0); // not yet working
634 TChain *chain = dynamic_cast<TChain*>(tree);
636 // Initialize locally all tasks
638 AliAnalysisTask *task;
639 while ((task=(AliAnalysisTask*)next())) {
647 AliAnalysisTask *task;
648 // Call CreateOutputObjects for all tasks
649 while ((task=(AliAnalysisTask*)next())) {
650 TDirectory *curdir = gDirectory;
651 task->CreateOutputObjects();
652 if (curdir) curdir->cd();
658 // Run tree-based analysis via AliAnalysisSelector
659 gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
660 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
661 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
662 gROOT->ProcessLine(line);
663 sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
664 gROOT->ProcessLine(line);
667 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
668 printf("StartAnalysis: no PROOF!!!\n");
671 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
672 gROOT->ProcessLine(line);
675 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
676 chain->Process(gSystem->ExpandPathName("$ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+"));
678 printf("StartAnalysis: no chain\n");
683 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
687 //______________________________________________________________________________
688 void AliAnalysisManager::ExecAnalysis(Option_t *option)
692 Error("ExecAnalysis", "Analysis manager was not initialized !");
695 AliAnalysisTask *task;
696 // Check if the top tree is active.
699 printf("AliAnalysisManager::ExecAnalysis\n");
702 // De-activate all tasks
703 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
704 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
706 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
709 cont->SetData(fTree); // This will notify all consumers
710 TIter next1(cont->GetConsumers());
711 while ((task=(AliAnalysisTask*)next1())) {
712 // task->SetActive(kTRUE);
714 cout << " Executing task " << task->GetName() << endl;
716 task->ExecuteTask(option);
720 // The event loop is not controlled by TSelector
721 TIter next2(fTopTasks);
722 while ((task=(AliAnalysisTask*)next2())) {
723 task->SetActive(kTRUE);
725 cout << " Executing task " << task->GetName() << endl;
727 task->ExecuteTask(option);
731 //______________________________________________________________________________
732 void AliAnalysisManager::FinishAnalysis()