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