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