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