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;
64 //______________________________________________________________________________
65 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
69 fMode(kLocalAnalysis),
79 // Default constructor.
80 fgAnalysisManager = this;
81 fTasks = new TObjArray();
82 fTopTasks = new TObjArray();
83 fZombies = new TObjArray();
84 fContainers = new TObjArray();
85 fInputs = new TObjArray();
86 fOutputs = new TObjArray();
89 //______________________________________________________________________________
90 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
95 fInitOK(other.fInitOK),
105 fTasks = new TObjArray(*other.fTasks);
106 fTopTasks = new TObjArray(*other.fTopTasks);
107 fZombies = new TObjArray(*other.fZombies);
108 fContainers = new TObjArray(*other.fContainers);
109 fInputs = new TObjArray(*other.fInputs);
110 fOutputs = new TObjArray(*other.fOutputs);
111 fgAnalysisManager = this;
114 //______________________________________________________________________________
115 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
118 if (&other != this) {
119 TNamed::operator=(other);
123 fInitOK = other.fInitOK;
124 fDebug = other.fDebug;
125 fTasks = new TObjArray(*other.fTasks);
126 fTopTasks = new TObjArray(*other.fTopTasks);
127 fZombies = new TObjArray(*other.fZombies);
128 fContainers = new TObjArray(*other.fContainers);
129 fInputs = new TObjArray(*other.fInputs);
130 fOutputs = new TObjArray(*other.fOutputs);
131 fgAnalysisManager = this;
136 //______________________________________________________________________________
137 AliAnalysisManager::~AliAnalysisManager()
140 if (fTasks) {fTasks->Delete(); delete fTasks;}
141 if (fTopTasks) delete fTopTasks;
142 if (fZombies) delete fZombies;
143 if (fContainers) {fContainers->Delete(); delete fContainers;}
144 if (fInputs) delete fInputs;
145 if (fOutputs) delete fOutputs;
146 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
149 //______________________________________________________________________________
150 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
152 // Read one entry of the tree or a whole branch.
154 cout << "== AliAnalysisManager::GetEntry()" << endl;
156 fCurrentEntry = entry;
157 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
160 //______________________________________________________________________________
161 void AliAnalysisManager::Init(TTree *tree)
163 // The Init() function is called when the selector needs to initialize
164 // a new tree or chain. Typically here the branch addresses of the tree
165 // will be set. It is normaly not necessary to make changes to the
166 // generated code, but the routine can be extended by the user if needed.
167 // Init() will be called many times when running with PROOF.
170 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
172 if (!fInitOK) InitAnalysis();
173 if (!fInitOK) return;
175 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
177 cout<<"Error: No top input container !" <<endl;
182 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
186 //______________________________________________________________________________
187 void AliAnalysisManager::Begin(TTree *tree)
189 // The Begin() function is called at the start of the query.
190 // When running with PROOF Begin() is only called on the client.
191 // The tree argument is deprecated (on PROOF 0 is passed).
193 cout << "AliAnalysisManager::Begin()" << endl;
198 //______________________________________________________________________________
199 void AliAnalysisManager::SlaveBegin(TTree *tree)
201 // The SlaveBegin() function is called after the Begin() function.
202 // When running with PROOF SlaveBegin() is called on each slave server.
203 // The tree argument is deprecated (on PROOF 0 is passed).
205 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
209 AliAnalysisTask *task;
210 // Call CreateOutputObjects for all tasks
211 while ((task=(AliAnalysisTask*)next()))
212 task->CreateOutputObjects();
213 if (fMode == kLocalAnalysis) Init(tree);
215 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
219 //______________________________________________________________________________
220 Bool_t AliAnalysisManager::Notify()
222 // The Notify() function is called when a new file is opened. This
223 // can be either for a new TTree in a TChain or when when a new TTree
224 // is started when using PROOF. It is normaly not necessary to make changes
225 // to the generated code, but the routine can be extended by the
226 // user if needed. The return value is currently not used.
228 TFile *curfile = fTree->GetCurrentFile();
229 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
231 AliAnalysisTask *task;
232 // Call Notify for all tasks
233 while ((task=(AliAnalysisTask*)next()))
239 //______________________________________________________________________________
240 Bool_t AliAnalysisManager::Process(Long64_t entry)
242 // The Process() function is called for each entry in the tree (or possibly
243 // keyed object in the case of PROOF) to be processed. The entry argument
244 // specifies which entry in the currently loaded tree is to be processed.
245 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
246 // to read either all or the required parts of the data. When processing
247 // keyed objects with PROOF, the object is already loaded and is available
248 // via the fObject pointer.
250 // This function should contain the "body" of the analysis. It can contain
251 // simple or elaborate selection criteria, run algorithms on the data
252 // of the event and typically fill histograms.
254 // WARNING when a selector is used with a TChain, you must use
255 // the pointer to the current TTree to call GetEntry(entry).
256 // The entry is always the local entry number in the current tree.
257 // Assuming that fChain is the pointer to the TChain being processed,
258 // use fChain->GetTree()->GetEntry(entry).
260 cout << "->AliAnalysisManager::Process()" << endl;
265 cout << "<-AliAnalysisManager::Process()" << endl;
270 //______________________________________________________________________________
271 void AliAnalysisManager::PackOutput(TList *target)
273 // Pack all output data containers in the output list. Called at SlaveTerminate
274 // stage in PROOF case for each slave.
276 cout << "->AliAnalysisManager::PackOutput()" << endl;
279 Error("PackOutput", "No target. Aborting.");
283 if (fMode == kProofAnalysis) {
284 TIter next(fOutputs);
285 AliAnalysisDataContainer *output;
286 while ((output=(AliAnalysisDataContainer*)next())) {
287 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
288 if (output->GetData()) target->Add(output->ExportData());
292 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
296 //______________________________________________________________________________
297 void AliAnalysisManager::ImportWrappers(TList *source)
299 // Import data in output containers from wrappers coming in source.
301 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
303 TIter next(fOutputs);
304 AliAnalysisDataContainer *cont;
305 AliAnalysisDataWrapper *wrap;
307 while ((cont=(AliAnalysisDataContainer*)next())) {
308 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
309 if (!wrap && fDebug>1) {
310 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
314 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
315 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
316 cont->ImportData(wrap);
319 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
323 //______________________________________________________________________________
324 void AliAnalysisManager::UnpackOutput(TList *source)
326 // Called by AliAnalysisSelector::Terminate. Output containers should
327 // be in source in the same order as in fOutputs.
329 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
332 Error("UnpackOutput", "No target. Aborting.");
336 printf(" Source list contains %d containers\n", source->GetSize());
339 if (fMode == kProofAnalysis) ImportWrappers(source);
341 TIter next(fOutputs);
342 AliAnalysisDataContainer *output;
343 while ((output=(AliAnalysisDataContainer*)next())) {
344 if (!output->GetData()) continue;
345 // Check if the output need to be written to a file.
346 const char *filename = output->GetFileName();
347 if (!filename || !strlen(filename)) continue;
348 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
349 if (file) file->cd();
350 else file = new TFile(filename, "RECREATE");
351 if (file->IsZombie()) continue;
352 // Reparent data to this file
354 if (output->GetData()->IsA())
355 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
356 if (callEnv.IsValid()) {
357 callEnv.SetParam((Long_t) file);
358 callEnv.Execute(output->GetData());
360 output->GetData()->Write();
363 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
367 //______________________________________________________________________________
368 void AliAnalysisManager::Terminate()
370 // The Terminate() function is the last function to be called during
371 // a query. It always runs on the client, it can be used to present
372 // the results graphically.
374 cout << "->AliAnalysisManager::Terminate()" << endl;
376 AliAnalysisTask *task;
378 // Call Terminate() for tasks
379 while ((task=(AliAnalysisTask*)next())) task->Terminate();
381 cout << "<-AliAnalysisManager::Terminate()" << endl;
385 //______________________________________________________________________________
386 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
388 // Adds a user task to the global list of tasks.
389 task->SetActive(kFALSE);
393 //______________________________________________________________________________
394 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
396 // Retreive task by name.
397 if (!fTasks) return NULL;
398 return (AliAnalysisTask*)fTasks->FindObject(name);
401 //______________________________________________________________________________
402 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
403 TClass *datatype, EAliAnalysisContType type, const char *filename)
405 // Create a data container of a certain type. Types can be:
406 // kExchangeContainer = 0, used to exchange date between tasks
407 // kInputContainer = 1, used to store input data
408 // kOutputContainer = 2, used for posting results
409 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
410 fContainers->Add(cont);
412 case kInputContainer:
415 case kOutputContainer:
416 if (fOutputs->FindObject(name)) printf("CreateContainer: warning: a container named %s existing !\n",name);
418 if (filename && strlen(filename)) cont->SetFileName(filename);
420 case kExchangeContainer:
426 //______________________________________________________________________________
427 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
428 AliAnalysisDataContainer *cont)
430 // Connect input of an existing task to a data container.
431 if (!fTasks->FindObject(task)) {
433 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
435 Bool_t connected = task->ConnectInput(islot, cont);
439 //______________________________________________________________________________
440 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
441 AliAnalysisDataContainer *cont)
443 // Connect output of an existing task to a data container.
444 if (!fTasks->FindObject(task)) {
446 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
448 Bool_t connected = task->ConnectOutput(islot, cont);
452 //______________________________________________________________________________
453 void AliAnalysisManager::CleanContainers()
455 // Clean data from all containers that have already finished all client tasks.
456 TIter next(fContainers);
457 AliAnalysisDataContainer *cont;
458 while ((cont=(AliAnalysisDataContainer *)next())) {
459 if (cont->IsOwnedData() &&
460 cont->IsDataReady() &&
461 cont->ClientsExecuted()) cont->DeleteData();
465 //______________________________________________________________________________
466 Bool_t AliAnalysisManager::InitAnalysis()
468 // Initialization of analysis chain of tasks. Should be called after all tasks
469 // and data containers are properly connected
470 // Check for input/output containers
472 // Check for top tasks (depending only on input data containers)
473 if (!fTasks->First()) {
474 Error("InitAnalysis", "Analysis has no tasks !");
478 AliAnalysisTask *task;
479 AliAnalysisDataContainer *cont;
482 Bool_t iszombie = kFALSE;
483 Bool_t istop = kTRUE;
485 while ((task=(AliAnalysisTask*)next())) {
488 Int_t ninputs = task->GetNinputs();
489 for (i=0; i<ninputs; i++) {
490 cont = task->GetInputSlot(i)->GetContainer();
498 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
501 if (iszombie) continue;
502 // Check if cont is an input container
503 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
504 // Connect to parent task
508 fTopTasks->Add(task);
512 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
515 // Check now if there are orphan tasks
516 for (i=0; i<ntop; i++) {
517 task = (AliAnalysisTask*)fTopTasks->At(i);
522 while ((task=(AliAnalysisTask*)next())) {
523 if (!task->IsUsed()) {
525 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
528 // Check the task hierarchy (no parent task should depend on data provided
529 // by a daughter task)
530 for (i=0; i<ntop; i++) {
531 task = (AliAnalysisTask*)fTopTasks->At(i);
532 if (task->CheckCircularDeps()) {
533 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
542 //______________________________________________________________________________
543 void AliAnalysisManager::PrintStatus(Option_t *option) const
545 // Print task hierarchy.
546 TIter next(fTopTasks);
547 AliAnalysisTask *task;
548 while ((task=(AliAnalysisTask*)next()))
549 task->PrintTask(option);
552 //______________________________________________________________________________
553 void AliAnalysisManager::ResetAnalysis()
555 // Reset all execution flags and clean containers.
559 //______________________________________________________________________________
560 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
562 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
564 Error("StartAnalysis","Analysis manager was not initialized !");
568 cout << "StartAnalysis: " << GetName() << endl;
570 TString anaType = type;
572 fMode = kLocalAnalysis;
574 if (anaType.Contains("proof")) fMode = kProofAnalysis;
575 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
577 if (fMode == kGridAnalysis) {
578 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
579 fMode = kLocalAnalysis;
582 // Disable by default all branches
583 if (tree) tree->SetBranchStatus("*",0);
584 TChain *chain = dynamic_cast<TChain*>(tree);
586 // Initialize locally all tasks
588 AliAnalysisTask *task;
589 while ((task=(AliAnalysisTask*)next())) task->LocalInit();
595 AliAnalysisTask *task;
596 // Call CreateOutputObjects for all tasks
597 while ((task=(AliAnalysisTask*)next())) task->CreateOutputObjects();
602 // Run tree-based analysis via AliAnalysisSelector
603 gROOT->ProcessLine(".L AliAnalysisSelector.cxx+");
604 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
605 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
606 gROOT->ProcessLine(line);
607 sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
608 gROOT->ProcessLine(line);
611 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
612 printf("StartAnalysis: no PROOF!!!\n");
615 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
616 gROOT->ProcessLine(line);
619 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
620 chain->Process(gSystem->ExpandPathName("$ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+"));
622 printf("StartAnalysis: no chain\n");
627 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
631 //______________________________________________________________________________
632 void AliAnalysisManager::ExecAnalysis(Option_t *option)
636 Error("ExecAnalysis", "Analysis manager was not initialized !");
639 AliAnalysisTask *task;
640 // Check if the top tree is active.
643 printf("AliAnalysisManager::ExecAnalysis\n");
646 // De-activate all tasks
647 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
648 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
650 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
653 cont->SetData(fTree); // This will notify all consumers
654 TIter next1(cont->GetConsumers());
655 while ((task=(AliAnalysisTask*)next1())) {
656 // task->SetActive(kTRUE);
658 cout << " Executing task " << task->GetName() << endl;
660 task->ExecuteTask(option);
664 // The event loop is not controlled by TSelector
665 TIter next2(fTopTasks);
666 while ((task=(AliAnalysisTask*)next2())) {
667 task->SetActive(kTRUE);
669 cout << " Executing task " << task->GetName() << endl;
671 task->ExecuteTask(option);
675 //______________________________________________________________________________
676 void AliAnalysisManager::FinishAnalysis()