- Introduced merging via TProofFile waiting for a full API to handle the output file.
[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 #include <TCanvas.h>
37
38 #include "AliAnalysisTask.h"
39 #include "AliAnalysisDataContainer.h"
40 #include "AliAnalysisDataSlot.h"
41 #include "AliVEventHandler.h"
42 #include "AliSysInfo.h"
43 #include "AliAnalysisManager.h"
44
45 ClassImp(AliAnalysisManager)
46
47 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
48
49 //______________________________________________________________________________
50 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
51                    :TNamed(name,title),
52                     fTree(NULL),
53                     fInputEventHandler(NULL),
54                     fOutputEventHandler(NULL),
55                     fMCtruthEventHandler(NULL),
56                     fCurrentEntry(-1),
57                     fNSysInfo(0),
58                     fMode(kLocalAnalysis),
59                     fInitOK(kFALSE),
60                     fDebug(0),
61                     fTasks(NULL),
62                     fTopTasks(NULL),
63                     fZombies(NULL),
64                     fContainers(NULL),
65                     fInputs(NULL),
66                     fOutputs(NULL)
67 {
68 // Default constructor.
69    fgAnalysisManager = this;
70    fTasks      = new TObjArray();
71    fTopTasks   = new TObjArray();
72    fZombies    = new TObjArray();
73    fContainers = new TObjArray();
74    fInputs     = new TObjArray();
75    fOutputs    = new TObjArray();
76    SetEventLoop(kTRUE);
77 }
78
79 //______________________________________________________________________________
80 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
81                    :TNamed(other),
82                     fTree(NULL),
83                     fInputEventHandler(NULL),
84                     fOutputEventHandler(NULL),
85                     fMCtruthEventHandler(NULL),
86                     fCurrentEntry(-1),
87                     fNSysInfo(0),
88                     fMode(other.fMode),
89                     fInitOK(other.fInitOK),
90                     fDebug(other.fDebug),
91                     fTasks(NULL),
92                     fTopTasks(NULL),
93                     fZombies(NULL),
94                     fContainers(NULL),
95                     fInputs(NULL),
96                     fOutputs(NULL)
97 {
98 // Copy constructor.
99    fTasks      = new TObjArray(*other.fTasks);
100    fTopTasks   = new TObjArray(*other.fTopTasks);
101    fZombies    = new TObjArray(*other.fZombies);
102    fContainers = new TObjArray(*other.fContainers);
103    fInputs     = new TObjArray(*other.fInputs);
104    fOutputs    = new TObjArray(*other.fOutputs);
105    fgAnalysisManager = this;
106 }
107    
108 //______________________________________________________________________________
109 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
110 {
111 // Assignment
112    if (&other != this) {
113       TNamed::operator=(other);
114       fInputEventHandler   = other.fInputEventHandler;
115       fOutputEventHandler  = other.fOutputEventHandler;
116       fMCtruthEventHandler = other.fMCtruthEventHandler;
117       fTree       = NULL;
118       fCurrentEntry = -1;
119       fNSysInfo   = other.fNSysInfo;
120       fMode       = other.fMode;
121       fInitOK     = other.fInitOK;
122       fDebug      = other.fDebug;
123       fTasks      = new TObjArray(*other.fTasks);
124       fTopTasks   = new TObjArray(*other.fTopTasks);
125       fZombies    = new TObjArray(*other.fZombies);
126       fContainers = new TObjArray(*other.fContainers);
127       fInputs     = new TObjArray(*other.fInputs);
128       fOutputs    = new TObjArray(*other.fOutputs);
129       fgAnalysisManager = this;
130    }
131    return *this;
132 }
133
134 //______________________________________________________________________________
135 AliAnalysisManager::~AliAnalysisManager()
136 {
137 // Destructor.
138    if (fTasks) {fTasks->Delete(); delete fTasks;}
139    if (fTopTasks) delete fTopTasks;
140    if (fZombies) delete fZombies;
141    if (fContainers) {fContainers->Delete(); delete fContainers;}
142    if (fInputs) delete fInputs;
143    if (fOutputs) delete fOutputs;
144    if (fgAnalysisManager==this) fgAnalysisManager = NULL;
145 }
146
147 //______________________________________________________________________________
148 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
149 {
150 // Read one entry of the tree or a whole branch.
151    if (fDebug > 1) {
152       cout << "== AliAnalysisManager::GetEntry()" << endl;
153    }   
154    fCurrentEntry = entry;
155    return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
156 }
157    
158 //______________________________________________________________________________
159 void AliAnalysisManager::Init(TTree *tree)
160 {
161   // The Init() function is called when the selector needs to initialize
162   // a new tree or chain. Typically here the branch addresses of the tree
163   // will be set. It is normaly not necessary to make changes to the
164   // generated code, but the routine can be extended by the user if needed.
165   // Init() will be called many times when running with PROOF.
166    if (!tree) return;
167    if (fDebug > 1) {
168       printf("->AliAnalysisManager::InitTree(%s)\n", tree->GetName());
169    }
170
171    // Call InitTree of EventHandler
172    if (fOutputEventHandler) {
173       if (fMode == kProofAnalysis) {
174          fOutputEventHandler->Init(0x0, "proof");
175       } else {
176          fOutputEventHandler->Init(0x0, "local");
177       }
178    }
179
180    if (fInputEventHandler) {
181       if (fMode == kProofAnalysis) {
182          fInputEventHandler->Init(tree, "proof");
183       } else {
184          fInputEventHandler->Init(tree, "local");
185       }
186    } else {
187       // If no input event handler we need to get the tree once
188       // for the chain
189       if(!tree->GetTree()) tree->LoadTree(0);
190    }
191    
192
193    if (fMCtruthEventHandler) {
194       if (fMode == kProofAnalysis) {
195          fMCtruthEventHandler->Init(0x0, "proof");
196       } else {
197          fMCtruthEventHandler->Init(0x0, "local");
198       }
199    }
200
201    if (!fInitOK) InitAnalysis();
202    if (!fInitOK) return;
203    fTree = tree;
204    AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
205    if (!top) {
206       cout<<"Error: No top input container !" <<endl;
207       return;
208    }
209    top->SetData(tree);
210    if (fDebug > 1) {
211       printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
212    }
213 }
214
215 //______________________________________________________________________________
216 void AliAnalysisManager::SlaveBegin(TTree *tree)
217 {
218   // The SlaveBegin() function is called after the Begin() function.
219   // When running with PROOF SlaveBegin() is called on each slave server.
220   // The tree argument is deprecated (on PROOF 0 is passed).
221    if (fDebug > 1) {
222       cout << "->AliAnalysisManager::SlaveBegin()" << endl;
223    }
224
225    // Call Init of EventHandler
226    if (fOutputEventHandler) {
227       if (fMode == kProofAnalysis) {
228          fOutputEventHandler->Init("proof");
229       } else {
230          fOutputEventHandler->Init("local");
231       }
232    }
233
234    if (fInputEventHandler) {
235       fInputEventHandler->SetInputTree(tree);
236       if (fMode == kProofAnalysis) {
237          fInputEventHandler->Init("proof");
238       } else {
239          fInputEventHandler->Init("local");
240       }
241    }
242
243    if (fMCtruthEventHandler) {
244       if (fMode == kProofAnalysis) {
245          fMCtruthEventHandler->Init("proof");
246       } else {
247          fMCtruthEventHandler->Init("local");
248       }
249    }
250
251    TIter next(fTasks);
252    AliAnalysisTask *task;
253    // Call CreateOutputObjects for all tasks
254    while ((task=(AliAnalysisTask*)next())) {
255       TDirectory *curdir = gDirectory;
256       task->CreateOutputObjects();
257       if (curdir) curdir->cd();
258    }
259
260    if (fDebug > 1) {
261       cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
262    }
263 }
264
265 //______________________________________________________________________________
266 Bool_t AliAnalysisManager::Notify()
267 {
268    // The Notify() function is called when a new file is opened. This
269    // can be either for a new TTree in a TChain or when when a new TTree
270    // is started when using PROOF. It is normaly not necessary to make changes
271    // to the generated code, but the routine can be extended by the
272    // user if needed. The return value is currently not used.
273     if (fTree) {
274         TFile *curfile = fTree->GetCurrentFile();
275         if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
276         TIter next(fTasks);
277         AliAnalysisTask *task;
278         // Call Notify for all tasks
279         while ((task=(AliAnalysisTask*)next())) 
280             task->Notify();
281         
282         // Call Notify of the event handlers
283         if (fInputEventHandler) {
284             fInputEventHandler->Notify(curfile->GetName());
285         }
286
287         if (fOutputEventHandler) {
288             fOutputEventHandler->Notify(curfile->GetName());
289         }
290
291         if (fMCtruthEventHandler) {
292             fMCtruthEventHandler->Notify(curfile->GetName());
293         }
294
295     }
296     return kTRUE;
297 }    
298
299 //______________________________________________________________________________
300 Bool_t AliAnalysisManager::Process(Long64_t entry)
301 {
302   // The Process() function is called for each entry in the tree (or possibly
303   // keyed object in the case of PROOF) to be processed. The entry argument
304   // specifies which entry in the currently loaded tree is to be processed.
305   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
306   // to read either all or the required parts of the data. When processing
307   // keyed objects with PROOF, the object is already loaded and is available
308   // via the fObject pointer.
309   //
310   // This function should contain the "body" of the analysis. It can contain
311   // simple or elaborate selection criteria, run algorithms on the data
312   // of the event and typically fill histograms.
313
314   // WARNING when a selector is used with a TChain, you must use
315   //  the pointer to the current TTree to call GetEntry(entry).
316   //  The entry is always the local entry number in the current tree.
317   //  Assuming that fChain is the pointer to the TChain being processed,
318   //  use fChain->GetTree()->GetEntry(entry).
319    if (fDebug > 1) {
320       cout << "->AliAnalysisManager::Process(" << entry << ")" << endl;
321    }
322    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
323    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
324    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
325    
326    GetEntry(entry);
327    ExecAnalysis();
328    if (fDebug > 1) {
329       cout << "<-AliAnalysisManager::Process()" << endl;
330    }
331    return kTRUE;
332 }
333
334 //______________________________________________________________________________
335 void AliAnalysisManager::PackOutput(TList *target)
336 {
337   // Pack all output data containers in the output list. Called at SlaveTerminate
338   // stage in PROOF case for each slave.
339    if (fDebug > 1) {
340       cout << "->AliAnalysisManager::PackOutput()" << endl;
341    }   
342    if (!target) {
343       Error("PackOutput", "No target. Aborting.");
344       return;
345    }
346    if (fInputEventHandler)   fInputEventHandler  ->Terminate();
347    if (fOutputEventHandler)  fOutputEventHandler ->Terminate();
348    if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
349    
350    if (fMode == kProofAnalysis) {
351       TIter next(fOutputs);
352       AliAnalysisDataContainer *output;
353       while ((output=(AliAnalysisDataContainer*)next())) {
354          if (output->GetData() && !output->IsSpecialOutput()) {
355             if (output->GetProducer()->IsPostEventLoop()) continue;
356
357             const char *filename = output->GetFileName();
358             if (!(strcmp(filename, "default"))) {
359                if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
360             }      
361             if (strlen(filename)) {
362                TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
363                TDirectory *opwd = gDirectory;
364                if (file) file->cd();
365                else      file = new TFile(filename, "RECREATE");
366                if (file->IsZombie()) continue;
367                // Clear file list to release object ownership to user.
368                // Save data to file, then close.
369                file->Clear();
370                output->GetData()->Write();
371                file->Close();
372                // Set null directory to histograms and trees.
373                TMethodCall callEnv;
374                if (output->GetData()->IsA())
375                   callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
376                if (callEnv.IsValid()) {
377                   callEnv.SetParam(Long_t(0));
378                   callEnv.Execute(output->GetData());
379                }
380                // Restore current directory
381                if (opwd) opwd->cd();
382             }   
383             AliAnalysisDataWrapper *wrap = output->ExportData();
384             // Output wrappers must delete data after merging (AG 13/11/07)
385             wrap->SetDeleteData(kTRUE);
386             if (fDebug > 1) printf("   Packing container %s...\n", output->GetName());
387             target->Add(wrap);
388          }   
389          // Special outputs files are closed and copied on the remote location
390          if (output->IsSpecialOutput() && strlen(output->GetFileName())) {
391             TDirectory *opwd = gDirectory;
392             TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(output->GetFileName());
393             if (!file) continue;
394             file->cd();
395             if (output->GetData()) output->GetData()->Write();
396             file->Close();
397             if (opwd) opwd->cd();
398             if (strlen(fSpecialOutputLocation.Data())) {
399                TString remote = fSpecialOutputLocation;
400                remote += "/";
401                Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
402                remote += Form("%s_%d_", gSystem->HostName(), gid);
403                remote += output->GetFileName();
404                TFile::Cp(output->GetFileName(), remote.Data());
405             } else {
406             // No special location specified-> use TProofFile as merging utility
407                char line[256];
408                sprintf(line, "((TList*)0x%lx)->Add(new TProofFile(\"%s\"));",
409                        (ULong_t)target, output->GetFileName());
410                gROOT->ProcessLine(line);
411             }
412          }      
413       }
414    } 
415    if (fDebug > 1) {
416       printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
417    }
418 }
419
420 //______________________________________________________________________________
421 void AliAnalysisManager::ImportWrappers(TList *source)
422 {
423 // Import data in output containers from wrappers coming in source.
424    if (fDebug > 1) {
425       cout << "->AliAnalysisManager::ImportWrappers()" << endl;
426    }   
427    TIter next(fOutputs);
428    AliAnalysisDataContainer *cont;
429    AliAnalysisDataWrapper   *wrap;
430    Int_t icont = 0;
431    while ((cont=(AliAnalysisDataContainer*)next())) {
432       if (cont->GetProducer()->IsPostEventLoop() ||
433           cont->IsSpecialOutput()) continue;
434       wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
435       if (!wrap && fDebug>1) {
436          printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
437          continue;
438       }
439       icont++;
440       if (fDebug > 1) printf("   Importing data for container %s\n", wrap->GetName());
441       if (cont->GetFileName()) printf("    -> %s\n", cont->GetFileName());
442       cont->ImportData(wrap);
443    }         
444    if (fDebug > 1) {
445       cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
446    }   
447 }
448
449 //______________________________________________________________________________
450 void AliAnalysisManager::UnpackOutput(TList *source)
451 {
452   // Called by AliAnalysisSelector::Terminate only on the client.
453    if (fDebug > 1) {
454       cout << "->AliAnalysisManager::UnpackOutput()" << endl;
455    }   
456    if (!source) {
457       Error("UnpackOutput", "No target. Aborting.");
458       return;
459    }
460    if (fDebug > 1) {
461       printf("   Source list contains %d containers\n", source->GetSize());
462    }   
463
464    if (fMode == kProofAnalysis) ImportWrappers(source);
465
466    TIter next(fOutputs);
467    AliAnalysisDataContainer *output;
468    while ((output=(AliAnalysisDataContainer*)next())) {
469       if (!output->GetData()) continue;
470       // Check if there are client tasks that run post event loop
471       if (output->HasConsumers()) {
472          // Disable event loop semaphore
473          output->SetPostEventLoop(kTRUE);
474          TObjArray *list = output->GetConsumers();
475          Int_t ncons = list->GetEntriesFast();
476          for (Int_t i=0; i<ncons; i++) {
477             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
478             task->CheckNotify(kTRUE);
479             // If task is active, execute it
480             if (task->IsPostEventLoop() && task->IsActive()) {
481                if (fDebug > 1) {
482                   cout << "== Executing post event loop task " << task->GetName() << endl;
483                }                  
484                task->ExecuteTask();
485             }   
486          }
487       }   
488       // Check if the output need to be written to a file.
489       const char *filename = output->GetFileName();
490       if (!(strcmp(filename, "default"))) {
491           if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
492       }
493       
494       if (!filename || !strlen(filename)) continue;
495       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
496       TDirectory *opwd = gDirectory;
497       if (file) file->cd();
498       else      file = new TFile(filename, "RECREATE");
499       if (file->IsZombie()) continue;
500       // Reparent data to this file
501       TMethodCall callEnv;
502       if (output->GetData()->IsA())
503          callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
504       if (callEnv.IsValid()) {
505          callEnv.SetParam((Long_t) file);
506          callEnv.Execute(output->GetData());
507       }
508       output->GetData()->Write();
509       if (opwd) opwd->cd();
510    }
511    if (fDebug > 1) {
512       cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
513    }   
514 }
515
516 //______________________________________________________________________________
517 void AliAnalysisManager::Terminate()
518 {
519   // The Terminate() function is the last function to be called during
520   // a query. It always runs on the client, it can be used to present
521   // the results graphically.
522    if (fDebug > 1) {
523       cout << "->AliAnalysisManager::Terminate()" << endl;
524    }   
525    AliAnalysisTask *task;
526    TIter next(fTasks);
527    // Call Terminate() for tasks
528    while ((task=(AliAnalysisTask*)next())) task->Terminate();
529    if (fDebug > 1) {
530       cout << "<-AliAnalysisManager::Terminate()" << endl;
531    }   
532    //
533    if (fInputEventHandler)   fInputEventHandler  ->TerminateIO();
534    if (fOutputEventHandler)  fOutputEventHandler ->TerminateIO();
535    if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
536    TIter next1(fOutputs);
537    AliAnalysisDataContainer *output;
538    while ((output=(AliAnalysisDataContainer*)next1())) {
539       // Close all files at output
540       const char *filename = output->GetFileName();
541       if (!(strcmp(filename, "default"))) {
542          if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
543       }
544       
545       if (!filename || !strlen(filename)) continue;
546       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
547       if (!file || file->IsZombie()) continue;
548       file->Close();
549    }   
550
551    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
552    if (getsysInfo) {
553       TDirectory *cdir = gDirectory;
554       TFile f("syswatch.root", "RECREATE");
555       if (!f.IsZombie()) {
556          TTree *tree = AliSysInfo::MakeTree("syswatch.log");
557          tree->SetMarkerStyle(kCircle);
558          tree->SetMarkerColor(kBlue);
559          tree->SetMarkerSize(0.5);
560          if (!gROOT->IsBatch()) {
561             tree->SetAlias("event", "id0");
562             tree->SetAlias("memUSED", "mi.fMemUsed");
563             tree->SetAlias("userCPU", "pI.fCpuUser");
564             TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
565             c->Divide(2,1,0.01,0.01);
566             c->cd(1);
567             tree->Draw("memUSED:event","","", 1234567890, 0);
568             c->cd(2);
569             tree->Draw("userCPU:event","","", 1234567890, 0);
570          }   
571          tree->Write();
572          f.Close();
573          delete tree;
574       }
575       if (cdir) cdir->cd();
576    }      
577 }
578
579 //______________________________________________________________________________
580 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
581 {
582 // Adds a user task to the global list of tasks.
583    task->SetActive(kFALSE);
584    fTasks->Add(task);
585 }  
586
587 //______________________________________________________________________________
588 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
589 {
590 // Retreive task by name.
591    if (!fTasks) return NULL;
592    return (AliAnalysisTask*)fTasks->FindObject(name);
593 }
594
595 //______________________________________________________________________________
596 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
597                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
598 {
599 // Create a data container of a certain type. Types can be:
600 //   kExchangeContainer  = 0, used to exchange date between tasks
601 //   kInputContainer   = 1, used to store input data
602 //   kOutputContainer  = 2, used for posting results
603    if (fContainers->FindObject(name)) {
604       Error("CreateContainer","A container named %s already defined !\n",name);
605       return NULL;
606    }   
607    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
608    fContainers->Add(cont);
609    switch (type) {
610       case kInputContainer:
611          fInputs->Add(cont);
612          break;
613       case kOutputContainer:
614          fOutputs->Add(cont);
615          if (filename && strlen(filename)) {
616             cont->SetFileName(filename);
617             cont->SetDataOwned(kFALSE);  // data owned by the file
618          }   
619          break;
620       case kExchangeContainer:
621          break;   
622    }
623    return cont;
624 }
625          
626 //______________________________________________________________________________
627 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
628                                         AliAnalysisDataContainer *cont)
629 {
630 // Connect input of an existing task to a data container.
631    if (!fTasks->FindObject(task)) {
632       AddTask(task);
633       Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
634    } 
635    Bool_t connected = task->ConnectInput(islot, cont);
636    return connected;
637 }   
638
639 //______________________________________________________________________________
640 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
641                                         AliAnalysisDataContainer *cont)
642 {
643 // Connect output of an existing task to a data container.
644    if (!fTasks->FindObject(task)) {
645       AddTask(task);
646       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
647    } 
648    Bool_t connected = task->ConnectOutput(islot, cont);
649    return connected;
650 }   
651                                
652 //______________________________________________________________________________
653 void AliAnalysisManager::CleanContainers()
654 {
655 // Clean data from all containers that have already finished all client tasks.
656    TIter next(fContainers);
657    AliAnalysisDataContainer *cont;
658    while ((cont=(AliAnalysisDataContainer *)next())) {
659       if (cont->IsOwnedData() && 
660           cont->IsDataReady() && 
661           cont->ClientsExecuted()) cont->DeleteData();
662    }
663 }
664
665 //______________________________________________________________________________
666 Bool_t AliAnalysisManager::InitAnalysis()
667 {
668 // Initialization of analysis chain of tasks. Should be called after all tasks
669 // and data containers are properly connected
670    // Check for input/output containers
671    fInitOK = kFALSE;
672    // Check for top tasks (depending only on input data containers)
673    if (!fTasks->First()) {
674       Error("InitAnalysis", "Analysis has no tasks !");
675       return kFALSE;
676    }   
677    TIter next(fTasks);
678    AliAnalysisTask *task;
679    AliAnalysisDataContainer *cont;
680    Int_t ntop = 0;
681    Int_t nzombies = 0;
682    Bool_t iszombie = kFALSE;
683    Bool_t istop = kTRUE;
684    Int_t i;
685    while ((task=(AliAnalysisTask*)next())) {
686       istop = kTRUE;
687       iszombie = kFALSE;
688       Int_t ninputs = task->GetNinputs();
689       for (i=0; i<ninputs; i++) {
690          cont = task->GetInputSlot(i)->GetContainer();
691          if (!cont) {
692             if (!iszombie) {
693                task->SetZombie();
694                fZombies->Add(task);
695                nzombies++;
696                iszombie = kTRUE;
697             }   
698             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
699                   i, task->GetName()); 
700          }
701          if (iszombie) continue;
702          // Check if cont is an input container
703          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
704          // Connect to parent task
705       }
706       if (istop) {
707          ntop++;
708          fTopTasks->Add(task);
709       }
710    }
711    if (!ntop) {
712       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
713       return kFALSE;
714    }                        
715    // Check now if there are orphan tasks
716    for (i=0; i<ntop; i++) {
717       task = (AliAnalysisTask*)fTopTasks->At(i);
718       task->SetUsed();
719    }
720    Int_t norphans = 0;
721    next.Reset();
722    while ((task=(AliAnalysisTask*)next())) {
723       if (!task->IsUsed()) {
724          norphans++;
725          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
726       }   
727    }          
728    // Check the task hierarchy (no parent task should depend on data provided
729    // by a daughter task)
730    for (i=0; i<ntop; i++) {
731       task = (AliAnalysisTask*)fTopTasks->At(i);
732       if (task->CheckCircularDeps()) {
733          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
734          PrintStatus("dep");
735          return kFALSE;
736       }   
737    }
738    // Check that all containers feeding post-event loop tasks are in the outputs list
739    TIter nextcont(fContainers); // loop over all containers
740    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
741       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
742          if (cont->HasConsumers()) {
743          // Check if one of the consumers is post event loop
744             TIter nextconsumer(cont->GetConsumers());
745             while ((task=(AliAnalysisTask*)nextconsumer())) {
746                if (task->IsPostEventLoop()) {
747                   fOutputs->Add(cont);
748                   break;
749                }
750             }
751          }
752       }
753    }   
754    fInitOK = kTRUE;
755    return kTRUE;
756 }   
757
758 //______________________________________________________________________________
759 void AliAnalysisManager::PrintStatus(Option_t *option) const
760 {
761 // Print task hierarchy.
762    if (!fInitOK) {
763       Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
764       return;
765    }   
766    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
767    if (getsysInfo)
768       Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
769    TIter next(fTopTasks);
770    AliAnalysisTask *task;
771    while ((task=(AliAnalysisTask*)next()))
772       task->PrintTask(option);
773 }
774
775 //______________________________________________________________________________
776 void AliAnalysisManager::ResetAnalysis()
777 {
778 // Reset all execution flags and clean containers.
779    CleanContainers();
780 }
781
782 //______________________________________________________________________________
783 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
784 {
785 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
786 // Process nentries starting from firstentry
787    if (!fInitOK) {
788       Error("StartAnalysis","Analysis manager was not initialized !");
789       return;
790    }
791    if (fDebug>1) {
792       cout << "StartAnalysis: " << GetName() << endl;   
793    }   
794    TString anaType = type;
795    anaType.ToLower();
796    fMode = kLocalAnalysis;
797    if (tree) {
798       if (anaType.Contains("proof"))     fMode = kProofAnalysis;
799       else if (anaType.Contains("grid")) fMode = kGridAnalysis;
800    }
801    if (fMode == kGridAnalysis) {
802       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
803       fMode = kLocalAnalysis;
804    }
805    char line[256];
806    SetEventLoop(kFALSE);
807    // Disable all branches if requested and set event loop mode
808    if (tree) {
809       if (TestBit(kDisableBranches)) {
810          printf("Disabling all branches...\n");
811 //         tree->SetBranchStatus("*",0); // not yet working
812       }   
813       SetEventLoop(kTRUE);
814    }   
815
816    TChain *chain = 0;
817    TString ttype = "TTree";
818    if (tree->IsA() == TChain::Class()) {
819       chain = (TChain*)tree;
820       if (!chain || !chain->GetListOfFiles()->First()) {
821          Error("StartAnalysis", "Cannot process null or empty chain...");
822          return;
823       }   
824       ttype = "TChain";
825    }   
826
827    // Initialize locally all tasks
828    TIter next(fTasks);
829    AliAnalysisTask *task;
830    while ((task=(AliAnalysisTask*)next())) {
831       task->LocalInit();
832    }
833    
834    switch (fMode) {
835       case kLocalAnalysis:
836          if (!tree) {
837             TIter next(fTasks);
838             AliAnalysisTask *task;
839             // Call CreateOutputObjects for all tasks
840             while ((task=(AliAnalysisTask*)next())) {
841                TDirectory *curdir = gDirectory;
842                task->CreateOutputObjects();
843                if (curdir) curdir->cd();
844             }   
845             ExecAnalysis();
846             Terminate();
847             return;
848          } 
849          // Run tree-based analysis via AliAnalysisSelector  
850          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
851          sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
852          gROOT->ProcessLine(line);
853          sprintf(line, "((%s*)0x%lx)->Process(selector, \"\",%lld, %lld);",ttype.Data(),(ULong_t)tree, nentries, firstentry);
854          gROOT->ProcessLine(line);
855          break;
856       case kProofAnalysis:
857          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
858             printf("StartAnalysis: no PROOF!!!\n");
859             return;
860          }   
861          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
862          gROOT->ProcessLine(line);
863          if (chain) {
864             chain->SetProof();
865             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
866             chain->Process("AliAnalysisSelector", "", nentries, firstentry);
867          } else {
868             printf("StartAnalysis: no chain\n");
869             return;
870          }      
871          break;
872       case kGridAnalysis:
873          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
874    }   
875 }   
876
877 //______________________________________________________________________________
878 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
879 {
880 // Start analysis for this manager on a given dataset. Analysis task can be: 
881 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
882    if (!fInitOK) {
883       Error("StartAnalysis","Analysis manager was not initialized !");
884       return;
885    }
886    if (fDebug>1) {
887       cout << "StartAnalysis: " << GetName() << endl;   
888    }   
889    TString anaType = type;
890    anaType.ToLower();
891    if (!anaType.Contains("proof")) {
892       Error("Cannot process datasets in %s mode. Try PROOF.", type);
893       return;
894    }   
895    fMode = kProofAnalysis;
896    char line[256];
897    SetEventLoop(kTRUE);
898    // Set the dataset flag
899    TObject::SetBit(kUseDataSet);
900    fTree = 0;
901
902    // Initialize locally all tasks
903    TIter next(fTasks);
904    AliAnalysisTask *task;
905    while ((task=(AliAnalysisTask*)next())) {
906       task->LocalInit();
907    }
908    
909    if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
910       printf("StartAnalysis: no PROOF!!!\n");
911       return;
912    }   
913    sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
914    gROOT->ProcessLine(line);
915    sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
916    if (!gROOT->ProcessLine(line)) {
917       Error("StartAnalysis", "Dataset %s not found", dataset);
918       return;
919    }   
920    sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
921            dataset, nentries, firstentry);
922    cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
923    gROOT->ProcessLine(line);
924 }   
925
926 //______________________________________________________________________________
927 void AliAnalysisManager::ExecAnalysis(Option_t *option)
928 {
929 // Execute analysis.
930    static Long64_t ncalls = 0;
931    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
932    if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
933    ncalls++;
934    if (!fInitOK) {
935      Error("ExecAnalysis", "Analysis manager was not initialized !");
936       return;
937    }   
938    AliAnalysisTask *task;
939    // Check if the top tree is active.
940    if (fTree) {
941       TIter next(fTasks);
942    // De-activate all tasks
943       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
944       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
945       if (!cont) {
946               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
947          return;
948       }   
949       cont->SetData(fTree); // This will notify all consumers
950       Long64_t entry = fTree->GetTree()->GetReadEntry();
951       
952 //
953 //    Call BeginEvent() for optional input/output and MC services 
954       if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
955       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
956       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
957 //
958 //    Execute the tasks
959 //      TIter next1(cont->GetConsumers());
960       TIter next1(fTopTasks);
961       while ((task=(AliAnalysisTask*)next1())) {
962          if (fDebug >1) {
963             cout << "    Executing task " << task->GetName() << endl;
964          }   
965          
966          task->ExecuteTask(option);
967       }
968 //
969 //    Call FinishEvent() for optional output and MC services 
970       if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
971       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
972       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
973       // Gather system information if requested
974       if (getsysInfo && ((ncalls%fNSysInfo)==0)) 
975          AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
976       return;
977    }   
978    // The event loop is not controlled by TSelector   
979 //
980 //  Call BeginEvent() for optional input/output and MC services 
981    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(-1);
982    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(-1);
983    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
984    TIter next2(fTopTasks);
985    while ((task=(AliAnalysisTask*)next2())) {
986       task->SetActive(kTRUE);
987       if (fDebug > 1) {
988          cout << "    Executing task " << task->GetName() << endl;
989       }   
990       task->ExecuteTask(option);
991    }   
992 //
993 // Call FinishEvent() for optional output and MC services 
994    if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
995    if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
996    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
997 }
998
999 //______________________________________________________________________________
1000 void AliAnalysisManager::FinishAnalysis()
1001 {
1002 // Finish analysis.
1003 }