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