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