]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
Fix length of water inside cooling pipes
[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         TString fileName(curfile->GetName());
262         fileName.ReplaceAll("AliESDs.root", "");
263         fMCtruthEventHandler->Notify(fileName.Data());
264       }
265    }
266    return kTRUE;
267 }    
268
269 //______________________________________________________________________________
270 Bool_t AliAnalysisManager::Process(Long64_t entry)
271 {
272   // The Process() function is called for each entry in the tree (or possibly
273   // keyed object in the case of PROOF) to be processed. The entry argument
274   // specifies which entry in the currently loaded tree is to be processed.
275   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
276   // to read either all or the required parts of the data. When processing
277   // keyed objects with PROOF, the object is already loaded and is available
278   // via the fObject pointer.
279   //
280   // This function should contain the "body" of the analysis. It can contain
281   // simple or elaborate selection criteria, run algorithms on the data
282   // of the event and typically fill histograms.
283
284   // WARNING when a selector is used with a TChain, you must use
285   //  the pointer to the current TTree to call GetEntry(entry).
286   //  The entry is always the local entry number in the current tree.
287   //  Assuming that fChain is the pointer to the TChain being processed,
288   //  use fChain->GetTree()->GetEntry(entry).
289    if (fDebug > 1) {
290       cout << "->AliAnalysisManager::Process()" << endl;
291    }
292
293    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent();
294    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
295    
296    GetEntry(entry);
297    ExecAnalysis();
298    if (fDebug > 1) {
299       cout << "<-AliAnalysisManager::Process()" << endl;
300    }
301    return kTRUE;
302 }
303
304 //______________________________________________________________________________
305 void AliAnalysisManager::PackOutput(TList *target)
306 {
307   // Pack all output data containers in the output list. Called at SlaveTerminate
308   // stage in PROOF case for each slave.
309    if (fDebug > 1) {
310       cout << "->AliAnalysisManager::PackOutput()" << endl;
311    }   
312    if (!target) {
313       Error("PackOutput", "No target. Aborting.");
314       return;
315    }
316
317    if (fOutputEventHandler)  fOutputEventHandler ->Terminate();
318    if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
319    
320    if (fMode == kProofAnalysis) {
321       TIter next(fOutputs);
322       AliAnalysisDataContainer *output;
323       while ((output=(AliAnalysisDataContainer*)next())) {
324          if (fDebug > 1) printf("   Packing container %s...\n", output->GetName());
325          if (output->GetData()) target->Add(output->ExportData());
326       }
327    } 
328    if (fDebug > 1) {
329       printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
330    }
331 }
332
333 //______________________________________________________________________________
334 void AliAnalysisManager::ImportWrappers(TList *source)
335 {
336 // Import data in output containers from wrappers coming in source.
337    if (fDebug > 1) {
338       cout << "->AliAnalysisManager::ImportWrappers()" << endl;
339    }   
340    TIter next(fOutputs);
341    AliAnalysisDataContainer *cont;
342    AliAnalysisDataWrapper   *wrap;
343    Int_t icont = 0;
344    while ((cont=(AliAnalysisDataContainer*)next())) {
345       wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
346       if (!wrap && fDebug>1) {
347          printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
348          continue;
349       }
350       icont++;
351       if (fDebug > 1) printf("   Importing data for container %s\n", wrap->GetName());
352       if (cont->GetFileName()) printf("    -> %s\n", cont->GetFileName());
353       cont->ImportData(wrap);
354    }         
355    if (fDebug > 1) {
356       cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
357    }   
358 }
359
360 //______________________________________________________________________________
361 void AliAnalysisManager::UnpackOutput(TList *source)
362 {
363   // Called by AliAnalysisSelector::Terminate. Output containers should
364   // be in source in the same order as in fOutputs.
365    if (fDebug > 1) {
366       cout << "->AliAnalysisManager::UnpackOutput()" << endl;
367    }   
368    if (!source) {
369       Error("UnpackOutput", "No target. Aborting.");
370       return;
371    }
372    if (fDebug > 1) {
373       printf("   Source list contains %d containers\n", source->GetSize());
374    }   
375
376    if (fMode == kProofAnalysis) ImportWrappers(source);
377
378    TIter next(fOutputs);
379    AliAnalysisDataContainer *output;
380    while ((output=(AliAnalysisDataContainer*)next())) {
381       if (!output->GetData()) continue;
382       // Check if there are client tasks that run post event loop
383       if (output->HasConsumers()) {
384          // Disable event loop semaphore
385          output->SetPostEventLoop(kTRUE);
386          TObjArray *list = output->GetConsumers();
387          Int_t ncons = list->GetEntriesFast();
388          for (Int_t i=0; i<ncons; i++) {
389             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
390             task->CheckNotify(kTRUE);
391             // If task is active, execute it
392             if (task->IsPostEventLoop() && task->IsActive()) {
393                if (fDebug > 1) {
394                   cout << "== Executing post event loop task " << task->GetName() << endl;
395                }                  
396                task->ExecuteTask();
397             }   
398          }
399       }   
400       // Check if the output need to be written to a file.
401       const char *filename = output->GetFileName();
402       if (!(strcmp(filename, "default"))) {
403           if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
404       }
405       
406       if (!filename || !strlen(filename)) continue;
407       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
408       if (file) file->cd();
409       else      file = new TFile(filename, "RECREATE");
410       if (file->IsZombie()) continue;
411       // Reparent data to this file
412       TMethodCall callEnv;
413       if (output->GetData()->IsA())
414          callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
415       if (callEnv.IsValid()) {
416          callEnv.SetParam((Long_t) file);
417          callEnv.Execute(output->GetData());
418       }
419       output->GetData()->Write();
420    }
421    if (fDebug > 1) {
422       cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
423    }   
424 }
425
426 //______________________________________________________________________________
427 void AliAnalysisManager::Terminate()
428 {
429   // The Terminate() function is the last function to be called during
430   // a query. It always runs on the client, it can be used to present
431   // the results graphically.
432    if (fDebug > 1) {
433       cout << "->AliAnalysisManager::Terminate()" << endl;
434    }   
435    AliAnalysisTask *task;
436    TIter next(fTasks);
437    // Call Terminate() for tasks
438    while ((task=(AliAnalysisTask*)next())) task->Terminate();
439    if (fDebug > 1) {
440       cout << "<-AliAnalysisManager::Terminate()" << endl;
441    }   
442    //
443    if (fOutputEventHandler) fOutputEventHandler->TerminateIO();
444 }
445
446 //______________________________________________________________________________
447 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
448 {
449 // Adds a user task to the global list of tasks.
450    task->SetActive(kFALSE);
451    fTasks->Add(task);
452 }  
453
454 //______________________________________________________________________________
455 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
456 {
457 // Retreive task by name.
458    if (!fTasks) return NULL;
459    return (AliAnalysisTask*)fTasks->FindObject(name);
460 }
461
462 //______________________________________________________________________________
463 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
464                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
465 {
466 // Create a data container of a certain type. Types can be:
467 //   kExchangeContainer  = 0, used to exchange date between tasks
468 //   kInputContainer   = 1, used to store input data
469 //   kOutputContainer  = 2, used for posting results
470    if (fContainers->FindObject(name)) {
471       Error("CreateContainer","A container named %s already defined !\n",name);
472       return NULL;
473    }   
474    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
475    fContainers->Add(cont);
476    switch (type) {
477       case kInputContainer:
478          fInputs->Add(cont);
479          break;
480       case kOutputContainer:
481          fOutputs->Add(cont);
482          if (filename && strlen(filename)) cont->SetFileName(filename);
483          break;
484       case kExchangeContainer:
485          break;   
486    }
487    return cont;
488 }
489          
490 //______________________________________________________________________________
491 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
492                                         AliAnalysisDataContainer *cont)
493 {
494 // Connect input of an existing task to a data container.
495    if (!fTasks->FindObject(task)) {
496       AddTask(task);
497       Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
498    } 
499    Bool_t connected = task->ConnectInput(islot, cont);
500    return connected;
501 }   
502
503 //______________________________________________________________________________
504 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
505                                         AliAnalysisDataContainer *cont)
506 {
507 // Connect output of an existing task to a data container.
508    if (!fTasks->FindObject(task)) {
509       AddTask(task);
510       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
511    } 
512    Bool_t connected = task->ConnectOutput(islot, cont);
513    return connected;
514 }   
515                                
516 //______________________________________________________________________________
517 void AliAnalysisManager::CleanContainers()
518 {
519 // Clean data from all containers that have already finished all client tasks.
520    TIter next(fContainers);
521    AliAnalysisDataContainer *cont;
522    while ((cont=(AliAnalysisDataContainer *)next())) {
523       if (cont->IsOwnedData() && 
524           cont->IsDataReady() && 
525           cont->ClientsExecuted()) cont->DeleteData();
526    }
527 }
528
529 //______________________________________________________________________________
530 Bool_t AliAnalysisManager::InitAnalysis()
531 {
532 // Initialization of analysis chain of tasks. Should be called after all tasks
533 // and data containers are properly connected
534    // Check for input/output containers
535    fInitOK = kFALSE;
536    // Check for top tasks (depending only on input data containers)
537    if (!fTasks->First()) {
538       Error("InitAnalysis", "Analysis has no tasks !");
539       return kFALSE;
540    }   
541    TIter next(fTasks);
542    AliAnalysisTask *task;
543    AliAnalysisDataContainer *cont;
544    Int_t ntop = 0;
545    Int_t nzombies = 0;
546    Bool_t iszombie = kFALSE;
547    Bool_t istop = kTRUE;
548    Int_t i;
549    while ((task=(AliAnalysisTask*)next())) {
550       istop = kTRUE;
551       iszombie = kFALSE;
552       Int_t ninputs = task->GetNinputs();
553       for (i=0; i<ninputs; i++) {
554          cont = task->GetInputSlot(i)->GetContainer();
555          if (!cont) {
556             if (!iszombie) {
557                task->SetZombie();
558                fZombies->Add(task);
559                nzombies++;
560                iszombie = kTRUE;
561             }   
562             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
563                   i, task->GetName()); 
564          }
565          if (iszombie) continue;
566          // Check if cont is an input container
567          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
568          // Connect to parent task
569       }
570       if (istop) {
571          ntop++;
572          fTopTasks->Add(task);
573       }
574    }
575    if (!ntop) {
576       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
577       return kFALSE;
578    }                        
579    // Check now if there are orphan tasks
580    for (i=0; i<ntop; i++) {
581       task = (AliAnalysisTask*)fTopTasks->At(i);
582       task->SetUsed();
583    }
584    Int_t norphans = 0;
585    next.Reset();
586    while ((task=(AliAnalysisTask*)next())) {
587       if (!task->IsUsed()) {
588          norphans++;
589          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
590       }   
591    }          
592    // Check the task hierarchy (no parent task should depend on data provided
593    // by a daughter task)
594    for (i=0; i<ntop; i++) {
595       task = (AliAnalysisTask*)fTopTasks->At(i);
596       if (task->CheckCircularDeps()) {
597          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
598          PrintStatus("dep");
599          return kFALSE;
600       }   
601    }
602    // Check that all containers feeding post-event loop tasks are in the outputs list
603    TIter nextcont(fContainers); // loop over all containers
604    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
605       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
606          if (cont->HasConsumers()) {
607          // Check if one of the consumers is post event loop
608             TIter nextconsumer(cont->GetConsumers());
609             while ((task=(AliAnalysisTask*)nextconsumer())) {
610                if (task->IsPostEventLoop()) {
611                   fOutputs->Add(cont);
612                   break;
613                }
614             }
615          }
616       }
617    }   
618    fInitOK = kTRUE;
619    return kTRUE;
620 }   
621
622 //______________________________________________________________________________
623 void AliAnalysisManager::PrintStatus(Option_t *option) const
624 {
625 // Print task hierarchy.
626    TIter next(fTopTasks);
627    AliAnalysisTask *task;
628    while ((task=(AliAnalysisTask*)next()))
629       task->PrintTask(option);
630 }
631
632 //______________________________________________________________________________
633 void AliAnalysisManager::ResetAnalysis()
634 {
635 // Reset all execution flags and clean containers.
636    CleanContainers();
637 }
638
639 //______________________________________________________________________________
640 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree)
641 {
642 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
643    if (!fInitOK) {
644       Error("StartAnalysis","Analysis manager was not initialized !");
645       return;
646    }
647    if (fDebug>1) {
648       cout << "StartAnalysis: " << GetName() << endl;   
649    }   
650    TString anaType = type;
651    anaType.ToLower();
652    fMode = kLocalAnalysis;
653    if (tree) {
654       if (anaType.Contains("proof"))     fMode = kProofAnalysis;
655       else if (anaType.Contains("grid")) fMode = kGridAnalysis;
656    }
657    if (fMode == kGridAnalysis) {
658       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
659       fMode = kLocalAnalysis;
660    }
661    char line[128];
662    SetEventLoop(kFALSE);
663    // Disable all branches if requested and set event loop mode
664    if (tree) {
665       if (TestBit(kDisableBranches)) {
666          printf("Disabling all branches...\n");
667 //         tree->SetBranchStatus("*",0); // not yet working
668       }   
669       SetEventLoop(kTRUE);
670    }   
671
672    TChain *chain = dynamic_cast<TChain*>(tree);
673
674    // Initialize locally all tasks
675    TIter next(fTasks);
676    AliAnalysisTask *task;
677    while ((task=(AliAnalysisTask*)next())) {
678       task->LocalInit();
679    }
680    
681    switch (fMode) {
682       case kLocalAnalysis:
683          if (!tree) {
684             TIter next(fTasks);
685             AliAnalysisTask *task;
686             // Call CreateOutputObjects for all tasks
687             while ((task=(AliAnalysisTask*)next())) {
688                TDirectory *curdir = gDirectory;
689                task->CreateOutputObjects();
690                if (curdir) curdir->cd();
691             }   
692             ExecAnalysis();
693             Terminate();
694             return;
695          } 
696          // Run tree-based analysis via AliAnalysisSelector  
697 //         gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+");
698          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
699          sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
700          gROOT->ProcessLine(line);
701          sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree);
702          gROOT->ProcessLine(line);
703          break;
704       case kProofAnalysis:
705          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
706             printf("StartAnalysis: no PROOF!!!\n");
707             return;
708          }   
709          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
710          gROOT->ProcessLine(line);
711          if (chain) {
712             chain->SetProof();
713             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
714             chain->Process("AliAnalysisSelector");
715          } else {
716             printf("StartAnalysis: no chain\n");
717             return;
718          }      
719          break;
720       case kGridAnalysis:
721          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
722    }   
723 }   
724
725 //______________________________________________________________________________
726 void AliAnalysisManager::ExecAnalysis(Option_t *option)
727 {
728 // Execute analysis.
729    if (!fInitOK) {
730      Error("ExecAnalysis", "Analysis manager was not initialized !");
731       return;
732    }   
733    AliAnalysisTask *task;
734    // Check if the top tree is active.
735    if (fTree) {
736       if (fDebug>1) {
737          printf("AliAnalysisManager::ExecAnalysis\n");
738       }   
739       TIter next(fTasks);
740    // De-activate all tasks
741       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
742       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
743       if (!cont) {
744               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
745          return;
746       }   
747       cont->SetData(fTree); // This will notify all consumers
748 //
749 //    Call BeginEvent() for optional output and MC services 
750       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent();
751       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
752 //
753 //    Execute the tasks
754       TIter next1(cont->GetConsumers());
755       while ((task=(AliAnalysisTask*)next1())) {
756          if (fDebug >1) {
757             cout << "    Executing task " << task->GetName() << endl;
758          }   
759          
760          task->ExecuteTask(option);
761       }
762 //
763 //    Call FinishEvent() for optional output and MC services 
764       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
765       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
766 //
767       return;
768    }   
769    // The event loop is not controlled by TSelector   
770 //
771 //  Call BeginEvent() for optional output and MC services 
772    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent();
773    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent();
774    TIter next2(fTopTasks);
775    while ((task=(AliAnalysisTask*)next2())) {
776       task->SetActive(kTRUE);
777       if (fDebug > 1) {
778          cout << "    Executing task " << task->GetName() << endl;
779       }   
780       task->ExecuteTask(option);
781    }   
782 //
783 // Call FinishEvent() for optional output and MC services 
784    if (fOutputEventHandler)  fOutputEventHandler->FinishEvent();
785    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
786 }
787
788 //______________________________________________________________________________
789 void AliAnalysisManager::FinishAnalysis()
790 {
791 // Finish analysis.
792 }