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()
51 fInputEventHandler(NULL),
52 fOutputEventHandler(NULL),
53 fMCtruthEventHandler(NULL),
55 fMode(kLocalAnalysis),
66 fgAnalysisManager = this;
70 //______________________________________________________________________________
71 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
74 fInputEventHandler(NULL),
75 fOutputEventHandler(NULL),
76 fMCtruthEventHandler(NULL),
78 fMode(kLocalAnalysis),
88 // Default constructor.
89 fgAnalysisManager = this;
90 fTasks = new TObjArray();
91 fTopTasks = new TObjArray();
92 fZombies = new TObjArray();
93 fContainers = new TObjArray();
94 fInputs = new TObjArray();
95 fOutputs = new TObjArray();
99 //______________________________________________________________________________
100 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
103 fInputEventHandler(NULL),
104 fOutputEventHandler(NULL),
105 fMCtruthEventHandler(NULL),
108 fInitOK(other.fInitOK),
109 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;
127 //______________________________________________________________________________
128 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
131 if (&other != this) {
132 TNamed::operator=(other);
133 fInputEventHandler = other.fInputEventHandler;
134 fOutputEventHandler = other.fOutputEventHandler;
135 fMCtruthEventHandler = other.fMCtruthEventHandler;
139 fInitOK = other.fInitOK;
140 fDebug = other.fDebug;
141 fTasks = new TObjArray(*other.fTasks);
142 fTopTasks = new TObjArray(*other.fTopTasks);
143 fZombies = new TObjArray(*other.fZombies);
144 fContainers = new TObjArray(*other.fContainers);
145 fInputs = new TObjArray(*other.fInputs);
146 fOutputs = new TObjArray(*other.fOutputs);
147 fgAnalysisManager = this;
152 //______________________________________________________________________________
153 AliAnalysisManager::~AliAnalysisManager()
156 if (fTasks) {fTasks->Delete(); delete fTasks;}
157 if (fTopTasks) delete fTopTasks;
158 if (fZombies) delete fZombies;
159 if (fContainers) {fContainers->Delete(); delete fContainers;}
160 if (fInputs) delete fInputs;
161 if (fOutputs) delete fOutputs;
162 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
165 //______________________________________________________________________________
166 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
168 // Read one entry of the tree or a whole branch.
170 cout << "== AliAnalysisManager::GetEntry()" << endl;
172 fCurrentEntry = entry;
173 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
176 //______________________________________________________________________________
177 void AliAnalysisManager::Init(TTree *tree)
179 // The Init() function is called when the selector needs to initialize
180 // a new tree or chain. Typically here the branch addresses of the tree
181 // will be set. It is normaly not necessary to make changes to the
182 // generated code, but the routine can be extended by the user if needed.
183 // Init() will be called many times when running with PROOF.
186 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
189 if (fInputEventHandler) {
190 fInputEventHandler->SetInputTree(tree);
191 fInputEventHandler->InitIO("proof");
194 if (!fInitOK) InitAnalysis();
195 if (!fInitOK) return;
197 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
199 cout<<"Error: No top input container !" <<endl;
204 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
208 //______________________________________________________________________________
209 void AliAnalysisManager::Begin(TTree *tree)
211 // The Begin() function is called at the start of the query.
212 // When running with PROOF Begin() is only called on the client.
213 // The tree argument is deprecated (on PROOF 0 is passed).
215 cout << "AliAnalysisManager::Begin()" << endl;
220 //______________________________________________________________________________
221 void AliAnalysisManager::SlaveBegin(TTree *tree)
223 // The SlaveBegin() function is called after the Begin() function.
224 // When running with PROOF SlaveBegin() is called on each slave server.
225 // The tree argument is deprecated (on PROOF 0 is passed).
227 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
229 // Call InitIO of EventHandler
230 if (fOutputEventHandler) {
231 if (fMode == kProofAnalysis) {
232 fOutputEventHandler->InitIO("proof");
234 fOutputEventHandler->InitIO("local");
237 if (fInputEventHandler) {
238 if (fMode == kProofAnalysis) {
239 fInputEventHandler->SetInputTree(tree);
240 fInputEventHandler->InitIO("proof");
242 fInputEventHandler->SetInputTree(tree);
243 fInputEventHandler->InitIO("local");
247 if (fMCtruthEventHandler) {
248 if (fMode == kProofAnalysis) {
249 fMCtruthEventHandler->InitIO("proof");
251 fMCtruthEventHandler->InitIO("local");
257 AliAnalysisTask *task;
258 // Call CreateOutputObjects for all tasks
259 while ((task=(AliAnalysisTask*)next())) {
260 TDirectory *curdir = gDirectory;
261 task->CreateOutputObjects();
262 if (curdir) curdir->cd();
264 if (fMode == kLocalAnalysis)
267 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
271 //______________________________________________________________________________
272 Bool_t AliAnalysisManager::Notify()
274 // The Notify() function is called when a new file is opened. This
275 // can be either for a new TTree in a TChain or when when a new TTree
276 // is started when using PROOF. It is normaly not necessary to make changes
277 // to the generated code, but the routine can be extended by the
278 // user if needed. The return value is currently not used.
280 TFile *curfile = fTree->GetCurrentFile();
281 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
283 AliAnalysisTask *task;
284 // Call Notify for all tasks
285 while ((task=(AliAnalysisTask*)next()))
288 // Call Notify of the event handlers
289 if (fInputEventHandler) {
290 fInputEventHandler->Notify(curfile->GetName());
293 if (fOutputEventHandler) {
294 fOutputEventHandler->Notify(curfile->GetName());
297 if (fMCtruthEventHandler) {
298 fMCtruthEventHandler->Notify(curfile->GetName());
305 //______________________________________________________________________________
306 Bool_t AliAnalysisManager::Process(Long64_t entry)
308 // The Process() function is called for each entry in the tree (or possibly
309 // keyed object in the case of PROOF) to be processed. The entry argument
310 // specifies which entry in the currently loaded tree is to be processed.
311 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
312 // to read either all or the required parts of the data. When processing
313 // keyed objects with PROOF, the object is already loaded and is available
314 // via the fObject pointer.
316 // This function should contain the "body" of the analysis. It can contain
317 // simple or elaborate selection criteria, run algorithms on the data
318 // of the event and typically fill histograms.
320 // WARNING when a selector is used with a TChain, you must use
321 // the pointer to the current TTree to call GetEntry(entry).
322 // The entry is always the local entry number in the current tree.
323 // Assuming that fChain is the pointer to the TChain being processed,
324 // use fChain->GetTree()->GetEntry(entry).
326 cout << "->AliAnalysisManager::Process()" << endl;
328 if (fInputEventHandler) fInputEventHandler ->BeginEvent();
329 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent();
330 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
335 cout << "<-AliAnalysisManager::Process()" << endl;
340 //______________________________________________________________________________
341 void AliAnalysisManager::PackOutput(TList *target)
343 // Pack all output data containers in the output list. Called at SlaveTerminate
344 // stage in PROOF case for each slave.
346 cout << "->AliAnalysisManager::PackOutput()" << endl;
349 Error("PackOutput", "No target. Aborting.");
352 if (fInputEventHandler) fInputEventHandler ->Terminate();
353 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
354 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
356 if (fMode == kProofAnalysis) {
357 TIter next(fOutputs);
358 AliAnalysisDataContainer *output;
359 while ((output=(AliAnalysisDataContainer*)next())) {
360 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
361 if (output->GetData()) target->Add(output->ExportData());
365 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
369 //______________________________________________________________________________
370 void AliAnalysisManager::ImportWrappers(TList *source)
372 // Import data in output containers from wrappers coming in source.
374 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
376 TIter next(fOutputs);
377 AliAnalysisDataContainer *cont;
378 AliAnalysisDataWrapper *wrap;
380 while ((cont=(AliAnalysisDataContainer*)next())) {
381 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
382 if (!wrap && fDebug>1) {
383 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
387 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
388 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
389 cont->ImportData(wrap);
392 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
396 //______________________________________________________________________________
397 void AliAnalysisManager::UnpackOutput(TList *source)
399 // Called by AliAnalysisSelector::Terminate. Output containers should
400 // be in source in the same order as in fOutputs.
402 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
405 Error("UnpackOutput", "No target. Aborting.");
409 printf(" Source list contains %d containers\n", source->GetSize());
412 if (fMode == kProofAnalysis) ImportWrappers(source);
414 TIter next(fOutputs);
415 AliAnalysisDataContainer *output;
416 while ((output=(AliAnalysisDataContainer*)next())) {
417 if (!output->GetData()) continue;
418 // Check if there are client tasks that run post event loop
419 if (output->HasConsumers()) {
420 // Disable event loop semaphore
421 output->SetPostEventLoop(kTRUE);
422 TObjArray *list = output->GetConsumers();
423 Int_t ncons = list->GetEntriesFast();
424 for (Int_t i=0; i<ncons; i++) {
425 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
426 task->CheckNotify(kTRUE);
427 // If task is active, execute it
428 if (task->IsPostEventLoop() && task->IsActive()) {
430 cout << "== Executing post event loop task " << task->GetName() << endl;
436 // Check if the output need to be written to a file.
437 const char *filename = output->GetFileName();
438 if (!(strcmp(filename, "default"))) {
439 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
442 if (!filename || !strlen(filename)) continue;
443 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
444 if (file) file->cd();
445 else file = new TFile(filename, "RECREATE");
446 if (file->IsZombie()) continue;
447 // Reparent data to this file
449 if (output->GetData()->IsA())
450 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
451 if (callEnv.IsValid()) {
452 callEnv.SetParam((Long_t) file);
453 callEnv.Execute(output->GetData());
455 output->GetData()->Write();
458 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
462 //______________________________________________________________________________
463 void AliAnalysisManager::Terminate()
465 // The Terminate() function is the last function to be called during
466 // a query. It always runs on the client, it can be used to present
467 // the results graphically.
469 cout << "->AliAnalysisManager::Terminate()" << endl;
471 AliAnalysisTask *task;
473 // Call Terminate() for tasks
474 while ((task=(AliAnalysisTask*)next())) task->Terminate();
476 cout << "<-AliAnalysisManager::Terminate()" << endl;
479 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
480 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
481 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
484 //______________________________________________________________________________
485 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
487 // Adds a user task to the global list of tasks.
488 task->SetActive(kFALSE);
492 //______________________________________________________________________________
493 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
495 // Retreive task by name.
496 if (!fTasks) return NULL;
497 return (AliAnalysisTask*)fTasks->FindObject(name);
500 //______________________________________________________________________________
501 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
502 TClass *datatype, EAliAnalysisContType type, const char *filename)
504 // Create a data container of a certain type. Types can be:
505 // kExchangeContainer = 0, used to exchange date between tasks
506 // kInputContainer = 1, used to store input data
507 // kOutputContainer = 2, used for posting results
508 if (fContainers->FindObject(name)) {
509 Error("CreateContainer","A container named %s already defined !\n",name);
512 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
513 fContainers->Add(cont);
515 case kInputContainer:
518 case kOutputContainer:
520 if (filename && strlen(filename)) cont->SetFileName(filename);
522 case kExchangeContainer:
528 //______________________________________________________________________________
529 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
530 AliAnalysisDataContainer *cont)
532 // Connect input of an existing task to a data container.
533 if (!fTasks->FindObject(task)) {
535 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
537 Bool_t connected = task->ConnectInput(islot, cont);
541 //______________________________________________________________________________
542 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
543 AliAnalysisDataContainer *cont)
545 // Connect output of an existing task to a data container.
546 if (!fTasks->FindObject(task)) {
548 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
550 Bool_t connected = task->ConnectOutput(islot, cont);
554 //______________________________________________________________________________
555 void AliAnalysisManager::CleanContainers()
557 // Clean data from all containers that have already finished all client tasks.
558 TIter next(fContainers);
559 AliAnalysisDataContainer *cont;
560 while ((cont=(AliAnalysisDataContainer *)next())) {
561 if (cont->IsOwnedData() &&
562 cont->IsDataReady() &&
563 cont->ClientsExecuted()) cont->DeleteData();
567 //______________________________________________________________________________
568 Bool_t AliAnalysisManager::InitAnalysis()
570 // Initialization of analysis chain of tasks. Should be called after all tasks
571 // and data containers are properly connected
572 // Check for input/output containers
574 // Check for top tasks (depending only on input data containers)
575 if (!fTasks->First()) {
576 Error("InitAnalysis", "Analysis has no tasks !");
580 AliAnalysisTask *task;
581 AliAnalysisDataContainer *cont;
584 Bool_t iszombie = kFALSE;
585 Bool_t istop = kTRUE;
587 while ((task=(AliAnalysisTask*)next())) {
590 Int_t ninputs = task->GetNinputs();
591 for (i=0; i<ninputs; i++) {
592 cont = task->GetInputSlot(i)->GetContainer();
600 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
603 if (iszombie) continue;
604 // Check if cont is an input container
605 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
606 // Connect to parent task
610 fTopTasks->Add(task);
614 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
617 // Check now if there are orphan tasks
618 for (i=0; i<ntop; i++) {
619 task = (AliAnalysisTask*)fTopTasks->At(i);
624 while ((task=(AliAnalysisTask*)next())) {
625 if (!task->IsUsed()) {
627 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
630 // Check the task hierarchy (no parent task should depend on data provided
631 // by a daughter task)
632 for (i=0; i<ntop; i++) {
633 task = (AliAnalysisTask*)fTopTasks->At(i);
634 if (task->CheckCircularDeps()) {
635 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
640 // Check that all containers feeding post-event loop tasks are in the outputs list
641 TIter nextcont(fContainers); // loop over all containers
642 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
643 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
644 if (cont->HasConsumers()) {
645 // Check if one of the consumers is post event loop
646 TIter nextconsumer(cont->GetConsumers());
647 while ((task=(AliAnalysisTask*)nextconsumer())) {
648 if (task->IsPostEventLoop()) {
660 //______________________________________________________________________________
661 void AliAnalysisManager::PrintStatus(Option_t *option) const
663 // Print task hierarchy.
664 TIter next(fTopTasks);
665 AliAnalysisTask *task;
666 while ((task=(AliAnalysisTask*)next()))
667 task->PrintTask(option);
670 //______________________________________________________________________________
671 void AliAnalysisManager::ResetAnalysis()
673 // Reset all execution flags and clean containers.
677 //______________________________________________________________________________
678 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
680 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
682 Error("StartAnalysis","Analysis manager was not initialized !");
686 cout << "StartAnalysis: " << GetName() << endl;
688 TString anaType = type;
690 fMode = kLocalAnalysis;
692 if (anaType.Contains("proof")) fMode = kProofAnalysis;
693 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
695 if (fMode == kGridAnalysis) {
696 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
697 fMode = kLocalAnalysis;
700 SetEventLoop(kFALSE);
701 // Disable all branches if requested and set event loop mode
703 if (TestBit(kDisableBranches)) {
704 printf("Disabling all branches...\n");
705 // tree->SetBranchStatus("*",0); // not yet working
710 TChain *chain = dynamic_cast<TChain*>(tree);
712 // Initialize locally all tasks
714 AliAnalysisTask *task;
715 while ((task=(AliAnalysisTask*)next())) {
723 AliAnalysisTask *task;
724 // Call CreateOutputObjects for all tasks
725 while ((task=(AliAnalysisTask*)next())) {
726 TDirectory *curdir = gDirectory;
727 task->CreateOutputObjects();
728 if (curdir) curdir->cd();
734 // Run tree-based analysis via AliAnalysisSelector
735 // gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
736 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
737 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
738 gROOT->ProcessLine(line);
739 sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
740 gROOT->ProcessLine(line);
743 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
744 printf("StartAnalysis: no PROOF!!!\n");
747 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
748 gROOT->ProcessLine(line);
751 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
752 chain->Process("AliAnalysisSelector");
754 printf("StartAnalysis: no chain\n");
759 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
763 //______________________________________________________________________________
764 void AliAnalysisManager::ExecAnalysis(Option_t *option)
768 Error("ExecAnalysis", "Analysis manager was not initialized !");
771 AliAnalysisTask *task;
772 // Check if the top tree is active.
775 // De-activate all tasks
776 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
777 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
779 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
782 cont->SetData(fTree); // This will notify all consumers
784 // Call BeginEvent() for optional input/output and MC services
785 if (fInputEventHandler) fInputEventHandler ->BeginEvent();
786 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent();
787 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
790 TIter next1(cont->GetConsumers());
791 while ((task=(AliAnalysisTask*)next1())) {
793 cout << " Executing task " << task->GetName() << endl;
796 task->ExecuteTask(option);
799 // Call FinishEvent() for optional output and MC services
800 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
801 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
802 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
806 // The event loop is not controlled by TSelector
808 // Call BeginEvent() for optional input/output and MC services
809 if (fInputEventHandler) fInputEventHandler ->BeginEvent();
810 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent();
811 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
812 TIter next2(fTopTasks);
813 while ((task=(AliAnalysisTask*)next2())) {
814 task->SetActive(kTRUE);
816 cout << " Executing task " << task->GetName() << endl;
818 task->ExecuteTask(option);
821 // Call FinishEvent() for optional output and MC services
822 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
823 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
824 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
827 //______________________________________________________________________________
828 void AliAnalysisManager::FinishAnalysis()