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