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