]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
975bc28a0026eb32abec6587cb07012f48274c8f
[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 #include <TCanvas.h>
37
38 #include "AliAnalysisTask.h"
39 #include "AliAnalysisDataContainer.h"
40 #include "AliAnalysisDataSlot.h"
41 #include "AliVEventHandler.h"
42 #include "AliSysInfo.h"
43 #include "AliAnalysisManager.h"
44
45 ClassImp(AliAnalysisManager)
46
47 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
48
49 //______________________________________________________________________________
50 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
51                    :TNamed(name,title),
52                     fTree(NULL),
53                     fInputEventHandler(NULL),
54                     fOutputEventHandler(NULL),
55                     fMCtruthEventHandler(NULL),
56                     fCurrentEntry(-1),
57                     fNSysInfo(0),
58                     fMode(kLocalAnalysis),
59                     fInitOK(kFALSE),
60                     fDebug(0),
61                     fTasks(NULL),
62                     fTopTasks(NULL),
63                     fZombies(NULL),
64                     fContainers(NULL),
65                     fInputs(NULL),
66                     fOutputs(NULL)
67 {
68 // Default constructor.
69    fgAnalysisManager = this;
70    fTasks      = new TObjArray();
71    fTopTasks   = new TObjArray();
72    fZombies    = new TObjArray();
73    fContainers = new TObjArray();
74    fInputs     = new TObjArray();
75    fOutputs    = new TObjArray();
76    SetEventLoop(kTRUE);
77 }
78
79 //______________________________________________________________________________
80 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
81                    :TNamed(other),
82                     fTree(NULL),
83                     fInputEventHandler(NULL),
84                     fOutputEventHandler(NULL),
85                     fMCtruthEventHandler(NULL),
86                     fCurrentEntry(-1),
87                     fNSysInfo(0),
88                     fMode(other.fMode),
89                     fInitOK(other.fInitOK),
90                     fDebug(other.fDebug),
91                     fTasks(NULL),
92                     fTopTasks(NULL),
93                     fZombies(NULL),
94                     fContainers(NULL),
95                     fInputs(NULL),
96                     fOutputs(NULL)
97 {
98 // Copy constructor.
99    fTasks      = new TObjArray(*other.fTasks);
100    fTopTasks   = new TObjArray(*other.fTopTasks);
101    fZombies    = new TObjArray(*other.fZombies);
102    fContainers = new TObjArray(*other.fContainers);
103    fInputs     = new TObjArray(*other.fInputs);
104    fOutputs    = new TObjArray(*other.fOutputs);
105    fgAnalysisManager = this;
106 }
107    
108 //______________________________________________________________________________
109 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
110 {
111 // Assignment
112    if (&other != this) {
113       TNamed::operator=(other);
114       fInputEventHandler   = other.fInputEventHandler;
115       fOutputEventHandler  = other.fOutputEventHandler;
116       fMCtruthEventHandler = other.fMCtruthEventHandler;
117       fTree       = NULL;
118       fCurrentEntry = -1;
119       fNSysInfo   = other.fNSysInfo;
120       fMode       = other.fMode;
121       fInitOK     = other.fInitOK;
122       fDebug      = other.fDebug;
123       fTasks      = new TObjArray(*other.fTasks);
124       fTopTasks   = new TObjArray(*other.fTopTasks);
125       fZombies    = new TObjArray(*other.fZombies);
126       fContainers = new TObjArray(*other.fContainers);
127       fInputs     = new TObjArray(*other.fInputs);
128       fOutputs    = new TObjArray(*other.fOutputs);
129       fgAnalysisManager = this;
130    }
131    return *this;
132 }
133
134 //______________________________________________________________________________
135 AliAnalysisManager::~AliAnalysisManager()
136 {
137 // Destructor.
138    if (fTasks) {fTasks->Delete(); delete fTasks;}
139    if (fTopTasks) delete fTopTasks;
140    if (fZombies) delete fZombies;
141    if (fContainers) {fContainers->Delete(); delete fContainers;}
142    if (fInputs) delete fInputs;
143    if (fOutputs) delete fOutputs;
144    if (fgAnalysisManager==this) fgAnalysisManager = NULL;
145 }
146
147 //______________________________________________________________________________
148 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
149 {
150 // Read one entry of the tree or a whole branch.
151    if (fDebug > 1) {
152       cout << "== AliAnalysisManager::GetEntry()" << endl;
153    }   
154    fCurrentEntry = entry;
155    return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
156 }
157    
158 //______________________________________________________________________________
159 void AliAnalysisManager::Init(TTree *tree)
160 {
161   // The Init() function is called when the selector needs to initialize
162   // a new tree or chain. Typically here the branch addresses of the tree
163   // will be set. It is normaly not necessary to make changes to the
164   // generated code, but the routine can be extended by the user if needed.
165   // Init() will be called many times when running with PROOF.
166    if (!tree) return;
167    if (fDebug > 1) {
168       printf("->AliAnalysisManager::InitTree(%s)\n", tree->GetName());
169    }
170
171    // Call InitTree of EventHandler
172    if (fOutputEventHandler) {
173       if (fMode == kProofAnalysis) {
174          fOutputEventHandler->Init(0x0, "proof");
175       } else {
176          fOutputEventHandler->Init(0x0, "local");
177       }
178    }
179
180    if (fInputEventHandler) {
181       if (fMode == kProofAnalysis) {
182          fInputEventHandler->Init(tree, "proof");
183       } else {
184          fInputEventHandler->Init(tree, "local");
185       }
186    } else {
187       // If no input event handler we need to get the tree once
188       // for the chain
189       if(!tree->GetTree()) tree->LoadTree(0);
190    }
191    
192
193    if (fMCtruthEventHandler) {
194       if (fMode == kProofAnalysis) {
195          fMCtruthEventHandler->Init(0x0, "proof");
196       } else {
197          fMCtruthEventHandler->Init(0x0, "local");
198       }
199    }
200
201    if (!fInitOK) InitAnalysis();
202    if (!fInitOK) return;
203    fTree = tree;
204    AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
205    if (!top) {
206       cout<<"Error: No top input container !" <<endl;
207       return;
208    }
209    top->SetData(tree);
210    if (fDebug > 1) {
211       printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
212    }
213 }
214
215 //______________________________________________________________________________
216 void AliAnalysisManager::SlaveBegin(TTree *tree)
217 {
218   // The SlaveBegin() function is called after the Begin() function.
219   // When running with PROOF SlaveBegin() is called on each slave server.
220   // The tree argument is deprecated (on PROOF 0 is passed).
221    if (fDebug > 1) {
222       cout << "->AliAnalysisManager::SlaveBegin()" << endl;
223    }
224
225    // Call Init of EventHandler
226    if (fOutputEventHandler) {
227       if (fMode == kProofAnalysis) {
228          fOutputEventHandler->Init("proof");
229       } else {
230          fOutputEventHandler->Init("local");
231       }
232    }
233
234    if (fInputEventHandler) {
235       fInputEventHandler->SetInputTree(tree);
236       if (fMode == kProofAnalysis) {
237          fInputEventHandler->Init("proof");
238       } else {
239          fInputEventHandler->Init("local");
240       }
241    }
242
243    if (fMCtruthEventHandler) {
244       if (fMode == kProofAnalysis) {
245          fMCtruthEventHandler->Init("proof");
246       } else {
247          fMCtruthEventHandler->Init("local");
248       }
249    }
250
251    TIter next(fTasks);
252    AliAnalysisTask *task;
253    // Call CreateOutputObjects for all tasks
254    while ((task=(AliAnalysisTask*)next())) {
255       TDirectory *curdir = gDirectory;
256       task->CreateOutputObjects();
257       if (curdir) curdir->cd();
258    }
259
260    if (fDebug > 1) {
261       cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
262    }
263 }
264
265 //______________________________________________________________________________
266 Bool_t AliAnalysisManager::Notify()
267 {
268    // The Notify() function is called when a new file is opened. This
269    // can be either for a new TTree in a TChain or when when a new TTree
270    // is started when using PROOF. It is normaly not necessary to make changes
271    // to the generated code, but the routine can be extended by the
272    // user if needed. The return value is currently not used.
273     if (fTree) {
274         TFile *curfile = fTree->GetCurrentFile();
275         if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
276         TIter next(fTasks);
277         AliAnalysisTask *task;
278         // Call Notify for all tasks
279         while ((task=(AliAnalysisTask*)next())) 
280             task->Notify();
281         
282         // Call Notify of the event handlers
283         if (fInputEventHandler) {
284             fInputEventHandler->Notify(curfile->GetName());
285         }
286
287         if (fOutputEventHandler) {
288             fOutputEventHandler->Notify(curfile->GetName());
289         }
290
291         if (fMCtruthEventHandler) {
292             fMCtruthEventHandler->Notify(curfile->GetName());
293         }
294
295     }
296     return kTRUE;
297 }    
298
299 //______________________________________________________________________________
300 Bool_t AliAnalysisManager::Process(Long64_t entry)
301 {
302   // The Process() function is called for each entry in the tree (or possibly
303   // keyed object in the case of PROOF) to be processed. The entry argument
304   // specifies which entry in the currently loaded tree is to be processed.
305   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
306   // to read either all or the required parts of the data. When processing
307   // keyed objects with PROOF, the object is already loaded and is available
308   // via the fObject pointer.
309   //
310   // This function should contain the "body" of the analysis. It can contain
311   // simple or elaborate selection criteria, run algorithms on the data
312   // of the event and typically fill histograms.
313
314   // WARNING when a selector is used with a TChain, you must use
315   //  the pointer to the current TTree to call GetEntry(entry).
316   //  The entry is always the local entry number in the current tree.
317   //  Assuming that fChain is the pointer to the TChain being processed,
318   //  use fChain->GetTree()->GetEntry(entry).
319    if (fDebug > 1) {
320       cout << "->AliAnalysisManager::Process(" << entry << ")" << endl;
321    }
322    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
323    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
324    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
325    
326    GetEntry(entry);
327    ExecAnalysis();
328    if (fDebug > 1) {
329       cout << "<-AliAnalysisManager::Process()" << endl;
330    }
331    return kTRUE;
332 }
333
334 //______________________________________________________________________________
335 void AliAnalysisManager::PackOutput(TList *target)
336 {
337   // Pack all output data containers in the output list. Called at SlaveTerminate
338   // stage in PROOF case for each slave.
339    if (fDebug > 1) {
340       cout << "->AliAnalysisManager::PackOutput()" << endl;
341    }   
342    if (!target) {
343       Error("PackOutput", "No target. Aborting.");
344       return;
345    }
346    if (fInputEventHandler)   fInputEventHandler  ->Terminate();
347    if (fOutputEventHandler)  fOutputEventHandler ->Terminate();
348    if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
349    
350    if (fMode == kProofAnalysis) {
351       TIter next(fOutputs);
352       AliAnalysisDataContainer *output;
353       while ((output=(AliAnalysisDataContainer*)next())) {
354          if (output->GetData()) {
355             if (output->GetProducer()->IsPostEventLoop()) continue;
356             AliAnalysisDataWrapper *wrap = output->ExportData();
357             // Output wrappers must delete data after merging (AG 13/11/07)
358             wrap->SetDeleteData(kTRUE);
359             if (fDebug > 1) printf("   Packing container %s...\n", output->GetName());
360             target->Add(wrap);
361          }   
362       }
363    } 
364    if (fDebug > 1) {
365       printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
366    }
367 }
368
369 //______________________________________________________________________________
370 void AliAnalysisManager::ImportWrappers(TList *source)
371 {
372 // Import data in output containers from wrappers coming in source.
373    if (fDebug > 1) {
374       cout << "->AliAnalysisManager::ImportWrappers()" << endl;
375    }   
376    TIter next(fOutputs);
377    AliAnalysisDataContainer *cont;
378    AliAnalysisDataWrapper   *wrap;
379    Int_t icont = 0;
380    while ((cont=(AliAnalysisDataContainer*)next())) {
381       if (cont->GetProducer()->IsPostEventLoop()) continue;
382       wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
383       if (!wrap && fDebug>1) {
384          printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
385          continue;
386       }
387       icont++;
388       if (fDebug > 1) printf("   Importing data for container %s\n", wrap->GetName());
389       if (cont->GetFileName()) printf("    -> %s\n", cont->GetFileName());
390       cont->ImportData(wrap);
391    }         
392    if (fDebug > 1) {
393       cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
394    }   
395 }
396
397 //______________________________________________________________________________
398 void AliAnalysisManager::UnpackOutput(TList *source)
399 {
400   // Called by AliAnalysisSelector::Terminate. Output containers should
401   // be in source in the same order as in fOutputs.
402    if (fDebug > 1) {
403       cout << "->AliAnalysisManager::UnpackOutput()" << endl;
404    }   
405    if (!source) {
406       Error("UnpackOutput", "No target. Aborting.");
407       return;
408    }
409    if (fDebug > 1) {
410       printf("   Source list contains %d containers\n", source->GetSize());
411    }   
412
413    if (fMode == kProofAnalysis) ImportWrappers(source);
414
415    TIter next(fOutputs);
416    AliAnalysisDataContainer *output;
417    while ((output=(AliAnalysisDataContainer*)next())) {
418       if (!output->GetData()) continue;
419       // Check if there are client tasks that run post event loop
420       if (output->HasConsumers()) {
421          // Disable event loop semaphore
422          output->SetPostEventLoop(kTRUE);
423          TObjArray *list = output->GetConsumers();
424          Int_t ncons = list->GetEntriesFast();
425          for (Int_t i=0; i<ncons; i++) {
426             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
427             task->CheckNotify(kTRUE);
428             // If task is active, execute it
429             if (task->IsPostEventLoop() && task->IsActive()) {
430                if (fDebug > 1) {
431                   cout << "== Executing post event loop task " << task->GetName() << endl;
432                }                  
433                task->ExecuteTask();
434             }   
435          }
436       }   
437       // Check if the output need to be written to a file.
438       const char *filename = output->GetFileName();
439       if (!(strcmp(filename, "default"))) {
440           if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
441       }
442       
443       if (!filename || !strlen(filename)) continue;
444       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
445       TDirectory *opwd = gDirectory;
446       if (file) file->cd();
447       else      file = new TFile(filename, "RECREATE");
448       if (file->IsZombie()) continue;
449       // Reparent data to this file
450       TMethodCall callEnv;
451       if (output->GetData()->IsA())
452          callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
453       if (callEnv.IsValid()) {
454          callEnv.SetParam((Long_t) file);
455          callEnv.Execute(output->GetData());
456       }
457       output->GetData()->Write();
458       if (opwd) opwd->cd();
459    }
460    if (fDebug > 1) {
461       cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
462    }   
463 }
464
465 //______________________________________________________________________________
466 void AliAnalysisManager::Terminate()
467 {
468   // The Terminate() function is the last function to be called during
469   // a query. It always runs on the client, it can be used to present
470   // the results graphically.
471    if (fDebug > 1) {
472       cout << "->AliAnalysisManager::Terminate()" << endl;
473    }   
474    AliAnalysisTask *task;
475    TIter next(fTasks);
476    // Call Terminate() for tasks
477    while ((task=(AliAnalysisTask*)next())) task->Terminate();
478    if (fDebug > 1) {
479       cout << "<-AliAnalysisManager::Terminate()" << endl;
480    }   
481    //
482    if (fInputEventHandler)   fInputEventHandler  ->TerminateIO();
483    if (fOutputEventHandler)  fOutputEventHandler ->TerminateIO();
484    if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
485    TIter next1(fOutputs);
486    AliAnalysisDataContainer *output;
487    while ((output=(AliAnalysisDataContainer*)next1())) {
488       // Close all files at output
489       const char *filename = output->GetFileName();
490       if (!(strcmp(filename, "default"))) {
491          if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
492       }
493       
494       if (!filename || !strlen(filename)) continue;
495       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
496       if (!file || file->IsZombie()) continue;
497       file->Close();
498    }   
499
500    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
501    if (getsysInfo) {
502       TDirectory *cdir = gDirectory;
503       TFile f("syswatch.root", "RECREATE");
504       if (!f.IsZombie()) {
505          TTree *tree = AliSysInfo::MakeTree("syswatch.log");
506          tree->SetMarkerStyle(kCircle);
507          tree->SetMarkerColor(kBlue);
508          tree->SetMarkerSize(0.5);
509          if (!gROOT->IsBatch()) {
510             tree->SetAlias("event", "id0");
511             tree->SetAlias("memUSED", "mi.fMemUsed");
512             tree->SetAlias("userCPU", "pI.fCpuUser");
513             TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
514             c->Divide(2,1,0.01,0.01);
515             c->cd(1);
516             tree->Draw("memUSED:event","","", 1234567890, 0);
517             c->cd(2);
518             tree->Draw("userCPU:event","","", 1234567890, 0);
519          }   
520          tree->Write();
521          f.Close();
522          delete tree;
523       }
524       if (cdir) cdir->cd();
525    }      
526 }
527
528 //______________________________________________________________________________
529 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
530 {
531 // Adds a user task to the global list of tasks.
532    task->SetActive(kFALSE);
533    fTasks->Add(task);
534 }  
535
536 //______________________________________________________________________________
537 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
538 {
539 // Retreive task by name.
540    if (!fTasks) return NULL;
541    return (AliAnalysisTask*)fTasks->FindObject(name);
542 }
543
544 //______________________________________________________________________________
545 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
546                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
547 {
548 // Create a data container of a certain type. Types can be:
549 //   kExchangeContainer  = 0, used to exchange date between tasks
550 //   kInputContainer   = 1, used to store input data
551 //   kOutputContainer  = 2, used for posting results
552    if (fContainers->FindObject(name)) {
553       Error("CreateContainer","A container named %s already defined !\n",name);
554       return NULL;
555    }   
556    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
557    fContainers->Add(cont);
558    switch (type) {
559       case kInputContainer:
560          fInputs->Add(cont);
561          break;
562       case kOutputContainer:
563          fOutputs->Add(cont);
564          if (filename && strlen(filename)) {
565             cont->SetFileName(filename);
566             cont->SetDataOwned(kFALSE);  // data owned by the file
567          }   
568          break;
569       case kExchangeContainer:
570          break;   
571    }
572    return cont;
573 }
574          
575 //______________________________________________________________________________
576 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
577                                         AliAnalysisDataContainer *cont)
578 {
579 // Connect input of an existing task to a data container.
580    if (!fTasks->FindObject(task)) {
581       AddTask(task);
582       Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
583    } 
584    Bool_t connected = task->ConnectInput(islot, cont);
585    return connected;
586 }   
587
588 //______________________________________________________________________________
589 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
590                                         AliAnalysisDataContainer *cont)
591 {
592 // Connect output of an existing task to a data container.
593    if (!fTasks->FindObject(task)) {
594       AddTask(task);
595       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
596    } 
597    Bool_t connected = task->ConnectOutput(islot, cont);
598    return connected;
599 }   
600                                
601 //______________________________________________________________________________
602 void AliAnalysisManager::CleanContainers()
603 {
604 // Clean data from all containers that have already finished all client tasks.
605    TIter next(fContainers);
606    AliAnalysisDataContainer *cont;
607    while ((cont=(AliAnalysisDataContainer *)next())) {
608       if (cont->IsOwnedData() && 
609           cont->IsDataReady() && 
610           cont->ClientsExecuted()) cont->DeleteData();
611    }
612 }
613
614 //______________________________________________________________________________
615 Bool_t AliAnalysisManager::InitAnalysis()
616 {
617 // Initialization of analysis chain of tasks. Should be called after all tasks
618 // and data containers are properly connected
619    // Check for input/output containers
620    fInitOK = kFALSE;
621    // Check for top tasks (depending only on input data containers)
622    if (!fTasks->First()) {
623       Error("InitAnalysis", "Analysis has no tasks !");
624       return kFALSE;
625    }   
626    TIter next(fTasks);
627    AliAnalysisTask *task;
628    AliAnalysisDataContainer *cont;
629    Int_t ntop = 0;
630    Int_t nzombies = 0;
631    Bool_t iszombie = kFALSE;
632    Bool_t istop = kTRUE;
633    Int_t i;
634    while ((task=(AliAnalysisTask*)next())) {
635       istop = kTRUE;
636       iszombie = kFALSE;
637       Int_t ninputs = task->GetNinputs();
638       for (i=0; i<ninputs; i++) {
639          cont = task->GetInputSlot(i)->GetContainer();
640          if (!cont) {
641             if (!iszombie) {
642                task->SetZombie();
643                fZombies->Add(task);
644                nzombies++;
645                iszombie = kTRUE;
646             }   
647             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
648                   i, task->GetName()); 
649          }
650          if (iszombie) continue;
651          // Check if cont is an input container
652          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
653          // Connect to parent task
654       }
655       if (istop) {
656          ntop++;
657          fTopTasks->Add(task);
658       }
659    }
660    if (!ntop) {
661       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
662       return kFALSE;
663    }                        
664    // Check now if there are orphan tasks
665    for (i=0; i<ntop; i++) {
666       task = (AliAnalysisTask*)fTopTasks->At(i);
667       task->SetUsed();
668    }
669    Int_t norphans = 0;
670    next.Reset();
671    while ((task=(AliAnalysisTask*)next())) {
672       if (!task->IsUsed()) {
673          norphans++;
674          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
675       }   
676    }          
677    // Check the task hierarchy (no parent task should depend on data provided
678    // by a daughter task)
679    for (i=0; i<ntop; i++) {
680       task = (AliAnalysisTask*)fTopTasks->At(i);
681       if (task->CheckCircularDeps()) {
682          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
683          PrintStatus("dep");
684          return kFALSE;
685       }   
686    }
687    // Check that all containers feeding post-event loop tasks are in the outputs list
688    TIter nextcont(fContainers); // loop over all containers
689    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
690       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
691          if (cont->HasConsumers()) {
692          // Check if one of the consumers is post event loop
693             TIter nextconsumer(cont->GetConsumers());
694             while ((task=(AliAnalysisTask*)nextconsumer())) {
695                if (task->IsPostEventLoop()) {
696                   fOutputs->Add(cont);
697                   break;
698                }
699             }
700          }
701       }
702    }   
703    fInitOK = kTRUE;
704    return kTRUE;
705 }   
706
707 //______________________________________________________________________________
708 void AliAnalysisManager::PrintStatus(Option_t *option) const
709 {
710 // Print task hierarchy.
711    if (!fInitOK) {
712       Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
713       return;
714    }   
715    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
716    if (getsysInfo)
717       Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
718    TIter next(fTopTasks);
719    AliAnalysisTask *task;
720    while ((task=(AliAnalysisTask*)next()))
721       task->PrintTask(option);
722 }
723
724 //______________________________________________________________________________
725 void AliAnalysisManager::ResetAnalysis()
726 {
727 // Reset all execution flags and clean containers.
728    CleanContainers();
729 }
730
731 //______________________________________________________________________________
732 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
733 {
734 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
735 // Process nentries starting from firstentry
736    if (!fInitOK) {
737       Error("StartAnalysis","Analysis manager was not initialized !");
738       return;
739    }
740    if (fDebug>1) {
741       cout << "StartAnalysis: " << GetName() << endl;   
742    }   
743    TString anaType = type;
744    anaType.ToLower();
745    fMode = kLocalAnalysis;
746    if (tree) {
747       if (anaType.Contains("proof"))     fMode = kProofAnalysis;
748       else if (anaType.Contains("grid")) fMode = kGridAnalysis;
749    }
750    if (fMode == kGridAnalysis) {
751       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
752       fMode = kLocalAnalysis;
753    }
754    char line[128];
755    SetEventLoop(kFALSE);
756    // Disable all branches if requested and set event loop mode
757    if (tree) {
758       if (TestBit(kDisableBranches)) {
759          printf("Disabling all branches...\n");
760 //         tree->SetBranchStatus("*",0); // not yet working
761       }   
762       SetEventLoop(kTRUE);
763    }   
764
765    TChain *chain = 0;
766    TString ttype = "TTree";
767    if (tree->IsA() == TChain::Class()) {
768       chain = (TChain*)tree;
769       ttype = "TChain";
770    }   
771
772    // Initialize locally all tasks
773    TIter next(fTasks);
774    AliAnalysisTask *task;
775    while ((task=(AliAnalysisTask*)next())) {
776       task->LocalInit();
777    }
778    
779    switch (fMode) {
780       case kLocalAnalysis:
781          if (!tree) {
782             TIter next(fTasks);
783             AliAnalysisTask *task;
784             // Call CreateOutputObjects for all tasks
785             while ((task=(AliAnalysisTask*)next())) {
786                TDirectory *curdir = gDirectory;
787                task->CreateOutputObjects();
788                if (curdir) curdir->cd();
789             }   
790             ExecAnalysis();
791             Terminate();
792             return;
793          } 
794          // Run tree-based analysis via AliAnalysisSelector  
795          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
796          sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
797          gROOT->ProcessLine(line);
798          sprintf(line, "((%s*)0x%lx)->Process(selector, \"\",%lld, %lld);",ttype.Data(),(ULong_t)tree, nentries, firstentry);
799          gROOT->ProcessLine(line);
800          break;
801       case kProofAnalysis:
802          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
803             printf("StartAnalysis: no PROOF!!!\n");
804             return;
805          }   
806          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
807          gROOT->ProcessLine(line);
808          if (chain) {
809             chain->SetProof();
810             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
811             chain->Process("AliAnalysisSelector", "", nentries, firstentry);
812          } else {
813             printf("StartAnalysis: no chain\n");
814             return;
815          }      
816          break;
817       case kGridAnalysis:
818          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
819    }   
820 }   
821
822 //______________________________________________________________________________
823 void AliAnalysisManager::ExecAnalysis(Option_t *option)
824 {
825 // Execute analysis.
826    static Long64_t ncalls = 0;
827    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
828    if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
829    ncalls++;
830    if (!fInitOK) {
831      Error("ExecAnalysis", "Analysis manager was not initialized !");
832       return;
833    }   
834    AliAnalysisTask *task;
835    // Check if the top tree is active.
836    if (fTree) {
837       TIter next(fTasks);
838    // De-activate all tasks
839       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
840       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
841       if (!cont) {
842               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
843          return;
844       }   
845       cont->SetData(fTree); // This will notify all consumers
846       Long64_t entry = fTree->GetTree()->GetReadEntry();
847       
848 //
849 //    Call BeginEvent() for optional input/output and MC services 
850       if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
851       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
852       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
853 //
854 //    Execute the tasks
855       TIter next1(cont->GetConsumers());
856       while ((task=(AliAnalysisTask*)next1())) {
857          if (fDebug >1) {
858             cout << "    Executing task " << task->GetName() << endl;
859          }   
860          
861          task->ExecuteTask(option);
862       }
863 //
864 //    Call FinishEvent() for optional output and MC services 
865       if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
866       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
867       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
868       // Gather system information if requested
869       if (getsysInfo && ((ncalls%fNSysInfo)==0)) 
870          AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
871       return;
872    }   
873    // The event loop is not controlled by TSelector   
874 //
875 //  Call BeginEvent() for optional input/output and MC services 
876    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(-1);
877    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(-1);
878    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
879    TIter next2(fTopTasks);
880    while ((task=(AliAnalysisTask*)next2())) {
881       task->SetActive(kTRUE);
882       if (fDebug > 1) {
883          cout << "    Executing task " << task->GetName() << endl;
884       }   
885       task->ExecuteTask(option);
886    }   
887 //
888 // Call FinishEvent() for optional output and MC services 
889    if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
890    if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
891    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
892 }
893
894 //______________________________________________________________________________
895 void AliAnalysisManager::FinishAnalysis()
896 {
897 // Finish analysis.
898 }