]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
new QA plot (nunmber of clusters/track/species) by AlexW
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisManager.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17 // Author: Andrei Gheata, 31/05/2006
18
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 //==============================================================================
25 //
26 //==============================================================================
27
28 #include <Riostream.h>
29
30 #include <TClass.h>
31 #include <TFile.h>
32 #include <TKey.h>
33 #include <TMethodCall.h>
34 #include <TChain.h>
35 #include <TSystem.h>
36 #include <TROOT.h>
37 #include <TCanvas.h>
38
39 #include "AliAnalysisSelector.h"
40 #include "AliAnalysisTask.h"
41 #include "AliAnalysisDataContainer.h"
42 #include "AliAnalysisDataSlot.h"
43 #include "AliVEventHandler.h"
44 #include "AliVEventPool.h"
45 #include "AliSysInfo.h"
46 #include "AliAnalysisManager.h"
47
48 ClassImp(AliAnalysisManager)
49
50 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
51
52 //______________________________________________________________________________
53 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
54                    :TNamed(name,title),
55                     fTree(NULL),
56                     fInputEventHandler(NULL),
57                     fOutputEventHandler(NULL),
58                     fMCtruthEventHandler(NULL),
59                     fEventPool(NULL),
60                     fCurrentEntry(-1),
61                     fNSysInfo(0),
62                     fMode(kLocalAnalysis),
63                     fInitOK(kFALSE),
64                     fDebug(0),
65                     fSpecialOutputLocation(""), 
66                     fTasks(NULL),
67                     fTopTasks(NULL),
68                     fZombies(NULL),
69                     fContainers(NULL),
70                     fInputs(NULL),
71                     fOutputs(NULL),
72                     fSelector(NULL)
73 {
74 // Default constructor.
75    fgAnalysisManager = this;
76    fTasks      = new TObjArray();
77    fTopTasks   = new TObjArray();
78    fZombies    = new TObjArray();
79    fContainers = new TObjArray();
80    fInputs     = new TObjArray();
81    fOutputs    = new TObjArray();
82    SetEventLoop(kTRUE);
83 }
84
85 //______________________________________________________________________________
86 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
87                    :TNamed(other),
88                     fTree(NULL),
89                     fInputEventHandler(NULL),
90                     fOutputEventHandler(NULL),
91                     fMCtruthEventHandler(NULL),
92                     fEventPool(NULL),
93                     fCurrentEntry(-1),
94                     fNSysInfo(0),
95                     fMode(other.fMode),
96                     fInitOK(other.fInitOK),
97                     fDebug(other.fDebug),
98                     fSpecialOutputLocation(""), 
99                     fTasks(NULL),
100                     fTopTasks(NULL),
101                     fZombies(NULL),
102                     fContainers(NULL),
103                     fInputs(NULL),
104                     fOutputs(NULL),
105                     fSelector(NULL)
106 {
107 // Copy constructor.
108    fTasks      = new TObjArray(*other.fTasks);
109    fTopTasks   = new TObjArray(*other.fTopTasks);
110    fZombies    = new TObjArray(*other.fZombies);
111    fContainers = new TObjArray(*other.fContainers);
112    fInputs     = new TObjArray(*other.fInputs);
113    fOutputs    = new TObjArray(*other.fOutputs);
114    fgAnalysisManager = this;
115 }
116    
117 //______________________________________________________________________________
118 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
119 {
120 // Assignment
121    if (&other != this) {
122       TNamed::operator=(other);
123       fInputEventHandler   = other.fInputEventHandler;
124       fOutputEventHandler  = other.fOutputEventHandler;
125       fMCtruthEventHandler = other.fMCtruthEventHandler;
126       fEventPool           = other.fEventPool;
127       fTree       = NULL;
128       fCurrentEntry = -1;
129       fNSysInfo   = other.fNSysInfo;
130       fMode       = other.fMode;
131       fInitOK     = other.fInitOK;
132       fDebug      = other.fDebug;
133       fTasks      = new TObjArray(*other.fTasks);
134       fTopTasks   = new TObjArray(*other.fTopTasks);
135       fZombies    = new TObjArray(*other.fZombies);
136       fContainers = new TObjArray(*other.fContainers);
137       fInputs     = new TObjArray(*other.fInputs);
138       fOutputs    = new TObjArray(*other.fOutputs);
139       fSelector   = NULL;
140       fgAnalysisManager = this;
141    }
142    return *this;
143 }
144
145 //______________________________________________________________________________
146 AliAnalysisManager::~AliAnalysisManager()
147 {
148 // Destructor.
149    if (fTasks) {fTasks->Delete(); delete fTasks;}
150    if (fTopTasks) delete fTopTasks;
151    if (fZombies) delete fZombies;
152    if (fContainers) {fContainers->Delete(); delete fContainers;}
153    if (fInputs) delete fInputs;
154    if (fOutputs) delete fOutputs;
155    if (fgAnalysisManager==this) fgAnalysisManager = NULL;
156 }
157
158 //______________________________________________________________________________
159 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
160 {
161 // Read one entry of the tree or a whole branch.
162    if (fDebug > 0) printf("== AliAnalysisManager::GetEntry(%lld)\n", entry);
163    fCurrentEntry = entry;
164    return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
165 }
166    
167 //______________________________________________________________________________
168 void AliAnalysisManager::Init(TTree *tree)
169 {
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.
175    if (!tree) return;
176    if (fDebug > 0) {
177       printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
178    }
179
180    // Call InitTree of EventHandler
181    if (fOutputEventHandler) {
182       if (fMode == kProofAnalysis) {
183          fOutputEventHandler->Init(0x0, "proof");
184       } else {
185          fOutputEventHandler->Init(0x0, "local");
186       }
187    }
188
189    if (fInputEventHandler) {
190       if (fMode == kProofAnalysis) {
191          fInputEventHandler->Init(tree, "proof");
192       } else {
193          fInputEventHandler->Init(tree, "local");
194       }
195    } else {
196       // If no input event handler we need to get the tree once
197       // for the chain
198       if(!tree->GetTree()) tree->LoadTree(0);
199    }
200    
201
202    if (fMCtruthEventHandler) {
203       if (fMode == kProofAnalysis) {
204          fMCtruthEventHandler->Init(0x0, "proof");
205       } else {
206          fMCtruthEventHandler->Init(0x0, "local");
207       }
208    }
209
210    if (!fInitOK) InitAnalysis();
211    if (!fInitOK) return;
212    fTree = tree;
213    AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
214    if (!top) {
215       Error("Init","No top input container !");
216       return;
217    }
218    top->SetData(tree);
219    if (fDebug > 0) {
220       printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
221    }
222 }
223
224 //______________________________________________________________________________
225 void AliAnalysisManager::SlaveBegin(TTree *tree)
226 {
227   // The SlaveBegin() function is called after the Begin() function.
228   // When running with PROOF SlaveBegin() is called on each slave server.
229   // The tree argument is deprecated (on PROOF 0 is passed).
230    if (fDebug > 0) printf("->AliAnalysisManager::SlaveBegin()\n");
231    static Bool_t isCalled = kFALSE;
232    TDirectory *curdir = gDirectory;
233    // Call SlaveBegin only once in case of mixing
234    if (isCalled && fMode==kMixingAnalysis) return;
235    // Call Init of EventHandler
236    if (fOutputEventHandler) {
237       if (fMode == kProofAnalysis) {
238          TIter nextout(fOutputs);
239          AliAnalysisDataContainer *c_aod;
240          while ((c_aod=(AliAnalysisDataContainer*)nextout())) if (!strcmp(c_aod->GetFileName(),"default")) break;
241          if (c_aod && c_aod->IsSpecialOutput()) {
242             // Merging via files
243             if (fDebug > 1) printf("   Initializing special output file %s...\n", fOutputEventHandler->GetOutputFileName());
244             OpenProofFile(fOutputEventHandler->GetOutputFileName(), "RECREATE");
245             c_aod->SetFile(gFile);
246             fOutputEventHandler->Init("proofspecial");
247          } else {
248             // Merging in memory
249             fOutputEventHandler->Init("proof");
250          }   
251       } else {
252          fOutputEventHandler->Init("local");
253       }
254    }
255
256    if (fInputEventHandler) {
257       fInputEventHandler->SetInputTree(tree);
258       if (fMode == kProofAnalysis) {
259          fInputEventHandler->Init("proof");
260       } else {
261          fInputEventHandler->Init("local");
262       }
263    }
264
265    if (fMCtruthEventHandler) {
266       if (fMode == kProofAnalysis) {
267          fMCtruthEventHandler->Init("proof");
268       } else {
269          fMCtruthEventHandler->Init("local");
270       }
271    }
272    if (curdir) curdir->cd();
273    
274    TIter next(fTasks);
275    AliAnalysisTask *task;
276    // Call CreateOutputObjects for all tasks
277    while ((task=(AliAnalysisTask*)next())) {
278       curdir = gDirectory;
279       task->CreateOutputObjects();
280       if (curdir) curdir->cd();
281    }
282    isCalled = kTRUE;
283
284    if (fDebug > 0) printf("<-AliAnalysisManager::SlaveBegin()\n");
285 }
286
287 //______________________________________________________________________________
288 Bool_t AliAnalysisManager::Notify()
289 {
290    // The Notify() function is called when a new file is opened. This
291    // can be either for a new TTree in a TChain or when when a new TTree
292    // is started when using PROOF. It is normaly not necessary to make changes
293    // to the generated code, but the routine can be extended by the
294    // user if needed. The return value is currently not used.
295    if (!fTree) {
296       Error("Notify","No current tree.");
297       return kFALSE;
298    }   
299    TFile *curfile = fTree->GetCurrentFile();
300    if (!curfile) {
301       Error("Notify","No current file");
302       return kFALSE;
303    }   
304    
305    if (fDebug > 0) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
306    TIter next(fTasks);
307    AliAnalysisTask *task;
308    // Call Notify for all tasks
309    while ((task=(AliAnalysisTask*)next())) 
310       task->Notify();
311         
312    // Call Notify of the event handlers
313    if (fInputEventHandler) {
314        fInputEventHandler->Notify(curfile->GetName());
315    }
316
317    if (fOutputEventHandler) {
318        fOutputEventHandler->Notify(curfile->GetName());
319    }
320
321    if (fMCtruthEventHandler) {
322        fMCtruthEventHandler->Notify(curfile->GetName());
323    }
324    if (fDebug > 0) printf("<-AliAnalysisManager::Notify()\n");
325    return kTRUE;
326 }    
327
328 //______________________________________________________________________________
329 Bool_t AliAnalysisManager::Process(Long64_t entry)
330 {
331   // The Process() function is called for each entry in the tree (or possibly
332   // keyed object in the case of PROOF) to be processed. The entry argument
333   // specifies which entry in the currently loaded tree is to be processed.
334   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
335   // to read either all or the required parts of the data. When processing
336   // keyed objects with PROOF, the object is already loaded and is available
337   // via the fObject pointer.
338   //
339   // This function should contain the "body" of the analysis. It can contain
340   // simple or elaborate selection criteria, run algorithms on the data
341   // of the event and typically fill histograms.
342
343   // WARNING when a selector is used with a TChain, you must use
344   //  the pointer to the current TTree to call GetEntry(entry).
345   //  The entry is always the local entry number in the current tree.
346   //  Assuming that fChain is the pointer to the TChain being processed,
347   //  use fChain->GetTree()->GetEntry(entry).
348    if (fDebug > 0) printf("->AliAnalysisManager::Process(%lld)\n", entry);
349
350    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
351    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
352    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
353    
354    GetEntry(entry);
355    ExecAnalysis();
356    if (fDebug > 0) printf("<-AliAnalysisManager::Process()\n");
357    return kTRUE;
358 }
359
360 //______________________________________________________________________________
361 void AliAnalysisManager::PackOutput(TList *target)
362 {
363   // Pack all output data containers in the output list. Called at SlaveTerminate
364   // stage in PROOF case for each slave.
365    if (fDebug > 0) printf("->AliAnalysisManager::PackOutput()\n");
366    if (!target) {
367       Error("PackOutput", "No target. Aborting.");
368       return;
369    }
370    if (fInputEventHandler)   fInputEventHandler  ->Terminate();
371    if (fOutputEventHandler)  fOutputEventHandler ->Terminate();
372    if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
373
374    // Call FinishTaskOutput() for each event loop task (not called for 
375    // post-event loop tasks - use Terminate() fo those)
376    TIter nexttask(fTasks);
377    AliAnalysisTask *task;
378    while ((task=(AliAnalysisTask*)nexttask())) {
379       if (!task->IsPostEventLoop()) {
380          if (fDebug > 0) printf("->FinishTaskOutput: task %s\n", task->GetName());
381          task->FinishTaskOutput();
382          if (fDebug > 0) printf("<-FinishTaskOutput: task %s\n", task->GetName());
383       }
384    }      
385    
386    if (fMode == kProofAnalysis) {
387       TIter next(fOutputs);
388       AliAnalysisDataContainer *output;
389       Bool_t isManagedByHandler = kFALSE;
390       while ((output=(AliAnalysisDataContainer*)next())) {
391          // Do not consider outputs of post event loop tasks
392          isManagedByHandler = kFALSE;
393          if (output->GetProducer()->IsPostEventLoop()) continue;
394          const char *filename = output->GetFileName();
395          if (!(strcmp(filename, "default")) && fOutputEventHandler) {
396             isManagedByHandler = kTRUE;
397             filename = fOutputEventHandler->GetOutputFileName();
398          }
399          // Check if data was posted to this container. If not, issue an error.
400          if (!output->GetData() && !isManagedByHandler) {
401             Error("PackOutput", "No data for output container %s. Forgot to PostData ?\n", output->GetName());
402             continue;
403          }   
404          if (!output->IsSpecialOutput()) {
405             // Normal outputs
406             if (strlen(filename) && !isManagedByHandler) {
407                // File resident outputs
408                TFile *file = output->GetFile();
409                // Backup current folder
410                TDirectory *opwd = gDirectory;
411                // Create file if not existing and register to container.
412                if (file) file->cd();
413                else      file = new TFile(filename, "RECREATE"); 
414                if (file->IsZombie()) {
415                   Fatal("PackOutput", "Could not recreate file %s\n", filename);
416                   return;
417                }   
418                output->SetFile(file);
419                // Clear file list to release object ownership to user.
420                file->Clear();
421                // Save data to file, then close.
422                if (output->GetData()->InheritsFrom(TCollection::Class())) {
423                   // If data is a collection, we set the name of the collection 
424                   // as the one of the container and we save as a single key.
425                   TCollection *coll = (TCollection*)output->GetData();
426                   coll->SetName(output->GetName());
427                   coll->Write(output->GetName(), TObject::kSingleKey);
428                } else {
429                   if (output->GetData()->InheritsFrom(TTree::Class())) {
430                      TTree *tree = (TTree*)output->GetData();
431                      tree->SetDirectory(file);
432                      tree->AutoSave();
433                   } else {
434                      output->GetData()->Write();
435                   }   
436                }      
437                if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
438                if (fDebug > 2) {
439                   printf("   file %s listing content:\n", filename);
440                   file->ls();
441                }   
442                file->Close();
443                // Restore current directory
444                if (opwd) opwd->cd();
445             } else {
446                // Memory-resident outputs   
447                if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
448             }   
449             AliAnalysisDataWrapper *wrap = 0;
450             if (isManagedByHandler) {
451                wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
452                wrap->SetName(output->GetName());
453             }   
454             else                    wrap =output->ExportData();
455             // Output wrappers must NOT delete data after merging - the user owns them
456             wrap->SetDeleteData(kFALSE);
457             target->Add(wrap);
458          } else {
459          // Special outputs
460             TDirectory *opwd = gDirectory;
461             TFile *file = output->GetFile();
462             if (fDebug > 1 && file) printf("PackOutput %s: file merge, special output\n", output->GetName());
463             if (isManagedByHandler) {
464                // Terminate IO for files managed by the output handler
465                if (file) file->Write();
466                if (file && fDebug > 2) {
467                   printf("   handled file %s listing content:\n", file->GetName());
468                   file->ls();
469                }   
470                fOutputEventHandler->TerminateIO();
471                continue;
472             }   
473             
474             if (!file) {
475                AliAnalysisTask *producer = output->GetProducer();
476                Error("PackOutput", 
477                      "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
478                      output->GetFileName(), output->GetName(), producer->ClassName());
479                continue;
480             }   
481             file->cd();
482             // Release object ownership to users after writing data to file
483             if (output->GetData()->InheritsFrom(TCollection::Class())) {
484                // If data is a collection, we set the name of the collection 
485                // as the one of the container and we save as a single key.
486                TCollection *coll = (TCollection*)output->GetData();
487                coll->SetName(output->GetName());
488                coll->Write(output->GetName(), TObject::kSingleKey);
489             } else {
490                if (output->GetData()->InheritsFrom(TTree::Class())) {
491                   TTree *tree = (TTree*)output->GetData();
492                   tree->SetDirectory(file);
493                   tree->AutoSave();
494                } else {
495                   output->GetData()->Write();
496                }   
497             }      
498             file->Clear();
499             if (fDebug > 2) {
500                printf("   file %s listing content:\n", output->GetFileName());
501                file->ls();
502             }
503             TString outFilename = file->GetName();
504             file->Close();
505             // Restore current directory
506             if (opwd) opwd->cd();
507             // Check if a special output location was provided or the output files have to be merged
508             if (strlen(fSpecialOutputLocation.Data())) {
509                TString remote = fSpecialOutputLocation;
510                remote += "/";
511                Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
512                remote += Form("%s_%d_", gSystem->HostName(), gid);
513                remote += output->GetFileName();
514                TFile::Cp ( outFilename.Data(), remote.Data() );
515             } else {
516             // No special location specified-> use TProofOutputFile as merging utility
517             // The file at this output slot must be opened in CreateOutputObjects
518                if (fDebug > 1) printf("   File %s to be merged...\n", output->GetFileName());
519             }
520          }      
521       }
522    } 
523    if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
524 }
525
526 //______________________________________________________________________________
527 void AliAnalysisManager::ImportWrappers(TList *source)
528 {
529 // Import data in output containers from wrappers coming in source.
530    if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
531    TIter next(fOutputs);
532    AliAnalysisDataContainer *cont;
533    AliAnalysisDataWrapper   *wrap;
534    Int_t icont = 0;
535    while ((cont=(AliAnalysisDataContainer*)next())) {
536       wrap = 0;
537       if (cont->GetProducer()->IsPostEventLoop()) continue;
538       const char *filename = cont->GetFileName();
539       Bool_t isManagedByHandler = kFALSE;
540       if (!(strcmp(filename, "default")) && fOutputEventHandler) {
541          isManagedByHandler = kTRUE;
542          filename = fOutputEventHandler->GetOutputFileName();
543       }
544       if (cont->IsSpecialOutput()) {
545          if (strlen(fSpecialOutputLocation.Data()) && !isManagedByHandler) continue;
546          // Copy merged file from PROOF scratch space
547          char full_path[512];
548          char ch_url[512];
549          TObject *pof =  source->FindObject(filename);
550          if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
551             Error("ImportWrappers", "TProofOutputFile object not found in output list for container %s", cont->GetName());
552             continue;
553          }
554          gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", ((TProofOutputFile*)0x%lx)->GetOutputFileName();)", full_path, pof));
555          gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", gProof->GetUrl();)", ch_url));
556          TString clientUrl(ch_url);
557          TString full_path_str(full_path);
558          if (clientUrl.Contains("localhost")){
559             TObjArray* array = full_path_str.Tokenize ( "//" );
560             TObjString *strobj = ( TObjString *)array->At(1);
561             full_path_str.ReplaceAll(strobj->GetString().Data(),"localhost:11094");
562             if (fDebug > 1) Info("ImportWrappers","Using tunnel from %s to %s",full_path_str.Data(),filename);
563          }
564          if (fDebug > 1) 
565             printf("   Copying file %s from PROOF scratch space\n", full_path_str.Data());
566          Bool_t gotit = TFile::Cp(full_path_str.Data(), filename); 
567          if (!gotit) {
568             Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
569             continue;
570          }
571          // Normally we should connect data from the copied file to the
572          // corresponding output container, but it is not obvious how to do this
573          // automatically if several objects in file...
574          TFile *f = TFile::Open(filename, "READ");
575          if (!f) {
576             Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
577             continue;
578          }   
579          TObject *obj = 0;
580          // Try to fetch first a list object having the container name.
581          obj = f->Get(cont->GetName());
582          if (!obj) {
583          // Fetch first object from file having the container type.
584             TIter nextkey(f->GetListOfKeys());
585             TKey *key;
586             while ((key=(TKey*)nextkey())) {
587                obj = f->Get(key->GetName());
588                if (obj && obj->IsA()->InheritsFrom(cont->GetType())) break;
589             }                     
590          }
591          if (!obj) {
592             Error("ImportWrappers", "Could not find object for container %s in file %s", cont->GetName(), filename);
593             continue;
594          }  
595          wrap = new AliAnalysisDataWrapper(obj);
596          wrap->SetDeleteData(kFALSE);
597       }   
598       if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
599       if (!wrap) {
600          Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
601          continue;
602       }
603       icont++;
604       if (fDebug > 1) {
605          printf("   Importing data for container %s", cont->GetName());
606          if (strlen(filename)) printf("    -> file %s\n", filename);
607          else printf("\n");
608       }   
609       cont->ImportData(wrap);
610    }         
611    if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
612 }
613
614 //______________________________________________________________________________
615 void AliAnalysisManager::UnpackOutput(TList *source)
616 {
617   // Called by AliAnalysisSelector::Terminate only on the client.
618    if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
619    if (!source) {
620       Error("UnpackOutput", "No target. Aborting.");
621       return;
622    }
623    if (fDebug > 1) printf("   Source list contains %d containers\n", source->GetSize());
624
625    if (fMode == kProofAnalysis) ImportWrappers(source);
626
627    TIter next(fOutputs);
628    AliAnalysisDataContainer *output;
629    while ((output=(AliAnalysisDataContainer*)next())) {
630       if (!output->GetData()) continue;
631       // Check if there are client tasks that run post event loop
632       if (output->HasConsumers()) {
633          // Disable event loop semaphore
634          output->SetPostEventLoop(kTRUE);
635          TObjArray *list = output->GetConsumers();
636          Int_t ncons = list->GetEntriesFast();
637          for (Int_t i=0; i<ncons; i++) {
638             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
639             task->CheckNotify(kTRUE);
640             // If task is active, execute it
641             if (task->IsPostEventLoop() && task->IsActive()) {
642                if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
643                task->ExecuteTask();
644             }   
645          }
646       }   
647    }
648    if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
649 }
650
651 //______________________________________________________________________________
652 void AliAnalysisManager::Terminate()
653 {
654   // The Terminate() function is the last function to be called during
655   // a query. It always runs on the client, it can be used to present
656   // the results graphically.
657    if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
658    AliAnalysisTask *task;
659    TIter next(fTasks);
660    // Call Terminate() for tasks
661    while ((task=(AliAnalysisTask*)next())) task->Terminate();
662    //
663    TIter next1(fOutputs);
664    AliAnalysisDataContainer *output;
665    while ((output=(AliAnalysisDataContainer*)next1())) {
666       // Special outputs have the files already closed and written.
667      if (output->IsSpecialOutput()&&(fMode == kProofAnalysis)) continue;
668       const char *filename = output->GetFileName();
669       if (!(strcmp(filename, "default"))) {
670          if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
671          TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
672          if (aodfile) {
673             if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
674             aodfile->Write();
675             continue;
676          }   
677       }      
678       if (!strlen(filename)) continue;
679       if (!output->GetData()) continue;
680       TFile *file = output->GetFile();
681       TDirectory *opwd = gDirectory;
682       file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
683       if (!file) file = new TFile(filename, "RECREATE");
684       if (file->IsZombie()) continue;
685       output->SetFile(file);
686       file->cd();
687       if (fDebug > 1) printf("   writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
688       if (output->GetData()->InheritsFrom(TCollection::Class())) {
689       // If data is a collection, we set the name of the collection 
690       // as the one of the container and we save as a single key.
691          TCollection *coll = (TCollection*)output->GetData();
692          coll->SetName(output->GetName());
693          coll->Write(output->GetName(), TObject::kSingleKey);
694       } else {
695          if (output->GetData()->InheritsFrom(TTree::Class())) {
696             TTree *tree = (TTree*)output->GetData();
697             tree->SetDirectory(file);
698             tree->AutoSave();
699          } else {
700             output->GetData()->Write();
701          }   
702       }      
703       if (opwd) opwd->cd();
704    }   
705    next1.Reset();
706    while ((output=(AliAnalysisDataContainer*)next1())) {
707       // Close all files at output
708       TDirectory *opwd = gDirectory;
709       if (output->GetFile()) output->GetFile()->Close();
710       if (opwd) opwd->cd();
711    }   
712
713    if (fInputEventHandler)   fInputEventHandler  ->TerminateIO();
714    if (fOutputEventHandler)  fOutputEventHandler ->TerminateIO();
715    if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
716
717    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
718    if (getsysInfo) {
719       TDirectory *cdir = gDirectory;
720       TFile f("syswatch.root", "RECREATE");
721       if (!f.IsZombie()) {
722          TTree *tree = AliSysInfo::MakeTree("syswatch.log");
723          tree->SetMarkerStyle(kCircle);
724          tree->SetMarkerColor(kBlue);
725          tree->SetMarkerSize(0.5);
726          if (!gROOT->IsBatch()) {
727             tree->SetAlias("event", "id0");
728             tree->SetAlias("memUSED", "pI.fMemVirtual");
729             tree->SetAlias("userCPU", "pI.fCpuUser");
730             TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
731             c->Divide(2,1,0.01,0.01);
732             c->cd(1);
733             tree->Draw("memUSED:event","","", 1234567890, 0);
734             c->cd(2);
735             tree->Draw("userCPU:event","","", 1234567890, 0);
736          }   
737          tree->Write();
738          f.Close();
739          delete tree;
740       }
741       if (cdir) cdir->cd();
742    }      
743    if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
744 }
745
746 //______________________________________________________________________________
747 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
748 {
749 // Adds a user task to the global list of tasks.
750    if (fTasks->FindObject(task)) {
751       Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
752       return;
753    }   
754    task->SetActive(kFALSE);
755    fTasks->Add(task);
756 }  
757
758 //______________________________________________________________________________
759 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
760 {
761 // Retreive task by name.
762    if (!fTasks) return NULL;
763    return (AliAnalysisTask*)fTasks->FindObject(name);
764 }
765
766 //______________________________________________________________________________
767 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
768                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
769 {
770 // Create a data container of a certain type. Types can be:
771 //   kExchangeContainer  = 0, used to exchange date between tasks
772 //   kInputContainer   = 1, used to store input data
773 //   kOutputContainer  = 2, used for posting results
774    if (fContainers->FindObject(name)) {
775       Error("CreateContainer","A container named %s already defined !\n",name);
776       return NULL;
777    }   
778    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
779    fContainers->Add(cont);
780    switch (type) {
781       case kInputContainer:
782          fInputs->Add(cont);
783          break;
784       case kOutputContainer:
785          fOutputs->Add(cont);
786          if (filename && strlen(filename)) {
787             cont->SetFileName(filename);
788             cont->SetDataOwned(kFALSE);  // data owned by the file
789          }   
790          break;
791       case kExchangeContainer:
792          break;   
793    }
794    return cont;
795 }
796          
797 //______________________________________________________________________________
798 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
799                                         AliAnalysisDataContainer *cont)
800 {
801 // Connect input of an existing task to a data container.
802    if (!fTasks->FindObject(task)) {
803       AddTask(task);
804       Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
805    } 
806    Bool_t connected = task->ConnectInput(islot, cont);
807    return connected;
808 }   
809
810 //______________________________________________________________________________
811 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
812                                         AliAnalysisDataContainer *cont)
813 {
814 // Connect output of an existing task to a data container.
815    if (!fTasks->FindObject(task)) {
816       AddTask(task);
817       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
818    } 
819    Bool_t connected = task->ConnectOutput(islot, cont);
820    return connected;
821 }   
822                                
823 //______________________________________________________________________________
824 void AliAnalysisManager::CleanContainers()
825 {
826 // Clean data from all containers that have already finished all client tasks.
827    TIter next(fContainers);
828    AliAnalysisDataContainer *cont;
829    while ((cont=(AliAnalysisDataContainer *)next())) {
830       if (cont->IsOwnedData() && 
831           cont->IsDataReady() && 
832           cont->ClientsExecuted()) cont->DeleteData();
833    }
834 }
835
836 //______________________________________________________________________________
837 Bool_t AliAnalysisManager::InitAnalysis()
838 {
839 // Initialization of analysis chain of tasks. Should be called after all tasks
840 // and data containers are properly connected
841    // Check for input/output containers
842    fInitOK = kFALSE;
843    // Check for top tasks (depending only on input data containers)
844    if (!fTasks->First()) {
845       Error("InitAnalysis", "Analysis has no tasks !");
846       return kFALSE;
847    }   
848    TIter next(fTasks);
849    AliAnalysisTask *task;
850    AliAnalysisDataContainer *cont;
851    Int_t ntop = 0;
852    Int_t nzombies = 0;
853    Bool_t iszombie = kFALSE;
854    Bool_t istop = kTRUE;
855    Int_t i;
856    while ((task=(AliAnalysisTask*)next())) {
857       istop = kTRUE;
858       iszombie = kFALSE;
859       Int_t ninputs = task->GetNinputs();
860       for (i=0; i<ninputs; i++) {
861          cont = task->GetInputSlot(i)->GetContainer();
862          if (!cont) {
863             if (!iszombie) {
864                task->SetZombie();
865                fZombies->Add(task);
866                nzombies++;
867                iszombie = kTRUE;
868             }   
869             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
870                   i, task->GetName()); 
871          }
872          if (iszombie) continue;
873          // Check if cont is an input container
874          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
875          // Connect to parent task
876       }
877       if (istop) {
878          ntop++;
879          fTopTasks->Add(task);
880       }
881    }
882    if (!ntop) {
883       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
884       return kFALSE;
885    }                        
886    // Check now if there are orphan tasks
887    for (i=0; i<ntop; i++) {
888       task = (AliAnalysisTask*)fTopTasks->At(i);
889       task->SetUsed();
890    }
891    Int_t norphans = 0;
892    next.Reset();
893    while ((task=(AliAnalysisTask*)next())) {
894       if (!task->IsUsed()) {
895          norphans++;
896          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
897       }   
898    }          
899    // Check the task hierarchy (no parent task should depend on data provided
900    // by a daughter task)
901    for (i=0; i<ntop; i++) {
902       task = (AliAnalysisTask*)fTopTasks->At(i);
903       if (task->CheckCircularDeps()) {
904          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
905          PrintStatus("dep");
906          return kFALSE;
907       }   
908    }
909    // Check that all containers feeding post-event loop tasks are in the outputs list
910    TIter nextcont(fContainers); // loop over all containers
911    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
912       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
913          if (cont->HasConsumers()) {
914          // Check if one of the consumers is post event loop
915             TIter nextconsumer(cont->GetConsumers());
916             while ((task=(AliAnalysisTask*)nextconsumer())) {
917                if (task->IsPostEventLoop()) {
918                   fOutputs->Add(cont);
919                   break;
920                }
921             }
922          }
923       }
924    }   
925    // Check if all special output containers have a file name provided
926    TIter nextout(fOutputs);
927    while ((cont=(AliAnalysisDataContainer*)nextout())) {
928       if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
929          Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
930          return kFALSE;
931       }
932    }      
933    fInitOK = kTRUE;
934    return kTRUE;
935 }   
936
937 //______________________________________________________________________________
938 void AliAnalysisManager::PrintStatus(Option_t *option) const
939 {
940 // Print task hierarchy.
941    if (!fInitOK) {
942       Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
943       return;
944    }   
945    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
946    if (getsysInfo)
947       Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
948    TIter next(fTopTasks);
949    AliAnalysisTask *task;
950    while ((task=(AliAnalysisTask*)next()))
951       task->PrintTask(option);
952 }
953
954 //______________________________________________________________________________
955 void AliAnalysisManager::ResetAnalysis()
956 {
957 // Reset all execution flags and clean containers.
958    CleanContainers();
959 }
960
961 //______________________________________________________________________________
962 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
963 {
964 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
965 // MIX. Process nentries starting from firstentry
966    if (!fInitOK) {
967       Error("StartAnalysis","Analysis manager was not initialized !");
968       return;
969    }
970    if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
971    TString anaType = type;
972    anaType.ToLower();
973    fMode = kLocalAnalysis;
974    if (anaType.Contains("proof"))     fMode = kProofAnalysis;
975    else if (anaType.Contains("grid")) fMode = kGridAnalysis;
976    else if (anaType.Contains("mix"))  fMode = kMixingAnalysis;
977
978    if (fMode == kGridAnalysis) {
979       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
980       fMode = kLocalAnalysis;
981    }
982    char line[256];
983    SetEventLoop(kFALSE);
984    // Enable event loop mode if a tree was provided
985    if (tree || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
986
987    TChain *chain = 0;
988    TString ttype = "TTree";
989    if (tree && tree->IsA() == TChain::Class()) {
990       chain = (TChain*)tree;
991       if (!chain || !chain->GetListOfFiles()->First()) {
992          Error("StartAnalysis", "Cannot process null or empty chain...");
993          return;
994       }   
995       ttype = "TChain";
996    }   
997
998    // Initialize locally all tasks (happens for all modes)
999    TIter next(fTasks);
1000    AliAnalysisTask *task;
1001    while ((task=(AliAnalysisTask*)next())) {
1002       task->LocalInit();
1003    }
1004    
1005    switch (fMode) {
1006       case kLocalAnalysis:
1007          if (!tree) {
1008             TIter nextT(fTasks);
1009             // Call CreateOutputObjects for all tasks
1010             while ((task=(AliAnalysisTask*)nextT())) {
1011                TDirectory *curdir = gDirectory;
1012                task->CreateOutputObjects();
1013                if (curdir) curdir->cd();
1014             }   
1015             ExecAnalysis();
1016             Terminate();
1017             return;
1018          } 
1019          // Run tree-based analysis via AliAnalysisSelector  
1020          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1021          fSelector = new AliAnalysisSelector(this);
1022          tree->Process(fSelector, "", nentries, firstentry);
1023          break;
1024       case kProofAnalysis:
1025          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1026             printf("StartAnalysis: no PROOF!!!\n");
1027             return;
1028          }   
1029          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1030          gROOT->ProcessLine(line);
1031          if (chain) {
1032             chain->SetProof();
1033             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1034             chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1035          } else {
1036             printf("StartAnalysis: no chain\n");
1037             return;
1038          }      
1039          break;
1040       case kGridAnalysis:
1041          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
1042          break;
1043       case kMixingAnalysis:   
1044          // Run event mixing analysis
1045          if (!fEventPool) {
1046             Error("StartAnalysis", "Cannot run event mixing without event pool");
1047             return;
1048          }
1049          cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1050          fSelector = new AliAnalysisSelector(this);
1051          while ((chain=fEventPool->GetNextChain())) {
1052             next.Reset();
1053             // Call NotifyBinChange for all tasks
1054             while ((task=(AliAnalysisTask*)next()))
1055                if (!task->IsPostEventLoop()) task->NotifyBinChange();
1056             chain->Process(fSelector);
1057          }
1058          PackOutput(fSelector->GetOutputList());
1059          Terminate();
1060    }   
1061 }   
1062
1063 //______________________________________________________________________________
1064 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1065 {
1066 // Start analysis for this manager on a given dataset. Analysis task can be: 
1067 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1068    if (!fInitOK) {
1069       Error("StartAnalysis","Analysis manager was not initialized !");
1070       return;
1071    }
1072    if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
1073    TString anaType = type;
1074    anaType.ToLower();
1075    if (!anaType.Contains("proof")) {
1076       Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1077       return;
1078    }   
1079    fMode = kProofAnalysis;
1080    char line[256];
1081    SetEventLoop(kTRUE);
1082    // Set the dataset flag
1083    TObject::SetBit(kUseDataSet);
1084    fTree = 0;
1085
1086    // Initialize locally all tasks
1087    TIter next(fTasks);
1088    AliAnalysisTask *task;
1089    while ((task=(AliAnalysisTask*)next())) {
1090       task->LocalInit();
1091    }
1092    
1093    if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1094       printf("StartAnalysis: no PROOF!!!\n");
1095       return;
1096    }   
1097    sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1098    gROOT->ProcessLine(line);
1099    sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
1100    if (!gROOT->ProcessLine(line)) {
1101       Error("StartAnalysis", "Dataset %s not found", dataset);
1102       return;
1103    }   
1104    sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1105            dataset, nentries, firstentry);
1106    cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1107    gROOT->ProcessLine(line);
1108 }   
1109
1110 //______________________________________________________________________________
1111 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
1112 {
1113 // Opens a special output file used in PROOF.
1114    char line[256];
1115    if (fMode!=kProofAnalysis || !fSelector) {
1116       Error("OpenProofFile","Cannot open PROOF file %s",filename);
1117       return NULL;
1118    }   
1119    sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
1120    if (fDebug > 1) printf("=== %s\n", line);
1121    gROOT->ProcessLine(line);
1122    sprintf(line, "pf->OpenFile(\"%s\");", option);
1123    gROOT->ProcessLine(line);
1124    if (fDebug > 1) {
1125       gROOT->ProcessLine("pf->Print()");
1126       printf(" == proof file name: %s\n", gFile->GetName());
1127    }   
1128    sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1129    if (fDebug > 1) printf("=== %s\n", line);
1130    gROOT->ProcessLine(line);
1131    return gFile;
1132 }   
1133
1134 //______________________________________________________________________________
1135 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1136 {
1137 // Execute analysis.
1138    static Long64_t ncalls = 0;
1139    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1140    if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1141    ncalls++;
1142    if (!fInitOK) {
1143      Error("ExecAnalysis", "Analysis manager was not initialized !");
1144       return;
1145    }   
1146    AliAnalysisTask *task;
1147    // Check if the top tree is active.
1148    if (fTree) {
1149       TIter next(fTasks);
1150    // De-activate all tasks
1151       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1152       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1153       if (!cont) {
1154               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1155          return;
1156       }   
1157       cont->SetData(fTree); // This will notify all consumers
1158       Long64_t entry = fTree->GetTree()->GetReadEntry();
1159       
1160 //
1161 //    Call BeginEvent() for optional input/output and MC services 
1162       if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
1163       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
1164       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1165 //
1166 //    Execute the tasks
1167 //      TIter next1(cont->GetConsumers());
1168       TIter next1(fTopTasks);
1169       while ((task=(AliAnalysisTask*)next1())) {
1170          if (fDebug >1) {
1171             cout << "    Executing task " << task->GetName() << endl;
1172          }   
1173          
1174          task->ExecuteTask(option);
1175       }
1176 //
1177 //    Call FinishEvent() for optional output and MC services 
1178       if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
1179       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
1180       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1181       // Gather system information if requested
1182       if (getsysInfo && ((ncalls%fNSysInfo)==0)) 
1183          AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1184       return;
1185    }   
1186    // The event loop is not controlled by TSelector   
1187 //
1188 //  Call BeginEvent() for optional input/output and MC services 
1189    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(-1);
1190    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(-1);
1191    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1192    TIter next2(fTopTasks);
1193    while ((task=(AliAnalysisTask*)next2())) {
1194       task->SetActive(kTRUE);
1195       if (fDebug > 1) {
1196          cout << "    Executing task " << task->GetName() << endl;
1197       }   
1198       task->ExecuteTask(option);
1199    }   
1200 //
1201 // Call FinishEvent() for optional output and MC services 
1202    if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
1203    if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
1204    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1205 }
1206
1207 //______________________________________________________________________________
1208 void AliAnalysisManager::FinishAnalysis()
1209 {
1210 // Finish analysis.
1211 }