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