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