]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
- AliHLTComponent: bugfix in buffer handling; overwrite check added to
[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
205    TIter next(fTasks);
206    AliAnalysisTask *task;
207    // Call CreateOutputObjects for all tasks
208    while ((task=(AliAnalysisTask*)next())) 
209       task->CreateOutputObjects();
210    if (fMode == kLocalAnalysis) Init(tree);   
211 }
212
213 //______________________________________________________________________________
214 Bool_t AliAnalysisManager::Notify()
215 {
216    // The Notify() function is called when a new file is opened. This
217    // can be either for a new TTree in a TChain or when when a new TTree
218    // is started when using PROOF. It is normaly not necessary to make changes
219    // to the generated code, but the routine can be extended by the
220    // user if needed. The return value is currently not used.
221    if (fTree) {
222       TFile *curfile = fTree->GetCurrentFile();
223       if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
224    }
225    return kTRUE;
226 }    
227
228 //______________________________________________________________________________
229 Bool_t AliAnalysisManager::Process(Long64_t entry)
230 {
231   // The Process() function is called for each entry in the tree (or possibly
232   // keyed object in the case of PROOF) to be processed. The entry argument
233   // specifies which entry in the currently loaded tree is to be processed.
234   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
235   // to read either all or the required parts of the data. When processing
236   // keyed objects with PROOF, the object is already loaded and is available
237   // via the fObject pointer.
238   //
239   // This function should contain the "body" of the analysis. It can contain
240   // simple or elaborate selection criteria, run algorithms on the data
241   // of the event and typically fill histograms.
242
243   // WARNING when a selector is used with a TChain, you must use
244   //  the pointer to the current TTree to call GetEntry(entry).
245   //  The entry is always the local entry number in the current tree.
246   //  Assuming that fChain is the pointer to the TChain being processed,
247   //  use fChain->GetTree()->GetEntry(entry).
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       AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
270       if (!top) {
271           cout<<"Error: No top input container !" <<endl;
272           return;
273       }
274       top->SetData(0);
275
276       TIter next(fOutputs);
277       AliAnalysisDataContainer *output;
278       while ((output=(AliAnalysisDataContainer*)next())) {
279          if (fDebug > 1) printf("   Packing container %s\n", output->GetName());
280          if (output->GetData()) target->Add(output);
281 //         output->SetDataOwned(kFALSE);
282       }
283    } 
284    fContainers->Clear();
285    if (fDebug > 1) {
286       printf("   ->output list contains %d containers\n", target->GetSize());
287    }
288 }
289
290 //______________________________________________________________________________
291 void AliAnalysisManager::ReplaceOutputContainers(TList *source)
292 {
293 // Replace all exising containers with the ones coming in source.
294    TIter next(fOutputs);
295    AliAnalysisDataContainer *cont, *output;
296    while ((cont=(AliAnalysisDataContainer*)next())) {
297       output = (AliAnalysisDataContainer*)source->FindObject(cont->GetName());
298       if (!output) {
299          printf("Error: container %s not found in analysis output !\n", cont->GetName());
300          continue;
301       }
302       if (fDebug > 1) printf("...Replacing output container %s\n", output->GetName());
303       if (cont->GetFileName()) printf("    -> %s\n", output->GetFileName());
304       Int_t ntasks = fTasks->GetEntries();
305       AliAnalysisTask *task;
306       AliAnalysisDataSlot *oslot;
307       for (Int_t i=0; i<ntasks; i++) {
308          task = (AliAnalysisTask*)fTasks->At(i);
309          Int_t nout = task->GetNoutputs();
310          for (Int_t iout=0; iout<nout; iout++) {
311             oslot = task->GetOutputSlot(iout);
312             if (oslot->GetContainer() == cont) oslot->ConnectContainer(output);
313          }
314       }
315 //      output->GetConsumers()->Delete();
316 //      if (output->GetProducer()) delete output->GetProducer();
317    }         
318 }
319
320 //______________________________________________________________________________
321 void AliAnalysisManager::UnpackOutput(TList *source)
322 {
323   // Called by AliAnalysisSelector::Terminate. Output containers should
324   // be in source in the same order as in fOutputs.
325
326    if (!source) {
327       Error("PackOutput", "No target. Aborting.");
328       return;
329    }
330    if (fDebug > 1) {
331       printf("AliAnalysisManager::UnpackOutput(): %d containers\n", source->GetSize());
332       printf("   Source list contains %d containers\n", source->GetSize());
333    }   
334
335    if (fMode == kProofAnalysis) {
336       ReplaceOutputContainers(source);
337       fOutputs->Clear();
338    }
339
340    TCollection *collection = source;
341    if (fMode == kLocalAnalysis) collection = fOutputs;
342    TIter next(collection);
343
344    AliAnalysisDataContainer *output;
345    while ((output=(AliAnalysisDataContainer*)next())) {
346       if (fMode == kProofAnalysis) {
347          output->SetDataOwned(kTRUE);
348          fOutputs->Add(output);
349       }
350       if (!output->GetData()) continue;
351       // Check if the output need to be written to a file.
352       const char *filename = output->GetFileName();
353       if (!filename || !strlen(filename)) continue;
354       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
355       if (file) file->cd();
356       else      file = new TFile(filename, "RECREATE");
357       if (file->IsZombie()) continue;
358       // Reparent data to this file
359       TMethodCall callEnv;
360       if (output->GetData()->IsA())
361          callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
362       if (callEnv.IsValid()) {
363          callEnv.SetParam((Long_t) file);
364          callEnv.Execute(output->GetData());
365       }
366       output->GetData()->Write();      
367    }
368 }
369
370 //______________________________________________________________________________
371 void AliAnalysisManager::Terminate()
372 {
373   // The Terminate() function is the last function to be called during
374   // a query. It always runs on the client, it can be used to present
375   // the results graphically.
376    if (fDebug > 1) {
377       cout << "AliAnalysisManager::Terminate()" << endl;
378    }   
379    AliAnalysisTask *task;
380    TIter next(fTasks);
381    // Call Terminate() for tasks
382    while ((task=(AliAnalysisTask*)next())) task->Terminate();
383 }
384
385 //______________________________________________________________________________
386 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
387 {
388 // Adds a user task to the global list of tasks.
389    task->SetActive(kFALSE);
390    fTasks->Add(task);
391 }  
392
393 //______________________________________________________________________________
394 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
395 {
396 // Retreive task by name.
397    if (!fTasks) return NULL;
398    return (AliAnalysisTask*)fTasks->FindObject(name);
399 }
400
401 //______________________________________________________________________________
402 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
403                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
404 {
405 // Create a data container of a certain type. Types can be:
406 //   kExchangeContainer  = 0, used to exchange date between tasks
407 //   kInputContainer   = 1, used to store input data
408 //   kOutputContainer  = 2, used for posting results
409    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
410    fContainers->Add(cont);
411    switch (type) {
412       case kInputContainer:
413          fInputs->Add(cont);
414          break;
415       case kOutputContainer:
416          if (fOutputs->FindObject(name)) printf("CreateContainer: warning: a container named %s existing !\n",name);
417          fOutputs->Add(cont);
418          if (filename && strlen(filename)) cont->SetFileName(filename);
419          break;
420       case kExchangeContainer:
421          break;   
422    }
423    return cont;
424 }
425          
426 //______________________________________________________________________________
427 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
428                                         AliAnalysisDataContainer *cont)
429 {
430 // Connect input of an existing task to a data container.
431    if (!fTasks->FindObject(task)) {
432       AddTask(task);
433       Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
434    } 
435    Bool_t connected = task->ConnectInput(islot, cont);
436    return connected;
437 }   
438
439 //______________________________________________________________________________
440 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
441                                         AliAnalysisDataContainer *cont)
442 {
443 // Connect output of an existing task to a data container.
444    if (!fTasks->FindObject(task)) {
445       AddTask(task);
446       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
447    } 
448    Bool_t connected = task->ConnectOutput(islot, cont);
449    return connected;
450 }   
451                                
452 //______________________________________________________________________________
453 void AliAnalysisManager::CleanContainers()
454 {
455 // Clean data from all containers that have already finished all client tasks.
456    TIter next(fContainers);
457    AliAnalysisDataContainer *cont;
458    while ((cont=(AliAnalysisDataContainer *)next())) {
459       if (cont->IsOwnedData() && 
460           cont->IsDataReady() && 
461           cont->ClientsExecuted()) cont->DeleteData();
462    }
463 }
464
465 //______________________________________________________________________________
466 Bool_t AliAnalysisManager::InitAnalysis()
467 {
468 // Initialization of analysis chain of tasks. Should be called after all tasks
469 // and data containers are properly connected
470    // Check for input/output containers
471    fInitOK = kFALSE;
472    // Check for top tasks (depending only on input data containers)
473    if (!fTasks->First()) {
474       Error("InitAnalysis", "Analysis has no tasks !");
475       return kFALSE;
476    }   
477    TIter next(fTasks);
478    AliAnalysisTask *task;
479    AliAnalysisDataContainer *cont;
480    Int_t ntop = 0;
481    Int_t nzombies = 0;
482    Bool_t iszombie = kFALSE;
483    Bool_t istop = kTRUE;
484    Int_t i;
485    while ((task=(AliAnalysisTask*)next())) {
486       istop = kTRUE;
487       iszombie = kFALSE;
488       Int_t ninputs = task->GetNinputs();
489       for (i=0; i<ninputs; i++) {
490          cont = task->GetInputSlot(i)->GetContainer();
491          if (!cont) {
492             if (!iszombie) {
493                task->SetZombie();
494                fZombies->Add(task);
495                nzombies++;
496                iszombie = kTRUE;
497             }   
498             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
499                   i, task->GetName()); 
500          }
501          if (iszombie) continue;
502          // Check if cont is an input container
503          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
504          // Connect to parent task
505       }
506       if (istop) {
507          ntop++;
508          fTopTasks->Add(task);
509       }
510    }
511    if (!ntop) {
512       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
513       return kFALSE;
514    }                        
515    // Check now if there are orphan tasks
516    for (i=0; i<ntop; i++) {
517       task = (AliAnalysisTask*)fTopTasks->At(i);
518       task->SetUsed();
519    }
520    Int_t norphans = 0;
521    next.Reset();
522    while ((task=(AliAnalysisTask*)next())) {
523       if (!task->IsUsed()) {
524          norphans++;
525          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
526       }   
527    }          
528    // Check the task hierarchy (no parent task should depend on data provided
529    // by a daughter task)
530    for (i=0; i<ntop; i++) {
531       task = (AliAnalysisTask*)fTopTasks->At(i);
532       if (task->CheckCircularDeps()) {
533          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
534          PrintStatus("dep");
535          return kFALSE;
536       }   
537    }
538    fInitOK = kTRUE;
539    return kTRUE;
540 }   
541
542 //______________________________________________________________________________
543 void AliAnalysisManager::PrintStatus(Option_t *option) const
544 {
545 // Print task hierarchy.
546    TIter next(fTopTasks);
547    AliAnalysisTask *task;
548    while ((task=(AliAnalysisTask*)next()))
549       task->PrintTask(option);
550 }
551
552 //______________________________________________________________________________
553 void AliAnalysisManager::ResetAnalysis()
554 {
555 // Reset all execution flags and clean containers.
556    CleanContainers();
557 }
558
559 //______________________________________________________________________________
560 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
561 {
562 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
563    if (!fInitOK) {
564       Error("StartAnalysis","Analysis manager was not initialized !");
565       return;
566    }
567    if (fDebug>1) {
568       cout << "StartAnalysis: " << GetName() << endl;   
569    }   
570    TString anaType = type;
571    anaType.ToLower();
572    fMode = kLocalAnalysis;
573    if (tree) {
574       if (anaType.Contains("proof"))     fMode = kProofAnalysis;
575       else if (anaType.Contains("grid")) fMode = kGridAnalysis;
576    }
577    if (fMode == kGridAnalysis) {
578       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
579       fMode = kLocalAnalysis; 
580    }     
581    char line[128];   
582    TChain *chain = dynamic_cast<TChain*>(tree);
583    switch (fMode) {
584       case kLocalAnalysis:
585          if (!tree) {
586             ExecAnalysis();
587             return;
588          } 
589          // Run tree-based analysis via AliAnalysisSelector  
590          gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
591          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
592          sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
593          gROOT->ProcessLine(line);
594          sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
595          gROOT->ProcessLine(line);
596          break;
597       case kProofAnalysis:
598          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
599             printf("StartAnalysis: no PROOF!!!\n");
600             return;
601          }   
602          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
603          gROOT->ProcessLine(line);
604          if (chain) {
605             chain->SetProof();
606             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
607             chain->Process(gSystem->ExpandPathName("$ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+"));
608          } else {
609             printf("StartAnalysis: no chain\n");
610             return;
611          }      
612          break;
613       case kGridAnalysis:
614          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
615    }   
616 }   
617
618 //______________________________________________________________________________
619 void AliAnalysisManager::ExecAnalysis(Option_t *option)
620 {
621 // Execute analysis.
622    if (!fInitOK) {
623      Error("ExecAnalysis", "Analysis manager was not initialized !");
624       return;
625    }   
626    AliAnalysisTask *task;
627    // Check if the top tree is active.
628    if (fTree) {
629       if (fDebug>1) {
630          printf("AliAnalysisManager::ExecAnalysis\n");
631       }   
632       TIter next(fTasks);
633    // De-activate all tasks
634       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
635       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
636       if (!cont) {
637               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
638          return;
639       }   
640       cont->SetData(fTree); // This will notify all consumers
641       TIter next1(cont->GetConsumers());
642       while ((task=(AliAnalysisTask*)next1())) {
643 //         task->SetActive(kTRUE);
644          if (fDebug >1) {
645             cout << "    Executing task " << task->GetName() << endl;
646          }   
647          task->ExecuteTask(option);
648       }
649       return;
650    }   
651    // The event loop is not controlled by TSelector   
652    TIter next2(fTopTasks);
653    while ((task=(AliAnalysisTask*)next2())) {
654       task->SetActive(kTRUE);
655       if (fDebug > 1) {
656          cout << "    Executing task " << task->GetName() << endl;
657       }   
658       task->ExecuteTask(option);
659    }   
660 }
661
662 //______________________________________________________________________________
663 void AliAnalysisManager::FinishAnalysis()
664 {
665 // Finish analysis.
666 }