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