]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
PROOF-aware version of the analysis framework (Andrei)
[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 <TMethodCall.h>
33 #include <TChain.h>
34 #include <TSystem.h>
35 #include <TROOT.h>
36
37 #include "AliAnalysisTask.h"
38 #include "AliAnalysisDataContainer.h"
39 #include "AliAnalysisDataSlot.h"
40 #include "AliAnalysisManager.h"
41
42 ClassImp(AliAnalysisManager)
43
44 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
45
46 //______________________________________________________________________________
47 AliAnalysisManager::AliAnalysisManager() 
48                    :TNamed(),
49                     fTree(NULL),
50                     fCurrentEntry(-1),
51                     fMode(kLocalAnalysis),
52                     fInitOK(kFALSE),
53                     fTasks(NULL),
54                     fTopTasks(NULL),
55                     fZombies(NULL),
56                     fContainers(NULL),
57                     fInputs(NULL),
58                     fOutputs(NULL)
59 {
60 // Dummy constructor.
61    fgAnalysisManager = this;
62 }
63
64 //______________________________________________________________________________
65 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
66                    :TNamed(name,title),
67                     fTree(NULL),
68                     fCurrentEntry(-1),
69                     fMode(kLocalAnalysis),
70                     fInitOK(kFALSE),
71                     fDebug(0),
72                     fTasks(NULL),
73                     fTopTasks(NULL),
74                     fZombies(NULL),
75                     fContainers(NULL),
76                     fInputs(NULL),
77                     fOutputs(NULL)
78 {
79 // Default constructor.
80    fgAnalysisManager = this;
81    fTasks      = new TObjArray();
82    fTopTasks   = new TObjArray();
83    fZombies    = new TObjArray();
84    fContainers = new TObjArray();
85    fInputs     = new TObjArray();
86    fOutputs    = new TObjArray();  
87 }
88
89 //______________________________________________________________________________
90 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
91                    :TNamed(other),
92                     fTree(NULL),
93                     fCurrentEntry(-1),
94                     fMode(other.fMode),
95                     fInitOK(other.fInitOK),
96                     fDebug(other.fDebug),
97                     fTasks(NULL),
98                     fTopTasks(NULL),
99                     fZombies(NULL),
100                     fContainers(NULL),
101                     fInputs(NULL),
102                     fOutputs(NULL)
103 {
104 // Copy constructor.
105    fTasks      = new TObjArray(*other.fTasks);
106    fTopTasks   = new TObjArray(*other.fTopTasks);
107    fZombies    = new TObjArray(*other.fZombies);
108    fContainers = new TObjArray(*other.fContainers);
109    fInputs     = new TObjArray(*other.fInputs);
110    fOutputs    = new TObjArray(*other.fOutputs);
111    fgAnalysisManager = this;
112 }
113    
114 //______________________________________________________________________________
115 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
116 {
117 // Assignment
118    if (&other != this) {
119       TNamed::operator=(other);
120       fTree       = NULL;
121       fCurrentEntry = -1;
122       fMode       = other.fMode;
123       fInitOK     = other.fInitOK;
124       fDebug      = other.fDebug;
125       fTasks      = new TObjArray(*other.fTasks);
126       fTopTasks   = new TObjArray(*other.fTopTasks);
127       fZombies    = new TObjArray(*other.fZombies);
128       fContainers = new TObjArray(*other.fContainers);
129       fInputs     = new TObjArray(*other.fInputs);
130       fOutputs    = new TObjArray(*other.fOutputs);
131       fgAnalysisManager = this;
132    }
133    return *this;
134 }
135
136 //______________________________________________________________________________
137 AliAnalysisManager::~AliAnalysisManager()
138 {
139 // Destructor.
140    if (fTasks) {fTasks->Delete(); delete fTasks;}
141    if (fTopTasks) delete fTopTasks;
142    if (fZombies) delete fZombies;
143    if (fContainers) {fContainers->Delete(); delete fContainers;}
144    if (fInputs) delete fInputs;
145    if (fOutputs) delete fOutputs;
146    if (fgAnalysisManager==this) fgAnalysisManager = NULL;
147 }
148
149 //______________________________________________________________________________
150 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
151 {
152 // Read one entry of the tree or a whole branch.
153    if (fDebug > 1) {
154       cout << "== AliAnalysisManager::GetEntry()" << endl;
155    }   
156    fCurrentEntry = entry;
157    return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
158 }
159    
160 //______________________________________________________________________________
161 void AliAnalysisManager::Init(TTree *tree)
162 {
163   // The Init() function is called when the selector needs to initialize
164   // a new tree or chain. Typically here the branch addresses of the tree
165   // will be set. It is normaly not necessary to make changes to the
166   // generated code, but the routine can be extended by the user if needed.
167   // Init() will be called many times when running with PROOF.
168    if (!tree) return;
169    if (fDebug > 1) {
170       printf("AliAnalysisManager::Init(%s)\n", tree->GetName());
171    }
172    if (!fInitOK) InitAnalysis();
173    if (!fInitOK) return;
174    fTree = tree;
175    AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
176    if (!top) {
177       cout<<"Error: No top input container !" <<endl;
178       return;
179    }   
180    top->SetData(tree);
181 }
182
183 //______________________________________________________________________________
184 void AliAnalysisManager::Begin(TTree *tree)
185 {
186   // The Begin() function is called at the start of the query.
187   // When running with PROOF Begin() is only called on the client.
188   // The tree argument is deprecated (on PROOF 0 is passed).
189    if (fDebug > 1) {
190       cout << "AliAnalysisManager::Begin()" << endl;
191    }   
192    Init(tree);
193 }
194
195 //______________________________________________________________________________
196 void AliAnalysisManager::SlaveBegin(TTree *tree)
197 {
198   // The SlaveBegin() function is called after the Begin() function.
199   // When running with PROOF SlaveBegin() is called on each slave server.
200   // The tree argument is deprecated (on PROOF 0 is passed).
201    if (fDebug > 1) {
202       cout << "AliAnalysisManager::SlaveBegin()" << endl;
203    }   
204    TIter next(fTasks);
205    AliAnalysisTask *task;
206    // Call CreateOutputObjects for all tasks
207    while ((task=(AliAnalysisTask*)next())) 
208       task->CreateOutputObjects();
209    if (fMode == kLocalAnalysis) Init(tree);   
210 }
211
212 //______________________________________________________________________________
213 Bool_t AliAnalysisManager::Notify()
214 {
215    // The Notify() function is called when a new file is opened. This
216    // can be either for a new TTree in a TChain or when when a new TTree
217    // is started when using PROOF. It is normaly not necessary to make changes
218    // to the generated code, but the routine can be extended by the
219    // user if needed. The return value is currently not used.
220    if (fTree) {
221       TFile *curfile = fTree->GetCurrentFile();
222       if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
223    }
224    return kTRUE;
225 }    
226
227 //______________________________________________________________________________
228 Bool_t AliAnalysisManager::Process(Long64_t entry)
229 {
230   // The Process() function is called for each entry in the tree (or possibly
231   // keyed object in the case of PROOF) to be processed. The entry argument
232   // specifies which entry in the currently loaded tree is to be processed.
233   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
234   // to read either all or the required parts of the data. When processing
235   // keyed objects with PROOF, the object is already loaded and is available
236   // via the fObject pointer.
237   //
238   // This function should contain the "body" of the analysis. It can contain
239   // simple or elaborate selection criteria, run algorithms on the data
240   // of the event and typically fill histograms.
241
242   // WARNING when a selector is used with a TChain, you must use
243   //  the pointer to the current TTree to call GetEntry(entry).
244   //  The entry is always the local entry number in the current tree.
245   //  Assuming that fChain is the pointer to the TChain being processed,
246   //  use fChain->GetTree()->GetEntry(entry).
247   
248    if (fDebug > 1) {
249       cout << "AliAnalysisManager::Process()" << endl;
250    }   
251    GetEntry(entry);
252    ExecAnalysis();
253    return kTRUE;
254 }
255
256 //______________________________________________________________________________
257 void AliAnalysisManager::PackOutput(TList *target)
258 {
259   // Pack all output data containers in the output list.
260    if (fDebug > 1) {
261       cout << "AliAnalysisManager::PackOutput()" << endl;
262    }   
263    if (!target) {
264       Error("PackOutput", "No target. Aborting.");
265       return;
266    }   
267
268    if (fMode == kProofAnalysis) {
269       TIter next(fOutputs);
270       AliAnalysisDataContainer *output;
271       while ((output=(AliAnalysisDataContainer*)next())) {
272          if (fDebug > 1) printf("   Packing container %s\n", output->GetName());
273          if (output->GetData()) target->Add(output);
274 //         output->SetDataOwned(kFALSE);
275       }
276    } 
277    fContainers->Clear();
278    if (fDebug > 1) {
279       printf("   ->output list contains %d containers\n", target->GetSize());  
280    }   
281 }
282
283 //______________________________________________________________________________
284 void AliAnalysisManager::ReplaceOutputContainers(TList *source)
285 {
286 // Replace all exising containers with the ones coming in source.
287    TIter next(fOutputs);
288    AliAnalysisDataContainer *cont, *output;
289    while ((cont=(AliAnalysisDataContainer*)next())) {
290       output = (AliAnalysisDataContainer*)source->FindObject(cont->GetName());
291       if (!output) {
292          printf("Error: container %s not found in analysis output !\n", cont->GetName());
293          continue;
294       }
295       if (fDebug > 1) printf("...Replacing output container %s\n", output->GetName());
296       if (cont->GetFileName()) printf("    -> %s\n", output->GetFileName());
297       Int_t ntasks = fTasks->GetEntries();
298       AliAnalysisTask *task;
299       AliAnalysisDataSlot *oslot;
300       for (Int_t i=0; i<ntasks; i++) {  
301          task = (AliAnalysisTask*)fTasks->At(i);
302          Int_t nout = task->GetNoutputs();
303          for (Int_t iout=0; iout<nout; iout++) {
304             oslot = task->GetOutputSlot(iout);
305             if (oslot->GetContainer() == cont) oslot->ConnectContainer(output);
306          }
307       }
308 //      output->GetConsumers()->Delete();
309 //      if (output->GetProducer()) delete output->GetProducer();
310    }         
311 }
312
313 //______________________________________________________________________________
314 void AliAnalysisManager::UnpackOutput(TList *source)
315 {
316   // Called by AliAnalysisSelector::Terminate. Output containers should
317   // be in source in the same order as in fOutputs.
318
319    if (!source) {
320       Error("PackOutput", "No target. Aborting.");
321       return;
322    }
323    if (fDebug > 1) {
324       printf("AliAnalysisManager::UnpackOutput(): %d containers\n", source->GetSize());
325       printf("   Source list contains %d containers\n", source->GetSize());
326    }   
327    TCollection *collection = source;
328    if (fMode == kLocalAnalysis) collection = fOutputs;
329    TIter next(collection);
330
331    if (fMode == kProofAnalysis) {
332       ReplaceOutputContainers(source);
333       fOutputs->Clear();
334    }   
335    AliAnalysisDataContainer *output;
336    while ((output=(AliAnalysisDataContainer*)next())) {
337       if (fMode == kProofAnalysis) {
338          output->SetDataOwned(kTRUE);
339          fOutputs->Add(output);
340       }   
341       if (!output->GetData()) continue;
342       // Check if the output need to be written to a file.
343       const char *filename = output->GetFileName();
344       if (!filename || !strlen(filename)) continue;
345       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
346       if (file) file->cd();
347       else      file = new TFile(filename, "RECREATE");
348       if (file->IsZombie()) continue;
349       // Reparent data to this file
350       TMethodCall callEnv;
351       if (output->GetData()->IsA())
352          callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
353       if (callEnv.IsValid()) {
354          callEnv.SetParam((Long_t) file);
355          callEnv.Execute(output->GetData());
356       }
357       output->GetData()->Write();      
358    }
359 }
360
361 //______________________________________________________________________________
362 void AliAnalysisManager::Terminate()
363 {
364   // The Terminate() function is the last function to be called during
365   // a query. It always runs on the client, it can be used to present
366   // the results graphically.
367    if (fDebug > 1) {
368       cout << "AliAnalysisManager::Terminate()" << endl;
369    }   
370    AliAnalysisTask *task;
371    TIter next(fTasks);
372    // Call Terminate() for tasks
373    while ((task=(AliAnalysisTask*)next())) task->Terminate();
374 }
375
376 //______________________________________________________________________________
377 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
378 {
379 // Adds a user task to the global list of tasks.
380    task->SetActive(kFALSE);
381    fTasks->Add(task);
382 }  
383
384 //______________________________________________________________________________
385 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
386 {
387 // Retreive task by name.
388    if (!fTasks) return NULL;
389    return (AliAnalysisTask*)fTasks->FindObject(name);
390 }
391
392 //______________________________________________________________________________
393 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
394                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
395 {
396 // Create a data container of a certain type. Types can be:
397 //   kExchangeContainer  = 0, used to exchange date between tasks
398 //   kInputContainer   = 1, used to store input data
399 //   kOutputContainer  = 2, used for posting results
400    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
401    fContainers->Add(cont);
402    switch (type) {
403       case kInputContainer:
404          fInputs->Add(cont);
405          break;
406       case kOutputContainer:
407          if (fOutputs->FindObject(name)) printf("CreateContainer: warning: a container named %s existing !\n",name);
408          fOutputs->Add(cont);
409          if (filename && strlen(filename)) cont->SetFileName(filename);
410          break;
411       case kExchangeContainer:
412          break;   
413    }
414    return cont;
415 }
416          
417 //______________________________________________________________________________
418 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
419                                         AliAnalysisDataContainer *cont)
420 {
421 // Connect input of an existing task to a data container.
422    if (!fTasks->FindObject(task)) {
423       AddTask(task);
424       Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
425    } 
426    Bool_t connected = task->ConnectInput(islot, cont);
427    return connected;
428 }   
429
430 //______________________________________________________________________________
431 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
432                                         AliAnalysisDataContainer *cont)
433 {
434 // Connect output of an existing task to a data container.
435    if (!fTasks->FindObject(task)) {
436       AddTask(task);
437       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
438    } 
439    Bool_t connected = task->ConnectOutput(islot, cont);
440    return connected;
441 }   
442                                
443 //______________________________________________________________________________
444 void AliAnalysisManager::CleanContainers()
445 {
446 // Clean data from all containers that have already finished all client tasks.
447    TIter next(fContainers);
448    AliAnalysisDataContainer *cont;
449    while ((cont=(AliAnalysisDataContainer *)next())) {
450       if (cont->IsOwnedData() && 
451           cont->IsDataReady() && 
452           cont->ClientsExecuted()) cont->DeleteData();
453    }
454 }
455
456 //______________________________________________________________________________
457 Bool_t AliAnalysisManager::InitAnalysis()
458 {
459 // Initialization of analysis chain of tasks. Should be called after all tasks
460 // and data containers are properly connected
461    // Check for input/output containers
462    fInitOK = kFALSE;
463    // Check for top tasks (depending only on input data containers)
464    if (!fTasks->First()) {
465       Error("InitAnalysis", "Analysis has no tasks !");
466       return kFALSE;
467    }   
468    TIter next(fTasks);
469    AliAnalysisTask *task;
470    AliAnalysisDataContainer *cont;
471    Int_t ntop = 0;
472    Int_t nzombies = 0;
473    Bool_t iszombie = kFALSE;
474    Bool_t istop = kTRUE;
475    Int_t i;
476    while ((task=(AliAnalysisTask*)next())) {
477       istop = kTRUE;
478       iszombie = kFALSE;
479       Int_t ninputs = task->GetNinputs();
480       for (i=0; i<ninputs; i++) {
481          cont = task->GetInputSlot(i)->GetContainer();
482          if (!cont) {
483             if (!iszombie) {
484                task->SetZombie();
485                fZombies->Add(task);
486                nzombies++;
487                iszombie = kTRUE;
488             }   
489             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
490                   i, task->GetName()); 
491          }
492          if (iszombie) continue;
493          // Check if cont is an input container
494          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
495          // Connect to parent task
496       }
497       if (istop) {
498          ntop++;
499          fTopTasks->Add(task);
500       }
501    }
502    if (!ntop) {
503       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
504       return kFALSE;
505    }                        
506    // Check now if there are orphan tasks
507    for (i=0; i<ntop; i++) {
508       task = (AliAnalysisTask*)fTopTasks->At(i);
509       task->SetUsed();
510    }
511    Int_t norphans = 0;
512    next.Reset();
513    while ((task=(AliAnalysisTask*)next())) {
514       if (!task->IsUsed()) {
515          norphans++;
516          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
517       }   
518    }          
519    // Check the task hierarchy (no parent task should depend on data provided
520    // by a daughter task)
521    for (i=0; i<ntop; i++) {
522       task = (AliAnalysisTask*)fTopTasks->At(i);
523       if (task->CheckCircularDeps()) {
524          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
525          PrintStatus("dep");
526          return kFALSE;
527       }   
528    }
529    fInitOK = kTRUE;
530    return kTRUE;
531 }   
532
533 //______________________________________________________________________________
534 void AliAnalysisManager::PrintStatus(Option_t *option) const
535 {
536 // Print task hierarchy.
537    TIter next(fTopTasks);
538    AliAnalysisTask *task;
539    while ((task=(AliAnalysisTask*)next()))
540       task->PrintTask(option);
541 }
542
543 //______________________________________________________________________________
544 void AliAnalysisManager::ResetAnalysis()
545 {
546 // Reset all execution flags and clean containers.
547    CleanContainers();
548 }
549
550 //______________________________________________________________________________
551 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
552 {
553 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
554    if (!fInitOK) {
555       Error("StartAnalysis","Analysis manager was not initialized !");
556       return;
557    }
558    if (fDebug>1) {
559       cout << "StartAnalysis: " << GetName() << endl;   
560    }   
561    TString anaType = type;
562    anaType.ToLower();
563    fMode = kLocalAnalysis;
564    if (tree) {
565       if (anaType.Contains("proof"))     fMode = kProofAnalysis;
566       else if (anaType.Contains("grid")) fMode = kGridAnalysis;
567    }
568    if (fMode == kGridAnalysis) {
569       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
570       fMode = kLocalAnalysis; 
571    }     
572    char line[128];   
573    TChain *chain = dynamic_cast<TChain*>(tree);
574    switch (fMode) {
575       case kLocalAnalysis:
576          if (!tree) {
577             ExecAnalysis();
578             return;
579          } 
580          // Run tree-based analysis via AliAnalysisSelector  
581          gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
582          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
583          sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
584          gROOT->ProcessLine(line);
585          sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
586          gROOT->ProcessLine(line);
587          break;
588       case kProofAnalysis:
589          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
590             printf("StartAnalysis: no PROOF!!!\n");
591             return;
592          }   
593          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
594          gROOT->ProcessLine(line);
595          if (chain) {
596             chain->SetProof();
597             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
598             chain->Process(gSystem->ExpandPathName("$ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+"));
599          } else {
600             printf("StartAnalysis: no chain\n");
601             return;
602          }      
603          break;
604       case kGridAnalysis:
605          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
606    }   
607 }   
608
609 //______________________________________________________________________________
610 void AliAnalysisManager::ExecAnalysis(Option_t *option)
611 {
612 // Execute analysis.
613    if (!fInitOK) {
614      Error("ExecAnalysis", "Analysis manager was not initialized !");
615       return;
616    }   
617    AliAnalysisTask *task;
618    // Check if the top tree is active.
619    if (fTree) {
620       if (fDebug>1) {
621          printf("AliAnalysisManager::ExecAnalysis\n");
622       }   
623       TIter next(fTasks);
624    // De-activate all tasks
625       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
626       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
627       if (!cont) {
628               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
629          return;
630       }   
631       cont->SetData(fTree); // This will notify all consumers
632       TIter next1(cont->GetConsumers());
633       while ((task=(AliAnalysisTask*)next1())) {
634 //         task->SetActive(kTRUE);
635          if (fDebug >1) {
636             cout << "    Executing task " << task->GetName() << endl;
637          }   
638          task->ExecuteTask(option);
639       }
640       return;
641    }   
642    // The event loop is not controlled by TSelector   
643    TIter next2(fTopTasks);
644    while ((task=(AliAnalysisTask*)next2())) {
645       task->SetActive(kTRUE);
646       if (fDebug > 1) {
647          cout << "    Executing task " << task->GetName() << endl;
648       }   
649       task->ExecuteTask(option);
650    }   
651 }
652
653 //______________________________________________________________________________
654 void AliAnalysisManager::FinishAnalysis()
655 {
656 // Finish analysis.
657 }