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(const char *name, const char *title)
51 fInputEventHandler(NULL),
52 fOutputEventHandler(NULL),
53 fMCtruthEventHandler(NULL),
55 fMode(kLocalAnalysis),
65 // Default constructor.
66 fgAnalysisManager = this;
67 fTasks = new TObjArray();
68 fTopTasks = new TObjArray();
69 fZombies = new TObjArray();
70 fContainers = new TObjArray();
71 fInputs = new TObjArray();
72 fOutputs = new TObjArray();
76 //______________________________________________________________________________
77 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
80 fInputEventHandler(NULL),
81 fOutputEventHandler(NULL),
82 fMCtruthEventHandler(NULL),
85 fInitOK(other.fInitOK),
95 fTasks = new TObjArray(*other.fTasks);
96 fTopTasks = new TObjArray(*other.fTopTasks);
97 fZombies = new TObjArray(*other.fZombies);
98 fContainers = new TObjArray(*other.fContainers);
99 fInputs = new TObjArray(*other.fInputs);
100 fOutputs = new TObjArray(*other.fOutputs);
101 fgAnalysisManager = this;
104 //______________________________________________________________________________
105 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
108 if (&other != this) {
109 TNamed::operator=(other);
110 fInputEventHandler = other.fInputEventHandler;
111 fOutputEventHandler = other.fOutputEventHandler;
112 fMCtruthEventHandler = other.fMCtruthEventHandler;
116 fInitOK = other.fInitOK;
117 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;
129 //______________________________________________________________________________
130 AliAnalysisManager::~AliAnalysisManager()
133 if (fTasks) {fTasks->Delete(); delete fTasks;}
134 if (fTopTasks) delete fTopTasks;
135 if (fZombies) delete fZombies;
136 if (fContainers) {fContainers->Delete(); delete fContainers;}
137 if (fInputs) delete fInputs;
138 if (fOutputs) delete fOutputs;
139 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
142 //______________________________________________________________________________
143 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
145 // Read one entry of the tree or a whole branch.
147 cout << "== AliAnalysisManager::GetEntry()" << endl;
149 fCurrentEntry = entry;
150 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
153 //______________________________________________________________________________
154 void AliAnalysisManager::Init(TTree *tree)
156 // The Init() function is called when the selector needs to initialize
157 // a new tree or chain. Typically here the branch addresses of the tree
158 // will be set. It is normaly not necessary to make changes to the
159 // generated code, but the routine can be extended by the user if needed.
160 // Init() will be called many times when running with PROOF.
163 printf("->AliAnalysisManager::InitTree(%s)\n", tree->GetName());
166 // Call InitTree of EventHandler
167 if (fOutputEventHandler) {
168 if (fMode == kProofAnalysis) {
169 fOutputEventHandler->Init(0x0, "proof");
171 fOutputEventHandler->Init(0x0, "local");
175 if (fInputEventHandler) {
176 if (fMode == kProofAnalysis) {
177 fInputEventHandler->Init(tree, "proof");
179 fInputEventHandler->Init(tree, "local");
182 // If no input event handler we need to get the tree once
184 if(!tree->GetTree()) tree->LoadTree(0);
188 if (fMCtruthEventHandler) {
189 if (fMode == kProofAnalysis) {
190 fMCtruthEventHandler->Init(0x0, "proof");
192 fMCtruthEventHandler->Init(0x0, "local");
196 if (!fInitOK) InitAnalysis();
197 if (!fInitOK) return;
199 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
201 cout<<"Error: No top input container !" <<endl;
206 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
210 //______________________________________________________________________________
211 void AliAnalysisManager::SlaveBegin(TTree *tree)
213 // The SlaveBegin() function is called after the Begin() function.
214 // When running with PROOF SlaveBegin() is called on each slave server.
215 // The tree argument is deprecated (on PROOF 0 is passed).
217 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
220 // Call Init of EventHandler
221 if (fOutputEventHandler) {
222 if (fMode == kProofAnalysis) {
223 fOutputEventHandler->Init("proof");
225 fOutputEventHandler->Init("local");
229 if (fInputEventHandler) {
230 fInputEventHandler->SetInputTree(tree);
231 if (fMode == kProofAnalysis) {
232 fInputEventHandler->Init("proof");
234 fInputEventHandler->Init("local");
238 if (fMCtruthEventHandler) {
239 if (fMode == kProofAnalysis) {
240 fMCtruthEventHandler->Init("proof");
242 fMCtruthEventHandler->Init("local");
247 AliAnalysisTask *task;
248 // Call CreateOutputObjects for all tasks
249 while ((task=(AliAnalysisTask*)next())) {
250 TDirectory *curdir = gDirectory;
251 task->CreateOutputObjects();
252 if (curdir) curdir->cd();
256 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
260 //______________________________________________________________________________
261 Bool_t AliAnalysisManager::Notify()
263 // The Notify() function is called when a new file is opened. This
264 // can be either for a new TTree in a TChain or when when a new TTree
265 // is started when using PROOF. It is normaly not necessary to make changes
266 // to the generated code, but the routine can be extended by the
267 // user if needed. The return value is currently not used.
269 TFile *curfile = fTree->GetCurrentFile();
270 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
272 AliAnalysisTask *task;
273 // Call Notify for all tasks
274 while ((task=(AliAnalysisTask*)next()))
277 // Call Notify of the event handlers
278 if (fInputEventHandler) {
279 fInputEventHandler->Notify(curfile->GetName());
282 if (fOutputEventHandler) {
283 fOutputEventHandler->Notify(curfile->GetName());
286 if (fMCtruthEventHandler) {
287 fMCtruthEventHandler->Notify(curfile->GetName());
294 //______________________________________________________________________________
295 Bool_t AliAnalysisManager::Process(Long64_t entry)
297 // The Process() function is called for each entry in the tree (or possibly
298 // keyed object in the case of PROOF) to be processed. The entry argument
299 // specifies which entry in the currently loaded tree is to be processed.
300 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
301 // to read either all or the required parts of the data. When processing
302 // keyed objects with PROOF, the object is already loaded and is available
303 // via the fObject pointer.
305 // This function should contain the "body" of the analysis. It can contain
306 // simple or elaborate selection criteria, run algorithms on the data
307 // of the event and typically fill histograms.
309 // WARNING when a selector is used with a TChain, you must use
310 // the pointer to the current TTree to call GetEntry(entry).
311 // The entry is always the local entry number in the current tree.
312 // Assuming that fChain is the pointer to the TChain being processed,
313 // use fChain->GetTree()->GetEntry(entry).
315 cout << "->AliAnalysisManager::Process(" << entry << ")" << endl;
317 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
318 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
319 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
324 cout << "<-AliAnalysisManager::Process()" << endl;
329 //______________________________________________________________________________
330 void AliAnalysisManager::PackOutput(TList *target)
332 // Pack all output data containers in the output list. Called at SlaveTerminate
333 // stage in PROOF case for each slave.
335 cout << "->AliAnalysisManager::PackOutput()" << endl;
338 Error("PackOutput", "No target. Aborting.");
341 if (fInputEventHandler) fInputEventHandler ->Terminate();
342 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
343 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
345 if (fMode == kProofAnalysis) {
346 TIter next(fOutputs);
347 AliAnalysisDataContainer *output;
348 while ((output=(AliAnalysisDataContainer*)next())) {
349 if (output->GetData()) {
350 if (output->GetProducer()->IsPostEventLoop()) continue;
351 AliAnalysisDataWrapper *wrap = output->ExportData();
352 // Output wrappers must delete data after merging (AG 13/11/07)
353 wrap->SetDeleteData(kTRUE);
354 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
360 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
364 //______________________________________________________________________________
365 void AliAnalysisManager::ImportWrappers(TList *source)
367 // Import data in output containers from wrappers coming in source.
369 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
371 TIter next(fOutputs);
372 AliAnalysisDataContainer *cont;
373 AliAnalysisDataWrapper *wrap;
375 while ((cont=(AliAnalysisDataContainer*)next())) {
376 if (cont->GetProducer()->IsPostEventLoop()) continue;
377 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
378 if (!wrap && fDebug>1) {
379 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
383 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
384 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
385 cont->ImportData(wrap);
388 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
392 //______________________________________________________________________________
393 void AliAnalysisManager::UnpackOutput(TList *source)
395 // Called by AliAnalysisSelector::Terminate. Output containers should
396 // be in source in the same order as in fOutputs.
398 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
401 Error("UnpackOutput", "No target. Aborting.");
405 printf(" Source list contains %d containers\n", source->GetSize());
408 if (fMode == kProofAnalysis) ImportWrappers(source);
410 TIter next(fOutputs);
411 AliAnalysisDataContainer *output;
412 while ((output=(AliAnalysisDataContainer*)next())) {
413 if (!output->GetData()) continue;
414 // Check if there are client tasks that run post event loop
415 if (output->HasConsumers()) {
416 // Disable event loop semaphore
417 output->SetPostEventLoop(kTRUE);
418 TObjArray *list = output->GetConsumers();
419 Int_t ncons = list->GetEntriesFast();
420 for (Int_t i=0; i<ncons; i++) {
421 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
422 task->CheckNotify(kTRUE);
423 // If task is active, execute it
424 if (task->IsPostEventLoop() && task->IsActive()) {
426 cout << "== Executing post event loop task " << task->GetName() << endl;
432 // Check if the output need to be written to a file.
433 const char *filename = output->GetFileName();
434 if (!(strcmp(filename, "default"))) {
435 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
438 if (!filename || !strlen(filename)) continue;
439 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
440 TDirectory *opwd = gDirectory;
441 if (file) file->cd();
442 else file = new TFile(filename, "RECREATE");
443 if (file->IsZombie()) continue;
444 // Reparent data to this file
446 if (output->GetData()->IsA())
447 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
448 if (callEnv.IsValid()) {
449 callEnv.SetParam((Long_t) file);
450 callEnv.Execute(output->GetData());
452 output->GetData()->Write();
454 if (opwd) opwd->cd();
457 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
461 //______________________________________________________________________________
462 void AliAnalysisManager::Terminate()
464 // The Terminate() function is the last function to be called during
465 // a query. It always runs on the client, it can be used to present
466 // the results graphically.
468 cout << "->AliAnalysisManager::Terminate()" << endl;
470 AliAnalysisTask *task;
472 // Call Terminate() for tasks
473 while ((task=(AliAnalysisTask*)next())) task->Terminate();
475 cout << "<-AliAnalysisManager::Terminate()" << endl;
478 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
479 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
480 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
483 //______________________________________________________________________________
484 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
486 // Adds a user task to the global list of tasks.
487 task->SetActive(kFALSE);
491 //______________________________________________________________________________
492 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
494 // Retreive task by name.
495 if (!fTasks) return NULL;
496 return (AliAnalysisTask*)fTasks->FindObject(name);
499 //______________________________________________________________________________
500 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
501 TClass *datatype, EAliAnalysisContType type, const char *filename)
503 // Create a data container of a certain type. Types can be:
504 // kExchangeContainer = 0, used to exchange date between tasks
505 // kInputContainer = 1, used to store input data
506 // kOutputContainer = 2, used for posting results
507 if (fContainers->FindObject(name)) {
508 Error("CreateContainer","A container named %s already defined !\n",name);
511 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
512 fContainers->Add(cont);
514 case kInputContainer:
517 case kOutputContainer:
519 if (filename && strlen(filename)) cont->SetFileName(filename);
521 case kExchangeContainer:
527 //______________________________________________________________________________
528 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
529 AliAnalysisDataContainer *cont)
531 // Connect input of an existing task to a data container.
532 if (!fTasks->FindObject(task)) {
534 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
536 Bool_t connected = task->ConnectInput(islot, cont);
540 //______________________________________________________________________________
541 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
542 AliAnalysisDataContainer *cont)
544 // Connect output of an existing task to a data container.
545 if (!fTasks->FindObject(task)) {
547 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
549 Bool_t connected = task->ConnectOutput(islot, cont);
553 //______________________________________________________________________________
554 void AliAnalysisManager::CleanContainers()
556 // Clean data from all containers that have already finished all client tasks.
557 TIter next(fContainers);
558 AliAnalysisDataContainer *cont;
559 while ((cont=(AliAnalysisDataContainer *)next())) {
560 if (cont->IsOwnedData() &&
561 cont->IsDataReady() &&
562 cont->ClientsExecuted()) cont->DeleteData();
566 //______________________________________________________________________________
567 Bool_t AliAnalysisManager::InitAnalysis()
569 // Initialization of analysis chain of tasks. Should be called after all tasks
570 // and data containers are properly connected
571 // Check for input/output containers
573 // Check for top tasks (depending only on input data containers)
574 if (!fTasks->First()) {
575 Error("InitAnalysis", "Analysis has no tasks !");
579 AliAnalysisTask *task;
580 AliAnalysisDataContainer *cont;
583 Bool_t iszombie = kFALSE;
584 Bool_t istop = kTRUE;
586 while ((task=(AliAnalysisTask*)next())) {
589 Int_t ninputs = task->GetNinputs();
590 for (i=0; i<ninputs; i++) {
591 cont = task->GetInputSlot(i)->GetContainer();
599 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
602 if (iszombie) continue;
603 // Check if cont is an input container
604 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
605 // Connect to parent task
609 fTopTasks->Add(task);
613 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
616 // Check now if there are orphan tasks
617 for (i=0; i<ntop; i++) {
618 task = (AliAnalysisTask*)fTopTasks->At(i);
623 while ((task=(AliAnalysisTask*)next())) {
624 if (!task->IsUsed()) {
626 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
629 // Check the task hierarchy (no parent task should depend on data provided
630 // by a daughter task)
631 for (i=0; i<ntop; i++) {
632 task = (AliAnalysisTask*)fTopTasks->At(i);
633 if (task->CheckCircularDeps()) {
634 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
639 // Check that all containers feeding post-event loop tasks are in the outputs list
640 TIter nextcont(fContainers); // loop over all containers
641 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
642 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
643 if (cont->HasConsumers()) {
644 // Check if one of the consumers is post event loop
645 TIter nextconsumer(cont->GetConsumers());
646 while ((task=(AliAnalysisTask*)nextconsumer())) {
647 if (task->IsPostEventLoop()) {
659 //______________________________________________________________________________
660 void AliAnalysisManager::PrintStatus(Option_t *option) const
662 // Print task hierarchy.
663 TIter next(fTopTasks);
664 AliAnalysisTask *task;
665 while ((task=(AliAnalysisTask*)next()))
666 task->PrintTask(option);
669 //______________________________________________________________________________
670 void AliAnalysisManager::ResetAnalysis()
672 // Reset all execution flags and clean containers.
676 //______________________________________________________________________________
677 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
679 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
681 Error("StartAnalysis","Analysis manager was not initialized !");
685 cout << "StartAnalysis: " << GetName() << endl;
687 TString anaType = type;
689 fMode = kLocalAnalysis;
691 if (anaType.Contains("proof")) fMode = kProofAnalysis;
692 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
694 if (fMode == kGridAnalysis) {
695 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
696 fMode = kLocalAnalysis;
699 SetEventLoop(kFALSE);
700 // Disable all branches if requested and set event loop mode
702 if (TestBit(kDisableBranches)) {
703 printf("Disabling all branches...\n");
704 // tree->SetBranchStatus("*",0); // not yet working
709 TChain *chain = dynamic_cast<TChain*>(tree);
711 // Initialize locally all tasks
713 AliAnalysisTask *task;
714 while ((task=(AliAnalysisTask*)next())) {
722 AliAnalysisTask *task;
723 // Call CreateOutputObjects for all tasks
724 while ((task=(AliAnalysisTask*)next())) {
725 TDirectory *curdir = gDirectory;
726 task->CreateOutputObjects();
727 if (curdir) curdir->cd();
733 // Run tree-based analysis via AliAnalysisSelector
734 // gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
735 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
736 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
737 gROOT->ProcessLine(line);
738 sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
739 gROOT->ProcessLine(line);
742 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
743 printf("StartAnalysis: no PROOF!!!\n");
746 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
747 gROOT->ProcessLine(line);
750 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
751 chain->Process("AliAnalysisSelector");
753 printf("StartAnalysis: no chain\n");
758 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
762 //______________________________________________________________________________
763 void AliAnalysisManager::ExecAnalysis(Option_t *option)
767 Error("ExecAnalysis", "Analysis manager was not initialized !");
770 AliAnalysisTask *task;
771 // Check if the top tree is active.
774 // De-activate all tasks
775 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
776 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
778 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
781 cont->SetData(fTree); // This will notify all consumers
782 Long64_t entry = fTree->GetTree()->GetReadEntry();
785 // Call BeginEvent() for optional input/output and MC services
786 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
787 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
788 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
791 TIter next1(cont->GetConsumers());
792 while ((task=(AliAnalysisTask*)next1())) {
794 cout << " Executing task " << task->GetName() << endl;
797 task->ExecuteTask(option);
800 // Call FinishEvent() for optional output and MC services
801 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
802 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
803 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
807 // The event loop is not controlled by TSelector
809 // Call BeginEvent() for optional input/output and MC services
810 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
811 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
812 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
813 TIter next2(fTopTasks);
814 while ((task=(AliAnalysisTask*)next2())) {
815 task->SetActive(kTRUE);
817 cout << " Executing task " << task->GetName() << endl;
819 task->ExecuteTask(option);
822 // Call FinishEvent() for optional output and MC services
823 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
824 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
825 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
828 //______________________________________________________________________________
829 void AliAnalysisManager::FinishAnalysis()