]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
Fix for Coverity defect 21178
[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 //   AliAnalysisManager - 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 "AliAnalysisManager.h"
29
30 #include <cerrno>
31 #include <Riostream.h>
32 #include <TError.h>
33 #include <TMap.h>
34 #include <TClass.h>
35 #include <TFile.h>
36 #include <TTreeCache.h>
37 #include <TEnv.h>
38 #include <TMath.h>
39 #include <TH1.h>
40 #include <TMethodCall.h>
41 #include <TChain.h>
42 #include <TSystem.h>
43 #include <TROOT.h>
44 #include <TCanvas.h>
45 #include <TStopwatch.h>
46
47 #include "AliLog.h"
48 #include "AliAnalysisSelector.h"
49 #include "AliAnalysisGrid.h"
50 #include "AliAnalysisTask.h"
51 #include "AliAnalysisDataContainer.h"
52 #include "AliAnalysisDataSlot.h"
53 #include "AliVEventHandler.h"
54 #include "AliVEventPool.h"
55 #include "AliSysInfo.h"
56 #include "AliAnalysisStatistics.h"
57
58 using std::ofstream;
59 using std::ios;
60 using std::cout;
61 using std::endl;
62 ClassImp(AliAnalysisManager)
63
64 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
65 TString AliAnalysisManager::fgCommonFileName = "";
66 TString AliAnalysisManager::fgMacroNames = "";
67 Int_t AliAnalysisManager::fPBUpdateFreq = 1;
68
69 //______________________________________________________________________________
70 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
71                    :TNamed(name,title),
72                     fTree(0),
73                     fInputEventHandler(0),
74                     fOutputEventHandler(0),
75                     fMCtruthEventHandler(0),
76                     fEventPool(0),
77                     fCurrentEntry(-1),
78                     fNSysInfo(0),
79                     fMode(kLocalAnalysis),
80                     fInitOK(kFALSE),
81                     fMustClean(kFALSE),
82                     fIsRemote(kFALSE),
83                     fLocked(kFALSE),
84                     fMCLoop(kFALSE),
85                     fDebug(0),
86                     fSpecialOutputLocation(""), 
87                     fTasks(0),
88                     fTopTasks(0),
89                     fZombies(0),
90                     fContainers(0),
91                     fInputs(0),
92                     fOutputs(0),
93                     fParamCont(0),
94                     fDebugOptions(0),
95                     fFileDescriptors(new TObjArray()),
96                     fCurrentDescriptor(0),
97                     fCommonInput(0),
98                     fCommonOutput(0),
99                     fSelector(0),
100                     fGridHandler(0),
101                     fExtraFiles(),
102                     fFileInfoLog(),
103                     fAutoBranchHandling(kTRUE), 
104                     fAsyncReading(kFALSE), // default prefetching on
105                     fTable(),
106                     fRunFromPath(0),
107                     fNcalls(0),
108                     fMaxEntries(0),
109                     fCacheSize(100000000), // default 100 MB
110                     fStatisticsMsg(),
111                     fRequestedBranches(),
112                     fStatistics(0),
113                     fGlobals(0),
114                     fIOTimer(0),
115                     fCPUTimer(0),
116                     fInitTimer(0),
117                     fIOTime(0),
118                     fCPUTime(0),
119                     fInitTime(0)
120 {
121 // Default constructor.
122    fgAnalysisManager = this;
123    fgCommonFileName  = "AnalysisResults.root";
124    if (TClass::IsCallingNew() != TClass::kDummyNew) {
125      fTasks      = new TObjArray();
126      fTopTasks   = new TObjArray();
127      fZombies    = new TObjArray();
128      fContainers = new TObjArray();
129      fInputs     = new TObjArray();
130      fOutputs    = new TObjArray();
131      fParamCont  = new TObjArray();
132      fGlobals    = new TMap();
133    }
134    fIOTimer = new TStopwatch();
135    fCPUTimer = new TStopwatch();
136    fInitTimer = new TStopwatch();
137    SetEventLoop(kTRUE);
138 }
139
140 //______________________________________________________________________________
141 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
142                    :TNamed(other),
143                     fTree(NULL),
144                     fInputEventHandler(NULL),
145                     fOutputEventHandler(NULL),
146                     fMCtruthEventHandler(NULL),
147                     fEventPool(NULL),
148                     fCurrentEntry(-1),
149                     fNSysInfo(0),
150                     fMode(other.fMode),
151                     fInitOK(other.fInitOK),
152                     fMustClean(other.fMustClean),
153                     fIsRemote(other.fIsRemote),
154                     fLocked(other.fLocked),
155                     fMCLoop(other.fMCLoop),
156                     fDebug(other.fDebug),
157                     fSpecialOutputLocation(""), 
158                     fTasks(NULL),
159                     fTopTasks(NULL),
160                     fZombies(NULL),
161                     fContainers(NULL),
162                     fInputs(NULL),
163                     fOutputs(NULL),
164                     fParamCont(NULL),
165                     fDebugOptions(NULL),
166                     fFileDescriptors(new TObjArray()),
167                     fCurrentDescriptor(0),
168                     fCommonInput(NULL),
169                     fCommonOutput(NULL),
170                     fSelector(NULL),
171                     fGridHandler(NULL),
172                     fExtraFiles(other.fExtraFiles),
173                     fFileInfoLog(other.fFileInfoLog),
174                     fAutoBranchHandling(other.fAutoBranchHandling), 
175                     fAsyncReading(other.fAsyncReading),
176                     fTable(),
177                     fRunFromPath(0),
178                     fNcalls(other.fNcalls),
179                     fMaxEntries(other.fMaxEntries),
180                     fCacheSize(other.fCacheSize),
181                     fStatisticsMsg(other.fStatisticsMsg),
182                     fRequestedBranches(other.fRequestedBranches),
183                     fStatistics(other.fStatistics),
184                     fGlobals(other.fGlobals),
185                     fIOTimer(new TStopwatch()),
186                     fCPUTimer(new TStopwatch()),
187                     fInitTimer(new TStopwatch()),
188                     fIOTime(0),
189                     fCPUTime(0),
190                     fInitTime(0)
191 {
192 // Copy constructor.
193    fTasks      = new TObjArray(*other.fTasks);
194    fTopTasks   = new TObjArray(*other.fTopTasks);
195    fZombies    = new TObjArray(*other.fZombies);
196    fContainers = new TObjArray(*other.fContainers);
197    fInputs     = new TObjArray(*other.fInputs);
198    fOutputs    = new TObjArray(*other.fOutputs);
199    fParamCont  = new TObjArray(*other.fParamCont);
200    fgCommonFileName  = "AnalysisResults.root";
201    fgAnalysisManager = this;
202 }
203    
204 //______________________________________________________________________________
205 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
206 {
207 // Assignment
208    if (&other != this) {
209       TNamed::operator=(other);
210       fInputEventHandler   = other.fInputEventHandler;
211       fOutputEventHandler  = other.fOutputEventHandler;
212       fMCtruthEventHandler = other.fMCtruthEventHandler;
213       fEventPool           = other.fEventPool;
214       fTree       = NULL;
215       fCurrentEntry = -1;
216       fNSysInfo   = other.fNSysInfo;
217       fMode       = other.fMode;
218       fInitOK     = other.fInitOK;
219       fIsRemote   = other.fIsRemote;
220       fLocked     = other.fLocked;
221       fMCLoop     = other.fMCLoop;
222       fDebug      = other.fDebug;
223       fTasks      = new TObjArray(*other.fTasks);
224       fTopTasks   = new TObjArray(*other.fTopTasks);
225       fZombies    = new TObjArray(*other.fZombies);
226       fContainers = new TObjArray(*other.fContainers);
227       fInputs     = new TObjArray(*other.fInputs);
228       fOutputs    = new TObjArray(*other.fOutputs);
229       fParamCont  = new TObjArray(*other.fParamCont);
230       fDebugOptions = NULL;
231       fFileDescriptors = new TObjArray();
232       fCurrentDescriptor = 0;
233       fCommonInput = NULL;
234       fCommonOutput = NULL;
235       fSelector   = NULL;
236       fGridHandler = NULL;
237       fExtraFiles = other.fExtraFiles;
238       fFileInfoLog = other.fFileInfoLog;
239       fgCommonFileName = "AnalysisResults.root";
240       fgAnalysisManager = this;
241       fAutoBranchHandling = other.fAutoBranchHandling;
242       fAsyncReading = other.fAsyncReading;
243       fTable.Clear("nodelete");
244       fRunFromPath = other.fRunFromPath;
245       fNcalls     = other. fNcalls;
246       fMaxEntries = other.fMaxEntries;
247       fCacheSize = other.fCacheSize;
248       fStatisticsMsg = other.fStatisticsMsg;
249       fRequestedBranches = other.fRequestedBranches;
250       fStatistics = other.fStatistics;
251       fGlobals = new TMap();
252       fIOTimer = new TStopwatch();
253       fCPUTimer = new TStopwatch();
254       fInitTimer = new TStopwatch();
255       fIOTime = 0.;
256       fCPUTime = 0.;
257       fInitTime = 0.;
258    }
259    return *this;
260 }
261
262 //______________________________________________________________________________
263 AliAnalysisManager::~AliAnalysisManager()
264 {
265 // Destructor.
266    if (fTasks) {fTasks->Delete(); delete fTasks;}
267    if (fTopTasks) delete fTopTasks;
268    if (fZombies) delete fZombies;
269    if (fContainers) {fContainers->Delete(); delete fContainers;}
270    if (fInputs) delete fInputs;
271    if (fOutputs) delete fOutputs;
272    if (fParamCont) delete fParamCont;
273    if (fDebugOptions) delete fDebugOptions;
274    if (fGridHandler) delete fGridHandler;
275    if (fInputEventHandler) delete fInputEventHandler;
276    if (fOutputEventHandler) delete fOutputEventHandler;
277    if (fMCtruthEventHandler) delete fMCtruthEventHandler;
278    if (fEventPool) delete fEventPool;
279    if (fgAnalysisManager==this) fgAnalysisManager = NULL;
280    if (fGlobals) {fGlobals->DeleteAll(); delete fGlobals;}
281    if (fFileDescriptors) {fFileDescriptors->Delete(); delete fFileDescriptors;}
282    delete fIOTimer;
283    delete fCPUTimer;
284    delete fInitTimer;
285 }
286
287 //______________________________________________________________________________
288 void AliAnalysisManager::CreateReadCache()
289 {
290 // Create cache for reading according fCacheSize and fAsyncReading.
291    if (!fTree || !fTree->GetCurrentFile()) {
292       Error("CreateReadCache","Current tree or tree file not yet defined");
293       return;
294    }   
295    if (!fCacheSize) {
296       if (fDebug) Info("CreateReadCache","=== Read caching disabled ===");
297       return;
298    }
299    gEnv->SetValue("TFile.AsyncPrefetching",(Int_t)fAsyncReading);
300    if (fAsyncReading) gEnv->SetValue("Cache.Directory",Form("file://%s/cache", gSystem->WorkingDirectory()));
301    if (fAsyncReading) gEnv->SetValue("TFile.AsyncReading",1);
302    fTree->SetCacheSize(fCacheSize);
303    TTreeCache::SetLearnEntries(1);  //<<< we can take the decision after 1 entry
304    fTree->AddBranchToCache("*",kTRUE);    //<<< add all branches to the cache
305    if (fDebug) {
306       Info("CreateReadCache","Read cache enabled %lld bytes with async reading=%d",fCacheSize, (Int_t)fAsyncReading);
307    }
308    return;
309 }   
310
311 //______________________________________________________________________________
312 Bool_t AliAnalysisManager::EventLoop(Long64_t nevents)
313 {
314 // Initialize an event loop where the data producer is the input handler
315 // The handler must implement MakeTree creating the tree of events (likely
316 // memory resident) and generate the current event in the method BeginEvent.
317 // If the tree is memory resident, the handler should never call TTree::Fill
318 // method.
319    cout << "===== RUNNING IN EVENT LOOP MODE: " << GetName() << endl;
320    if (!fInputEventHandler) {
321      Error("EventLoop", "No input handler: exiting");
322      return kFALSE;
323    }
324    TTree *tree = new TTree("DummyTree", "Dummy tree for AliAnalysisManager::EventLoop");
325    SetExternalLoop(kTRUE);
326    if (!Init(tree)) return kFALSE;
327    SlaveBegin(tree);
328    for (Long64_t iev=0; iev<nevents; iev++)
329       ExecAnalysis();
330    TList dummyList;
331    PackOutput(&dummyList);
332    fIsRemote = kTRUE;
333    Terminate();
334    return kTRUE;
335 }
336       
337 //______________________________________________________________________________
338 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
339 {
340 // Read one entry of the tree or a whole branch.
341    fCurrentEntry = entry;
342    if (!fAutoBranchHandling)
343      return 123456789;
344    if (!fTree || !fTree->GetTree()) return -1;
345    fIOTimer->Start(kTRUE); 
346    Long64_t readbytes = fTree->GetTree()->GetEntry(entry, getall);
347    fIOTimer->Stop();
348    fIOTime += fIOTimer->RealTime();
349    return (Int_t)readbytes;
350 }
351
352 //______________________________________________________________________________
353 Int_t AliAnalysisManager::GetRunFromAlienPath(const char *path)
354 {
355 // Attempt to extract run number from input data path. Works only for paths to
356 // alice data in alien.
357 //    sim:  /alice/sim/<production>/run_no/...
358 //    data: /alice/data/year/period/000run_no/... (ESD or AOD)
359    TString type = "unknown";
360    TString s(path);
361    if (s.Contains("/alice/data")) type = "real";
362    else if (s.Contains("/alice/sim")) type = "simulated";
363    TString srun;
364    Int_t ind1, ind2;
365    ind1 = s.Index("/00");
366    if (ind1>0) {
367       ind2 = s.Index("/",ind1+1);
368       if (ind2-ind1>8) srun = s(ind1+1, ind2-ind1-1);
369    }   
370    if (srun.IsNull()) {
371       ind1 = s.Index("/LHC");
372       if (ind1>0) {
373          ind1 = s.Index("/",ind1+1);
374          if (ind1>0) {
375             ind2 = s.Index("/",ind1+1);
376             if (ind2>0) srun = s(ind1+1, ind2-ind1-1);
377          }
378       }
379    }         
380    Int_t run = srun.Atoi();
381    if (run>0) printf("=== GetRunFromAlienPath: run %d of %s data ===\n", run, type.Data());
382    return run;
383 }   
384
385 //______________________________________________________________________________
386 Bool_t AliAnalysisManager::Init(TTree *tree)
387 {
388   // The Init() function is called when the selector needs to initialize
389   // a new tree or chain. Typically here the branch addresses of the tree
390   // will be set. It is normaly not necessary to make changes to the
391   // generated code, but the routine can be extended by the user if needed.
392   // Init() will be called many times when running with PROOF.
393    Bool_t init = kFALSE;
394    if (!tree) return kFALSE; // Should not happen - protected in selector caller
395    if (fDebug > 1) {
396       printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
397    }
398    // Call InitTree of EventHandler
399    if (fOutputEventHandler) {
400       if (fMode == kProofAnalysis) {
401          init = fOutputEventHandler->Init(0x0, "proof");
402       } else {
403          init = fOutputEventHandler->Init(0x0, "local");
404       }
405       if (!init) {
406          Error("Init", "Output event handler failed to initialize");
407          return kFALSE;
408       }         
409    }
410    
411    if (fInputEventHandler) {
412       if (fMode == kProofAnalysis) {
413          init = fInputEventHandler->Init(tree, "proof");
414       } else {
415          init = fInputEventHandler->Init(tree, "local");
416       }
417       if (!init) {
418          Error("Init", "Input event handler failed to initialize tree"); 
419          return kFALSE;
420       }         
421    } else {
422       // If no input event handler we need to get the tree once
423       // for the chain
424       if(!tree->GetTree()) {
425          Long64_t readEntry = tree->LoadTree(0);
426          if (readEntry == -2) {
427             Error("Init", "Input tree has no entry. Exiting");
428             return kFALSE;
429          }
430       }   
431    }
432
433    if (fMCtruthEventHandler) {
434       if (fMode == kProofAnalysis) {
435          init = fMCtruthEventHandler->Init(0x0, "proof");
436       } else {
437          init = fMCtruthEventHandler->Init(0x0, "local");
438       }
439       if (!init) {
440          Error("Init", "MC event handler failed to initialize"); 
441          return kFALSE;
442       }         
443    }
444
445    if (!fInitOK) InitAnalysis();
446    if (!fInitOK) return kFALSE;
447    fTree = tree;
448    if (fMode != kProofAnalysis) CreateReadCache();
449    else {
450      // cholm - here we should re-add to the table or branches 
451      fTable.Clear();
452    }
453    AliAnalysisDataContainer *top = fCommonInput;
454    if (!top) top = (AliAnalysisDataContainer*)fInputs->At(0);
455    if (!top) {
456       Error("Init","No top input container !");
457       return kFALSE;
458    }
459    top->SetData(tree);
460    CheckBranches(kFALSE);
461    fTable.Rehash(100);
462    if (fDebug > 1) {
463       printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
464    }
465    return kTRUE;
466 }
467
468 //______________________________________________________________________________
469 void AliAnalysisManager::SlaveBegin(TTree *tree)
470 {
471   // The SlaveBegin() function is called after the Begin() function.
472   // When running with PROOF SlaveBegin() is called on each slave server.
473   // The tree argument is deprecated (on PROOF 0 is passed).
474    if (fDebug > 1) printf("->AliAnalysisManager::SlaveBegin()\n");
475    // Init timer should be already started
476    // Apply debug options
477    ApplyDebugOptions();
478    if (fCacheSize && 
479        fMCtruthEventHandler &&
480        (fMode != kProofAnalysis)) fMCtruthEventHandler->SetCacheSize(fCacheSize);
481    if (!CheckTasks()) Fatal("SlaveBegin", "Not all needed libraries were loaded");
482    static Bool_t isCalled = kFALSE;
483    Bool_t init = kFALSE;
484    Bool_t initOK = kTRUE;
485    TString msg;
486    TDirectory *curdir = gDirectory;
487    // Call SlaveBegin only once in case of mixing
488    if (isCalled && fMode==kMixingAnalysis) return;
489    gROOT->cd();
490    // Call Init of EventHandler
491    if (fOutputEventHandler) {
492       if (fMode == kProofAnalysis) {
493          // Merging AOD's in PROOF via TProofOutputFile
494          if (fDebug > 1) printf("   Initializing AOD output file %s...\n", fOutputEventHandler->GetOutputFileName());
495          init = fOutputEventHandler->Init("proof");
496          if (!init) msg = "Failed to initialize output handler on worker";
497       } else {
498          init = fOutputEventHandler->Init("local");
499          if (!init) msg = "Failed to initialize output handler";
500       }
501       initOK &= init;
502       if (!fSelector) Error("SlaveBegin", "Selector not set");
503       else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
504    }
505    gROOT->cd();
506    if (fInputEventHandler) {
507       fInputEventHandler->SetInputTree(tree);
508       if (fMode == kProofAnalysis) {
509          init = fInputEventHandler->Init("proof");
510          if (!init) msg = "Failed to initialize input handler on worker";
511       } else {
512          init = fInputEventHandler->Init("local");
513          if (!init) msg = "Failed to initialize input handler";
514       }
515       initOK &= init;
516       if (!fSelector) Error("SlaveBegin", "Selector not set");      
517       else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
518    }
519    gROOT->cd();
520    if (fMCtruthEventHandler) {
521       if (fMode == kProofAnalysis) {
522          init = fMCtruthEventHandler->Init("proof");
523          if (!init) msg = "Failed to initialize MC handler on worker";
524       } else {
525          init = fMCtruthEventHandler->Init("local");
526          if (!init) msg = "Failed to initialize MC handler";
527       }
528       initOK &= init;
529       if (!fSelector) Error("SlaveBegin", "Selector not set");      
530       else if (!init) {fSelector->Abort(msg); fSelector->SetStatus(-1);}
531    }
532    if (curdir) curdir->cd();
533    isCalled = kTRUE;
534    if (!initOK) return;   
535    TIter next(fTasks);
536    AliAnalysisTask *task;
537    // Call CreateOutputObjects for all tasks
538    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
539    Bool_t dirStatus = TH1::AddDirectoryStatus();
540    Int_t itask = 0;
541    while ((task=(AliAnalysisTask*)next())) {
542       gROOT->cd();
543       // Start with memory as current dir and make sure by default histograms do not get attached to files.
544       TH1::AddDirectory(kFALSE);
545       task->CreateOutputObjects();
546       if (!task->CheckPostData()) {
547          Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
548                 Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
549                 ####### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ##########", task->GetName(), task->ClassName());
550       }
551       if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
552       itask++;
553    }
554    TH1::AddDirectory(dirStatus);
555    if (curdir) curdir->cd();
556    fInitTimer->Stop();
557    fInitTime += fInitTimer->RealTime();
558    fInitTimer->Continue();
559    printf("Initialization time: %g [sec]\n", fInitTime);
560    if (fDebug > 1) printf("<-AliAnalysisManager::SlaveBegin()\n");
561 }
562
563 //______________________________________________________________________________
564 Bool_t AliAnalysisManager::Notify()
565 {
566    // The Notify() function is called when a new file is opened. This
567    // can be either for a new TTree in a TChain or when when a new TTree
568    // is started when using PROOF. It is normaly not necessary to make changes
569    // to the generated code, but the routine can be extended by the
570    // user if needed. The return value is currently not used.
571    fIOTimer->Start(kTRUE); 
572    if (!fTree) return kFALSE;
573    if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) return kFALSE;
574
575    fTable.Clear("nodelete"); // clearing the hash table may not be needed -> C.L.
576    if (fMode == kProofAnalysis) fIsRemote = kTRUE;
577
578    TFile *curfile = fTree->GetCurrentFile();
579    if (!curfile) {
580       Error("Notify","No current file");
581       return kFALSE;
582    }  
583    if (IsCollectThroughput()) {
584       if (fCurrentDescriptor) fCurrentDescriptor->Done();
585       fCurrentDescriptor = new AliAnalysisFileDescriptor(curfile);
586       fFileDescriptors->Add(fCurrentDescriptor);
587    } 
588    
589    if (fDebug > 1) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
590    Int_t run = AliAnalysisManager::GetRunFromAlienPath(curfile->GetName());
591    if (run && (run != fRunFromPath)) {
592       fRunFromPath = run;
593       if (fDebug > 1) printf("   ### run found from path: %d\n", run);
594    }
595    TIter next(fTasks);
596    AliAnalysisTask *task;
597         
598    // Call Notify of the event handlers
599    if (fInputEventHandler) {
600        fInputEventHandler->Notify(curfile->GetName());
601    }
602
603    if (fOutputEventHandler) {
604        fOutputEventHandler->Notify(curfile->GetName());
605    }
606
607    if (fMCtruthEventHandler) {
608        fMCtruthEventHandler->Notify(curfile->GetName());
609    }
610
611    // Call Notify for all tasks
612    while ((task=(AliAnalysisTask*)next())) 
613       task->Notify();
614
615    if (fDebug > 1) printf("<-AliAnalysisManager::Notify()\n");
616    fIOTimer->Stop();
617    fIOTime += fIOTimer->RealTime();
618    return kTRUE;
619 }    
620
621 //______________________________________________________________________________
622 Bool_t AliAnalysisManager::Process(Long64_t)
623 {
624   // The Process() function is called for each entry in the tree (or possibly
625   // keyed object in the case of PROOF) to be processed. The entry argument
626   // specifies which entry in the currently loaded tree is to be processed.
627   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
628   // to read either all or the required parts of the data. When processing
629   // keyed objects with PROOF, the object is already loaded and is available
630   // via the fObject pointer.
631   //
632   // This function should contain the "body" of the analysis. It can contain
633   // simple or elaborate selection criteria, run algorithms on the data
634   // of the event and typically fill histograms.
635
636   // WARNING when a selector is used with a TChain, you must use
637   //  the pointer to the current TTree to call GetEntry(entry).
638   //  The entry is always the local entry number in the current tree.
639   //  Assuming that fChain is the pointer to the TChain being processed,
640   //  use fChain->GetTree()->GetEntry(entry).
641
642    // This method is obsolete. ExecAnalysis is called instead.
643    return kTRUE;
644 }
645
646 //______________________________________________________________________________
647 void AliAnalysisManager::PackOutput(TList *target)
648 {
649   // Pack all output data containers in the output list. Called at SlaveTerminate
650   // stage in PROOF case for each slave.
651    if (fDebug > 1) printf("->AliAnalysisManager::PackOutput()\n");
652    fIOTimer->Start(kTRUE);
653    std::ofstream out;
654    if (IsCollectThroughput()) {
655       if (fCurrentDescriptor) fCurrentDescriptor->Done();
656       fFileDescriptors->Print();
657       if (fFileInfoLog.IsNull()) fFileInfoLog = "fileinfo.log";
658       out.open(fFileInfoLog, std::ios::app);
659       if (out.bad()) Error("SavePrimitive", "Bad file name: %s", fFileInfoLog.Data());
660       else {
661          TIter nextflog(fFileDescriptors);
662          TObject *log;
663          while ((log=nextflog())) log->SavePrimitive(out,"");
664       }
665    }   
666    if (!target) {
667       Error("PackOutput", "No target. Exiting.");
668       return;
669    }
670    TDirectory *cdir = gDirectory;
671    gROOT->cd();
672    if (fInputEventHandler)   fInputEventHandler  ->Terminate();
673    if (fOutputEventHandler)  fOutputEventHandler ->Terminate();
674    if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
675    gROOT->cd();
676
677    // Call FinishTaskOutput() for each event loop task (not called for 
678    // post-event loop tasks - use Terminate() fo those)
679    TIter nexttask(fTasks);
680    AliAnalysisTask *task;
681    while ((task=(AliAnalysisTask*)nexttask())) {
682       if (!task->IsPostEventLoop()) {
683          if (fDebug > 1) printf("->FinishTaskOutput: task %s\n", task->GetName());
684          task->FinishTaskOutput();
685          gROOT->cd();
686          if (fDebug > 1) printf("<-FinishTaskOutput: task %s\n", task->GetName());
687       }
688    }
689    // Write statistics message on the workers.
690    if (fStatistics) WriteStatisticsMsg(fNcalls);
691    
692    if (fMode == kProofAnalysis) {
693       TIter next(fOutputs);
694       AliAnalysisDataContainer *output;
695       Bool_t isManagedByHandler = kFALSE;
696       TList filestmp;
697       filestmp.SetOwner();
698       while ((output=(AliAnalysisDataContainer*)next())) {
699          // Do not consider outputs of post event loop tasks
700          isManagedByHandler = kFALSE;
701          if (output->GetProducer() && output->GetProducer()->IsPostEventLoop()) continue;
702          const char *filename = output->GetFileName();
703          if (!(strcmp(filename, "default")) && fOutputEventHandler) {
704             isManagedByHandler = kTRUE;
705             printf("#### Handler output. Extra: %s\n", fExtraFiles.Data());
706             filename = fOutputEventHandler->GetOutputFileName();
707          }
708          // Check if data was posted to this container. If not, issue an error.
709          if (!output->GetData() && !isManagedByHandler) {
710             Error("PackOutput", "No data for output container %s. Forgot to PostData ?", output->GetName());
711             continue;
712          }   
713          if (!output->IsSpecialOutput()) {
714             // Normal outputs
715             if (strlen(filename) && !isManagedByHandler) {
716                // Backup current folder
717                TDirectory *opwd = gDirectory;
718                // File resident outputs. 
719                // Check first if the file exists.
720                TString openoption = "RECREATE";
721                Bool_t firsttime = kTRUE;
722                if (filestmp.FindObject(output->GetFileName())) {
723                   firsttime = kFALSE;
724                } else {   
725                   filestmp.Add(new TNamed(output->GetFileName(),""));
726                }   
727                if (!gSystem->AccessPathName(output->GetFileName()) && !firsttime) openoption = "UPDATE";
728 //               TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
729                // Save data to file, then close.
730                if (output->GetData()->InheritsFrom(TCollection::Class())) {
731                   // If data is a collection, we set the name of the collection 
732                   // as the one of the container and we save as a single key.
733                   TCollection *coll = (TCollection*)output->GetData();
734                   coll->SetName(output->GetName());
735 //                  coll->Write(output->GetName(), TObject::kSingleKey);
736                } else {
737                   if (output->GetData()->InheritsFrom(TTree::Class())) {
738                      TFile *file = AliAnalysisManager::OpenFile(output, openoption, kTRUE);
739                      // Save data to file, then close.
740                      TTree *tree = (TTree*)output->GetData();
741                      // Check if tree is in memory
742                      if (tree->GetDirectory()==gROOT) tree->SetDirectory(gDirectory);
743                      tree->AutoSave();
744                      file->Close();
745                   } else {
746 //                     output->GetData()->Write();
747                   }   
748                }      
749                if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
750 //               if (fDebug > 2) {
751 //                  printf("   file %s listing content:\n", filename);
752 //                  file->ls();
753 //               }   
754                // Clear file list to release object ownership to user.
755 //               file->Clear();
756 //               file->Close();
757                output->SetFile(NULL);
758                // Restore current directory
759                if (opwd) opwd->cd();
760             } else {
761                // Memory-resident outputs   
762                if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
763             }   
764             AliAnalysisDataWrapper *wrap = 0;
765             if (isManagedByHandler) {
766                wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
767                wrap->SetName(output->GetName());
768             }   
769             else                    wrap =output->ExportData();
770             // Output wrappers must NOT delete data after merging - the user owns them
771             wrap->SetDeleteData(kFALSE);
772             target->Add(wrap);
773          } else {
774          // Special outputs. The file must be opened and connected to the container.
775             TDirectory *opwd = gDirectory;
776             TFile *file = output->GetFile();
777             if (!file) {
778                AliAnalysisTask *producer = output->GetProducer();
779                Fatal("PackOutput", 
780                      "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
781                      output->GetFileName(), output->GetName(), producer->ClassName());
782                continue;
783             }   
784             TString outFilename = file->GetName();
785             if (fDebug > 1) printf("PackOutput %s: special output\n", output->GetName());
786             if (isManagedByHandler) {
787                // Terminate IO for files managed by the output handler
788                // file->Write() moved to AOD handler (A.G. 11.01.10)
789 //               if (file) file->Write();
790                if (file && fDebug > 2) {
791                   printf("   handled file %s listing content:\n", file->GetName());
792                   file->ls();
793                }   
794                fOutputEventHandler->TerminateIO();
795             } else {               
796                file->cd();
797                // Release object ownership to users after writing data to file
798                if (output->GetData()->InheritsFrom(TCollection::Class())) {
799                   // If data is a collection, we set the name of the collection 
800                   // as the one of the container and we save as a single key.
801                   TCollection *coll = (TCollection*)output->GetData();
802                   coll->SetName(output->GetName());
803                   coll->Write(output->GetName(), TObject::kSingleKey);
804                } else {
805                   if (output->GetData()->InheritsFrom(TTree::Class())) {
806                      TTree *tree = (TTree*)output->GetData();
807                      tree->SetDirectory(file);
808                      tree->AutoSave();
809                   } else {
810                      output->GetData()->Write();
811                   }   
812                }      
813                if (fDebug > 2) {
814                   printf("   file %s listing content:\n", output->GetFileName());
815                   file->ls();
816                }
817                // Clear file list to release object ownership to user.
818 //               file->Clear();
819                file->Close();
820                output->SetFile(NULL);
821             }
822             // Restore current directory
823             if (opwd) opwd->cd();
824             // Check if a special output location was provided or the output files have to be merged
825             if (strlen(fSpecialOutputLocation.Data())) {
826                TString remote = fSpecialOutputLocation;
827                remote += "/";
828                Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
829                if (remote.BeginsWith("alien:")) {
830                   gROOT->ProcessLine("TGrid::Connect(\"alien:\", gProofServ->GetUser());");
831                   remote += outFilename;
832                   remote.ReplaceAll(".root", Form("_%d.root", gid));
833                } else {   
834                   remote += Form("%s_%d_", gSystem->HostName(), gid);
835                   remote += outFilename;
836                }   
837                if (fDebug > 1) 
838                   Info("PackOutput", "Output file for container %s to be copied \n   at: %s. No merging.",
839                        output->GetName(), remote.Data());
840                TFile::Cp ( outFilename.Data(), remote.Data() );
841                // Copy extra outputs
842                if (fExtraFiles.Length() && isManagedByHandler) {
843                   TObjArray *arr = fExtraFiles.Tokenize(" ");
844                   TObjString *os;
845                   TIter nextfilename(arr);
846                   while ((os=(TObjString*)nextfilename())) {
847                      outFilename = os->GetString();
848                      remote = fSpecialOutputLocation;
849                      remote += "/";
850                      if (remote.BeginsWith("alien://")) {
851                         remote += outFilename;
852                         remote.ReplaceAll(".root", Form("_%d.root", gid));
853                      } else {   
854                         remote += Form("%s_%d_", gSystem->HostName(), gid);
855                         remote += outFilename;
856                      }   
857                      if (fDebug > 1) 
858                         Info("PackOutput", "Extra AOD file %s to be copied \n   at: %s. No merging.",
859                              outFilename.Data(), remote.Data());
860                      TFile::Cp ( outFilename.Data(), remote.Data() );
861                   }   
862                   delete arr;
863                }   
864             } else {
865             // No special location specified-> use TProofOutputFile as merging utility
866             // The file at this output slot must be opened in CreateOutputObjects
867                if (fDebug > 1) printf("   File for container %s to be merged via file merger...\n", output->GetName());
868             }
869          }      
870       }
871    } 
872    fIOTime += fIOTimer->RealTime();
873    if ((fDebug || IsCollectThroughput())) {
874       fInitTimer->Stop();
875       fInitTime = fInitTimer->RealTime()-fIOTime-fCPUTime;
876       printf("=Analysis %s= init time:       %g[sec]\
877             \n              I/O & data mng.: %g [sec]\
878             \n              task execution: %g [sec]\
879             \n              total time:     CPU=%g [sec]  REAL=%g[sec]\n",
880             GetName(), fInitTime, fIOTime, fCPUTime, fInitTimer->CpuTime(), fInitTimer->RealTime());
881       if (IsCollectThroughput()) {
882          out << "#summary#########################################################" << endl;
883          out << "train_name   " << GetName() << endl;
884          out << "root_time    " << fInitTimer->RealTime() << endl;
885          out << "root_cpu     " << fInitTimer->CpuTime() << endl;
886          out << "init_time    " << fInitTime << endl;
887          out << "io_mng_time  " << fIOTime << endl;
888          out << "exec_time    " << fCPUTime << endl;
889          TString aliensite = gSystem->Getenv("ALIEN_SITE");
890          out << "alien_site   " << aliensite << endl;
891          out << "host_name    ";
892          TString hostname = gSystem->Getenv("ALIEN_HOSTNAME");
893          if (hostname.IsNull()) {
894             out.close();
895             gSystem->Exec(Form("hostname -f >> %s", fFileInfoLog.Data()));
896          } else {
897             out << hostname << endl;
898          }   
899       }
900    }
901               
902    if (cdir) cdir->cd();
903    if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
904 }
905
906 //______________________________________________________________________________
907 void AliAnalysisManager::ImportWrappers(TList *source)
908 {
909 // Import data in output containers from wrappers coming in source.
910    if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
911    fIOTimer->Start(kTRUE);
912    TIter next(fOutputs);
913    AliAnalysisDataContainer *cont;
914    AliAnalysisDataWrapper   *wrap;
915    Int_t icont = 0;
916    Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
917    TDirectory *cdir = gDirectory;
918    while ((cont=(AliAnalysisDataContainer*)next())) {
919       wrap = 0;
920       if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
921       if (cont->IsRegisterDataset()) continue;
922       const char *filename = cont->GetFileName();
923       Bool_t isManagedByHandler = kFALSE;
924       if (!(strcmp(filename, "default")) && fOutputEventHandler) {
925          isManagedByHandler = kTRUE;
926          filename = fOutputEventHandler->GetOutputFileName();
927       }
928       if (cont->IsSpecialOutput() || inGrid) {
929          if (strlen(fSpecialOutputLocation.Data())) continue;
930          // Copy merged file from PROOF scratch space. 
931          // In case of grid the files are already in the current directory.
932          if (!inGrid) {
933             if (isManagedByHandler && fExtraFiles.Length()) {
934                // Copy extra registered dAOD files.
935                TObjArray *arr = fExtraFiles.Tokenize(" ");
936                TObjString *os;
937                TIter nextfilename(arr);
938                while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
939                delete arr;
940             }
941             if (!GetFileFromWrapper(filename, source)) continue;
942          }   
943          // Normally we should connect data from the copied file to the
944          // corresponding output container, but it is not obvious how to do this
945          // automatically if several objects in file...
946          TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
947          if (!f) f = TFile::Open(filename, "READ");
948          if (!f) {
949             Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
950             continue;
951          }   
952          f->cd();
953          TObject *obj = 0;
954          // Cd to the directory pointed by the container
955          TString folder = cont->GetFolderName();
956          if (!folder.IsNull()) f->cd(folder);
957          // Try to fetch first an object having the container name.
958          obj = gDirectory->Get(cont->GetName());
959          if (!obj) {
960             Warning("ImportWrappers", "Could not import object of type:%s for container %s in file %s:%s.\n Object will not be available in Terminate(). Try if possible to name the output object as the container (%s) or to embed it in a TList", 
961                     cont->GetType()->GetName(), cont->GetName(), filename, cont->GetFolderName(), cont->GetName());
962             continue;
963          }  
964          wrap = new AliAnalysisDataWrapper(obj);
965          wrap->SetDeleteData(kFALSE);
966       }   
967       if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
968       if (!wrap) {
969          Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
970          continue;
971       }
972       icont++;
973       if (fDebug > 1) {
974          printf("   Importing data for container %s\n", cont->GetName());
975          if (strlen(filename)) printf("    -> file %s\n", filename);
976          else printf("\n");
977       }   
978       cont->ImportData(wrap);
979    }
980    if (cdir) cdir->cd();
981    fIOTimer->Stop();
982    fIOTime += fIOTimer->RealTime();
983    if (fDebug > 1) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
984 }
985
986 //______________________________________________________________________________
987 void AliAnalysisManager::UnpackOutput(TList *source)
988 {
989   // Called by AliAnalysisSelector::Terminate only on the client.
990    fIOTimer->Start(kTRUE);
991    if (fDebug > 1) printf("->AliAnalysisManager::UnpackOutput()\n");
992    if (!source) {
993       Error("UnpackOutput", "No target. Exiting.");
994       return;
995    }
996    if (fDebug > 1) printf("   Source list contains %d containers\n", source->GetSize());
997
998    if (fMode == kProofAnalysis) ImportWrappers(source);
999
1000    TIter next(fOutputs);
1001    AliAnalysisDataContainer *output;
1002    while ((output=(AliAnalysisDataContainer*)next())) {
1003       if (!output->GetData()) continue;
1004       // Check if there are client tasks that run post event loop
1005       if (output->HasConsumers()) {
1006          // Disable event loop semaphore
1007          output->SetPostEventLoop(kTRUE);
1008          TObjArray *list = output->GetConsumers();
1009          Int_t ncons = list->GetEntriesFast();
1010          for (Int_t i=0; i<ncons; i++) {
1011             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
1012             task->CheckNotify(kTRUE);
1013             // If task is active, execute it
1014             if (task->IsPostEventLoop() && task->IsActive()) {
1015                if (fDebug > 1) printf("== Executing post event loop task %s\n", task->GetName());
1016                if (fStatistics) fStatistics->StartTimer(GetTaskIndex(task), task->GetName(), task->ClassName());
1017                task->ExecuteTask();
1018             }   
1019          }
1020          if (fStatistics) fStatistics->StopTimer();
1021       }   
1022    }
1023    fIOTimer->Stop();
1024    fIOTime += fIOTimer->RealTime();
1025    if (fDebug > 1) printf("<-AliAnalysisManager::UnpackOutput()\n");
1026 }
1027
1028 //______________________________________________________________________________
1029 void AliAnalysisManager::Terminate()
1030 {
1031   // The Terminate() function is the last function to be called during
1032   // a query. It always runs on the client, it can be used to present
1033   // the results graphically.
1034    if (fDebug > 1) printf("->AliAnalysisManager::Terminate()\n");
1035    fInitTimer->Start(kTRUE);
1036    TDirectory *cdir = gDirectory;
1037    gROOT->cd();
1038    AliAnalysisTask *task;
1039    AliAnalysisDataContainer *output;
1040    TIter next(fTasks);
1041    TStopwatch timer;
1042    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1043    // Call Terminate() for tasks
1044    Int_t itask = 0;
1045    while (!IsSkipTerminate() && (task=(AliAnalysisTask*)next())) {
1046       // Save all the canvases produced by the Terminate
1047       TString pictname = Form("%s_%s", task->GetName(), task->ClassName());
1048       task->Terminate();
1049       gROOT->cd();
1050       if (getsysInfo) 
1051          AliSysInfo::AddStamp(Form("%s_TERMINATE",task->ClassName()),0, itask, 2);
1052       itask++;   
1053       if (TObject::TestBit(kSaveCanvases)) {
1054          if (!gROOT->IsBatch()) {
1055             if (fDebug>1) printf("Waiting 5 sec for %s::Terminate() to finish drawing ...\n", task->ClassName());
1056             timer.Start(kTRUE);
1057             while (timer.RealTime()<5) {
1058                timer.Continue();
1059                gSystem->ProcessEvents();
1060             }
1061          }
1062          Int_t iend = gROOT->GetListOfCanvases()->GetEntries();
1063          if (iend==0) continue;
1064          TCanvas *canvas;
1065          for (Int_t ipict=0; ipict<iend; ipict++) {
1066             canvas = (TCanvas*)gROOT->GetListOfCanvases()->At(ipict);
1067             if (!canvas) continue;         
1068             canvas->SaveAs(Form("%s_%02d.gif", pictname.Data(),ipict));
1069          } 
1070          gROOT->GetListOfCanvases()->Delete(); 
1071       }
1072    }   
1073    //
1074    if (fInputEventHandler)   fInputEventHandler  ->TerminateIO();
1075    if (fOutputEventHandler)  fOutputEventHandler ->TerminateIO();
1076    if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
1077    gROOT->cd();
1078    TObjArray *allOutputs = new TObjArray();
1079    Int_t icont;
1080    for (icont=0; icont<fOutputs->GetEntriesFast(); icont++) allOutputs->Add(fOutputs->At(icont));
1081    if (!IsSkipTerminate())
1082       for (icont=0; icont<fParamCont->GetEntriesFast(); icont++) allOutputs->Add(fParamCont->At(icont));
1083    TIter next1(allOutputs);
1084    TString handlerFile = "";
1085    TString extraOutputs = "";
1086    if (fOutputEventHandler) {
1087       handlerFile = fOutputEventHandler->GetOutputFileName();
1088       extraOutputs = fOutputEventHandler->GetExtraOutputs();
1089    }
1090    icont = 0;
1091    TList filestmp;
1092    while ((output=(AliAnalysisDataContainer*)next1())) {
1093       // Special outputs or grid files have the files already closed and written.
1094       icont++;
1095       if (fMode == kGridAnalysis && icont<=fOutputs->GetEntriesFast()) continue;
1096       if (fMode == kProofAnalysis) {
1097         if (output->IsSpecialOutput() || output->IsRegisterDataset()) continue;
1098       }  
1099       const char *filename = output->GetFileName();
1100       TString openoption = "RECREATE";
1101       if (!(strcmp(filename, "default"))) continue;
1102       if (!strlen(filename)) continue;
1103       if (!output->GetData()) continue;
1104       TDirectory *opwd = gDirectory;
1105       TFile *file = output->GetFile();
1106       if (!file) file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1107       if (!file) {
1108               //if (handlerFile == filename && !gSystem->AccessPathName(filename)) openoption = "UPDATE";
1109          Bool_t firsttime = kTRUE;
1110          if (filestmp.FindObject(filename) || extraOutputs.Contains(filename)) {
1111             firsttime = kFALSE;
1112          } else {   
1113             filestmp.Add(new TNamed(filename,""));
1114          }   
1115          if (!gSystem->AccessPathName(filename) && !firsttime) openoption = "UPDATE";
1116               if (fDebug>1) printf("Opening file: %s  option=%s\n",filename, openoption.Data());
1117          file = new TFile(filename, openoption);
1118       } else {
1119          if (fDebug>1) printf("File <%s> already opened with option: <%s> \n", filename, file->GetOption());
1120          openoption = file->GetOption();
1121          if (openoption == "READ") {
1122             if (fDebug>1) printf("...reopening in UPDATE mode\n");
1123             file->ReOpen("UPDATE");            
1124          }
1125       }   
1126       if (file->IsZombie()) {
1127          Error("Terminate", "Cannot open output file %s", filename);
1128          continue;
1129       }   
1130       output->SetFile(file);
1131       file->cd();
1132       // Check for a folder request
1133       TString dir = output->GetFolderName();
1134       if (!dir.IsNull()) {
1135          if (!file->GetDirectory(dir)) file->mkdir(dir);
1136          file->cd(dir);
1137       }  
1138       if (fDebug > 1) printf("...writing container %s to file %s:%s\n", output->GetName(), file->GetName(), output->GetFolderName());
1139       if (output->GetData()->InheritsFrom(TCollection::Class())) {
1140       // If data is a collection, we set the name of the collection 
1141       // as the one of the container and we save as a single key.
1142          TCollection *coll = (TCollection*)output->GetData();
1143          coll->SetName(output->GetName());
1144          coll->Write(output->GetName(), TObject::kSingleKey);
1145       } else {
1146          if (output->GetData()->InheritsFrom(TTree::Class())) {
1147             TTree *tree = (TTree*)output->GetData();
1148             tree->SetDirectory(gDirectory);
1149             tree->AutoSave();
1150          } else {
1151             output->GetData()->Write();
1152          }   
1153       }      
1154       if (opwd) opwd->cd();
1155    }
1156    gROOT->cd();
1157    next1.Reset();
1158    TString copiedFiles;
1159    while ((output=(AliAnalysisDataContainer*)next1())) {
1160       // Close all files at output
1161       TDirectory *opwd = gDirectory;
1162       if (output->GetFile()) {
1163          // Clear file list to release object ownership to user.
1164 //         output->GetFile()->Clear();
1165          output->GetFile()->Close();
1166          // Copy merged outputs in alien if requested
1167          if (fSpecialOutputLocation.BeginsWith("alien://")) {
1168             if (copiedFiles.Contains(output->GetFile()->GetName())) {
1169                if (opwd) opwd->cd();
1170                output->SetFile(NULL);
1171                continue;
1172             } 
1173             Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data()); 
1174             gROOT->ProcessLine("if (!gGrid) TGrid::Connect(\"alien:\");");
1175             TFile::Cp(output->GetFile()->GetName(), 
1176                       Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
1177             copiedFiles += output->GetFile()->GetName();
1178          }             
1179          output->SetFile(NULL);
1180       }   
1181       if (opwd) opwd->cd();
1182    }   
1183    delete allOutputs;
1184    //Write statistics information on the client
1185    if (fStatistics) WriteStatisticsMsg(fNcalls);
1186    if (getsysInfo) {
1187       TDirectory *crtdir = gDirectory;
1188       TFile f("syswatch.root", "RECREATE");
1189       TH1 *hist;
1190       TString cut;
1191       if (!f.IsZombie()) {
1192          TTree *tree = AliSysInfo::MakeTree("syswatch.log");
1193          tree->SetName("syswatch");
1194          tree->SetMarkerStyle(kCircle);
1195          tree->SetMarkerColor(kBlue);
1196          tree->SetMarkerSize(0.5);
1197          if (!gROOT->IsBatch()) {
1198             tree->SetAlias("event", "id0");
1199             tree->SetAlias("task",  "id1");
1200             tree->SetAlias("stage", "id2");
1201             // Already defined aliases
1202             // tree->SetAlias("deltaT","stampSec-stampOldSec");
1203             // tree->SetAlias("T","stampSec-first");
1204             // tree->SetAlias("deltaVM","(pI.fMemVirtual-pIOld.fMemVirtual)");
1205             // tree->SetAlias("VM","pI.fMemVirtual");
1206             TCanvas *canvas = new TCanvas("SysInfo","SysInfo",10,10,1200,800);
1207             Int_t npads = 1 /*COO plot for all tasks*/ +
1208                           fTopTasks->GetEntries() /*Exec plot per task*/ +
1209                           1 /*Terminate plot for all tasks*/ +
1210                           1; /*vm plot*/
1211                           
1212             Int_t iopt = (Int_t)TMath::Sqrt((Double_t)npads);
1213             if (npads<iopt*(iopt+1))
1214                canvas->Divide(iopt, iopt+1, 0.01, 0.01);
1215             else
1216                canvas->Divide(iopt+1, iopt+1, 0.01, 0.01);
1217             Int_t ipad = 1;
1218             // draw the plot of deltaVM for Exec for each task
1219             for (itask=0; itask<fTopTasks->GetEntriesFast(); itask++) {
1220                task = (AliAnalysisTask*)fTopTasks->At(itask);
1221                canvas->cd(ipad++);
1222                cut = Form("task==%d && stage==1", itask);
1223                tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1224                hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");            
1225                if (hist) {
1226                   hist->SetTitle(Form("%s: Exec dVM[MB]/event", task->GetName()));
1227                   hist->GetYaxis()->SetTitle("deltaVM [MB]");
1228                }   
1229             }
1230             // Draw the plot of deltaVM for CreateOutputObjects for all tasks
1231             canvas->cd(ipad++);
1232             tree->SetMarkerStyle(kFullTriangleUp);
1233             tree->SetMarkerColor(kRed);
1234             tree->SetMarkerSize(0.8);
1235             cut = "task>=0 && task<1000 && stage==0";
1236             tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1237             hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");            
1238             if (hist) {
1239                hist->SetTitle("Memory in CreateOutputObjects()");
1240                hist->GetYaxis()->SetTitle("deltaVM [MB]");
1241                hist->GetXaxis()->SetTitle("task");
1242             }   
1243             // draw the plot of deltaVM for Terminate for all tasks
1244             canvas->cd(ipad++);
1245             tree->SetMarkerStyle(kOpenSquare);
1246             tree->SetMarkerColor(kMagenta);
1247             cut = "task>=0 && task<1000 && stage==2";
1248             tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1249             hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1250             if (hist) {
1251                hist->SetTitle("Memory in Terminate()");
1252                hist->GetYaxis()->SetTitle("deltaVM [MB]");
1253                hist->GetXaxis()->SetTitle("task");
1254             }   
1255             // Full VM profile
1256             canvas->cd(ipad++);
1257             tree->SetMarkerStyle(kFullCircle);
1258             tree->SetMarkerColor(kGreen);
1259             cut = Form("task==%d && stage==1",fTopTasks->GetEntriesFast()-1);            
1260             tree->Draw("VM:event",cut,"", 1234567890, 0);
1261             hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1262             if (hist) {
1263                hist->SetTitle("Virtual memory");
1264                hist->GetYaxis()->SetTitle("VM [MB]");
1265             }
1266             canvas->Modified();   
1267          }   
1268          tree->SetMarkerStyle(kCircle);
1269          tree->SetMarkerColor(kBlue);
1270          tree->SetMarkerSize(0.5);
1271          tree->Write();
1272          f.Close();
1273          delete tree;
1274       }
1275       if (crtdir) crtdir->cd();
1276    }
1277    // Validate the output files
1278    if (ValidateOutputFiles() && fIsRemote && fMode!=kProofAnalysis) {
1279       ofstream out;
1280       out.open("outputs_valid", ios::out);
1281       out.close();
1282    }
1283    if (cdir) cdir->cd();  
1284    fInitTimer->Stop();
1285    if (fDebug || IsCollectThroughput()) {
1286       printf("=Analysis %s= Terminate time:  %g[sec]\n", GetName(), fInitTimer->RealTime());
1287    }
1288    if (fDebug > 1) printf("<-AliAnalysisManager::Terminate()\n");
1289 }
1290 //______________________________________________________________________________
1291 void AliAnalysisManager::ProfileTask(Int_t itop, const char *option) const
1292 {
1293 // Profiles the task having the itop index in the list of top (first level) tasks.
1294    AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->At(itop);
1295    if (!task) {
1296       Error("ProfileTask", "There are only %d top tasks in the manager", fTopTasks->GetEntries());
1297       return;
1298    }
1299    ProfileTask(task->GetName(), option);
1300 }      
1301
1302 //______________________________________________________________________________
1303 void AliAnalysisManager::ProfileTask(const char *name, const char */*option*/) const
1304 {
1305 // Profile a managed task after the execution of the analysis in case NSysInfo
1306 // was used.
1307    if (gSystem->AccessPathName("syswatch.root")) {
1308       Error("ProfileTask", "No file syswatch.root found in the current directory");
1309       return;
1310    }
1311    if (gROOT->IsBatch()) return;
1312    AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->FindObject(name);
1313    if (!task) {
1314       Error("ProfileTask", "No top task named %s known by the manager.", name);
1315       return;
1316    }
1317    Int_t itop = fTopTasks->IndexOf(task);
1318    Int_t itask = fTasks->IndexOf(task);
1319    // Create canvas with 2 pads: first draw COO + Terminate, second Exec
1320    TDirectory *cdir = gDirectory;
1321    TFile f("syswatch.root");
1322    TTree *tree = (TTree*)f.Get("syswatch");
1323    if (!tree) {
1324       Error("ProfileTask", "No tree named <syswatch> found in file syswatch.root");
1325       return;
1326    }   
1327    if (fDebug > 1) printf("=== Profiling task %s (class %s)\n", name, task->ClassName());
1328    TCanvas *canvas = new TCanvas(Form("profile_%d",itop),Form("Profile of task %s (class %s)",name,task->ClassName()),10,10,800,600);
1329    canvas->Divide(2, 2, 0.01, 0.01);
1330    Int_t ipad = 1;
1331    TString cut;
1332    TH1 *hist;
1333    // VM profile for COO and Terminate methods
1334    canvas->cd(ipad++);
1335    cut = Form("task==%d && (stage==0 || stage==2)",itask);
1336    tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1337    hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1338    if (hist) {
1339       hist->SetTitle("Alocated VM[MB] for COO and Terminate");
1340       hist->GetYaxis()->SetTitle("deltaVM [MB]");
1341       hist->GetXaxis()->SetTitle("method");
1342    }   
1343    // CPU profile per event
1344    canvas->cd(ipad++);
1345    cut = Form("task==%d && stage==1",itop);
1346    tree->Draw("deltaT:event",cut,"", 1234567890, 0);
1347    hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1348    if (hist) {
1349       hist->SetTitle("Execution time per event");
1350       hist->GetYaxis()->SetTitle("CPU/event [s]");
1351    }   
1352    // VM profile for Exec
1353    canvas->cd(ipad++);
1354    cut = Form("task==%d && stage==1",itop);
1355    tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1356    hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1357    if (hist) {
1358       hist->SetTitle("Alocated VM[MB] per event");
1359       hist->GetYaxis()->SetTitle("deltaVM [MB]");
1360    }   
1361    canvas->Modified();
1362    delete tree;
1363    f.Close();
1364    if (cdir) cdir->cd();
1365 }     
1366
1367 //______________________________________________________________________________
1368 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1369 {
1370 // Adds a user task to the global list of tasks.
1371    if (fInitOK) {
1372       Error("AddTask", "Cannot add task %s since InitAnalysis was already called", task->GetName());
1373       return;
1374    }   
1375       
1376    if (fTasks->FindObject(task)) {
1377       Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1378       return;
1379    }   
1380    task->SetActive(kFALSE);
1381    fTasks->Add(task);
1382 }  
1383
1384 //______________________________________________________________________________
1385 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1386 {
1387 // Retreive task by name.
1388    if (!fTasks) return NULL;
1389    return (AliAnalysisTask*)fTasks->FindObject(name);
1390 }
1391
1392 //______________________________________________________________________________
1393 Int_t AliAnalysisManager::GetTaskIndex(const AliAnalysisTask *task) const
1394 {
1395 // Returns task inded in the manager's list, -1 if not registered.
1396    if (!fTasks) return -1;
1397    return fTasks->IndexOf(task);
1398 }  
1399
1400 //______________________________________________________________________________
1401 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
1402                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
1403 {
1404 // Create a data container of a certain type. Types can be:
1405 //   kExchangeContainer  = 0, used to exchange data between tasks
1406 //   kInputContainer   = 1, used to store input data
1407 //   kOutputContainer  = 2, used for writing result to a file
1408 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1409 // the output object to a folder inside the output file
1410    if (fContainers->FindObject(name)) {
1411       Error("CreateContainer","A container named %s already defined !",name);
1412       return NULL;
1413    }   
1414    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1415    fContainers->Add(cont);
1416    switch (type) {
1417       case kInputContainer:
1418          fInputs->Add(cont);
1419          break;
1420       case kOutputContainer:
1421          fOutputs->Add(cont);
1422          if (filename && strlen(filename)) {
1423             cont->SetFileName(filename);
1424             cont->SetDataOwned(kFALSE);  // data owned by the file
1425          }   
1426          break;
1427       case kParamContainer:
1428          fParamCont->Add(cont);
1429          if (filename && strlen(filename)) {
1430             cont->SetFileName(filename);
1431             cont->SetDataOwned(kFALSE);  // data owned by the file
1432          }   
1433          break;
1434       case kExchangeContainer:
1435          break;   
1436    }
1437    return cont;
1438 }
1439          
1440 //______________________________________________________________________________
1441 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1442                                         AliAnalysisDataContainer *cont)
1443 {
1444 // Connect input of an existing task to a data container.
1445    if (!task) {
1446       Error("ConnectInput", "Task pointer is NULL");
1447       return kFALSE;
1448    }   
1449    if (!fTasks->FindObject(task)) {
1450       AddTask(task);
1451       Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1452    } 
1453    Bool_t connected = task->ConnectInput(islot, cont);
1454    return connected;
1455 }   
1456
1457 //______________________________________________________________________________
1458 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1459                                         AliAnalysisDataContainer *cont)
1460 {
1461 // Connect output of an existing task to a data container.
1462    if (!task) {
1463       Error("ConnectOutput", "Task pointer is NULL");
1464       return kFALSE;
1465    }   
1466    if (!fTasks->FindObject(task)) {
1467       AddTask(task);
1468       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1469    } 
1470    Bool_t connected = task->ConnectOutput(islot, cont);
1471    return connected;
1472 }   
1473                                
1474 //______________________________________________________________________________
1475 void AliAnalysisManager::CleanContainers()
1476 {
1477 // Clean data from all containers that have already finished all client tasks.
1478    TIter next(fContainers);
1479    AliAnalysisDataContainer *cont;
1480    while ((cont=(AliAnalysisDataContainer *)next())) {
1481       if (cont->IsOwnedData() && 
1482           cont->IsDataReady() && 
1483           cont->ClientsExecuted()) cont->DeleteData();
1484    }
1485 }
1486
1487 //______________________________________________________________________________
1488 Bool_t AliAnalysisManager::InitAnalysis()
1489 {
1490 // Initialization of analysis chain of tasks. Should be called after all tasks
1491 // and data containers are properly connected
1492    // Reset flag and remove valid_outputs file if exists
1493    if (fInitOK) return kTRUE;
1494    if (!gSystem->AccessPathName("outputs_valid"))
1495       gSystem->Unlink("outputs_valid");
1496    // Check for top tasks (depending only on input data containers)
1497    if (!fTasks->First()) {
1498       Error("InitAnalysis", "Analysis has no tasks !");
1499       return kFALSE;
1500    }   
1501    TIter next(fTasks);
1502    AliAnalysisTask *task;
1503    AliAnalysisDataContainer *cont;
1504    Int_t ntop = 0;
1505    Int_t nzombies = 0;
1506    Bool_t iszombie = kFALSE;
1507    Bool_t istop = kTRUE;
1508    Int_t i;
1509    while ((task=(AliAnalysisTask*)next())) {
1510       istop = kTRUE;
1511       iszombie = kFALSE;
1512       Int_t ninputs = task->GetNinputs();
1513       for (i=0; i<ninputs; i++) {
1514          cont = task->GetInputSlot(i)->GetContainer();
1515          if (!cont) {
1516             if (!iszombie) {
1517                task->SetZombie();
1518                fZombies->Add(task);
1519                nzombies++;
1520                iszombie = kTRUE;
1521             }   
1522             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
1523                   i, task->GetName()); 
1524          }
1525          if (iszombie) continue;
1526          // Check if cont is an input container
1527          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1528          // Connect to parent task
1529       }
1530       if (istop) {
1531          ntop++;
1532          fTopTasks->Add(task);
1533       }
1534    }
1535    if (!ntop) {
1536       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1537       return kFALSE;
1538    }                        
1539    // Check now if there are orphan tasks
1540    for (i=0; i<ntop; i++) {
1541       task = (AliAnalysisTask*)fTopTasks->At(i);
1542       task->SetUsed();
1543    }
1544    Int_t norphans = 0;
1545    next.Reset();
1546    while ((task=(AliAnalysisTask*)next())) {
1547       if (!task->IsUsed()) {
1548          norphans++;
1549          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1550       }   
1551    }          
1552    // Check the task hierarchy (no parent task should depend on data provided
1553    // by a daughter task)
1554    for (i=0; i<ntop; i++) {
1555       task = (AliAnalysisTask*)fTopTasks->At(i);
1556       if (task->CheckCircularDeps()) {
1557          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1558          PrintStatus("dep");
1559          return kFALSE;
1560       }   
1561    }
1562    // Check that all containers feeding post-event loop tasks are in the outputs list
1563    TIter nextcont(fContainers); // loop over all containers
1564    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1565       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1566          if (cont->HasConsumers()) {
1567          // Check if one of the consumers is post event loop
1568             TIter nextconsumer(cont->GetConsumers());
1569             while ((task=(AliAnalysisTask*)nextconsumer())) {
1570                if (task->IsPostEventLoop()) {
1571                   fOutputs->Add(cont);
1572                   break;
1573                }
1574             }
1575          }
1576       }
1577    }   
1578    // Check if all special output containers have a file name provided
1579    TIter nextout(fOutputs);
1580    while ((cont=(AliAnalysisDataContainer*)nextout())) {
1581       if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1582          Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1583          return kFALSE;
1584       }
1585    }
1586    // Initialize requested branch list if needed
1587    if (!fAutoBranchHandling) {
1588       next.Reset();
1589       while ((task=(AliAnalysisTask*)next())) {
1590          if (!task->HasBranches()) {
1591             Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches.\nUse: fBranchNames = \"ESD:br1,br2,...,brN AOD:bra1,bra2,...,braM\"",
1592                   task->GetName(), task->ClassName());
1593             return kFALSE;
1594          }
1595          if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1596             Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1597             return kFALSE;
1598          }
1599          TString taskbranches;
1600          task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1601          if (taskbranches.IsNull()) {
1602             Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1603                   task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1604             return kFALSE;      
1605          }
1606          AddBranches(taskbranches);
1607       }         
1608    }
1609    fInitOK = kTRUE;
1610    return kTRUE;
1611 }   
1612
1613 //______________________________________________________________________________
1614 void AliAnalysisManager::AddBranches(const char *branches)
1615 {
1616 // Add branches to the existing fRequestedBranches.
1617    TString br(branches);
1618    TObjArray *arr = br.Tokenize(",");
1619    TIter next(arr);
1620    TObject *obj;
1621    while ((obj=next())) {
1622       if (!fRequestedBranches.Contains(obj->GetName())) {
1623          if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1624          fRequestedBranches += obj->GetName();
1625       }
1626    }
1627   delete arr;
1628 }   
1629
1630 //______________________________________________________________________________
1631 void AliAnalysisManager::CheckBranches(Bool_t load)
1632 {
1633 // The method checks the input branches to be loaded during the analysis.
1634    if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;   
1635    TObjArray *arr = fRequestedBranches.Tokenize(",");
1636    TIter next(arr);
1637    TObject *obj;
1638    while ((obj=next())) {
1639       TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1640       if (!br) {
1641          br = fTree->GetBranch(obj->GetName());
1642          if (!br) {
1643             Error("CheckBranches", "Could not find branch %s",obj->GetName());
1644             continue;
1645          }
1646          fTable.Add(br);
1647       }   
1648       if (load && br->GetReadEntry()!=GetCurrentEntry()) {
1649          br->GetEntry(GetCurrentEntry());
1650       }      
1651    }
1652   delete arr;
1653 }
1654
1655 //______________________________________________________________________________
1656 Bool_t AliAnalysisManager::CheckTasks() const
1657 {
1658 // Check consistency of tasks.
1659    Int_t ntasks = fTasks->GetEntries();
1660    if (!ntasks) {
1661       Error("CheckTasks", "No tasks connected to the manager. This may be due to forgetting to compile the task or to load their library.");
1662       return kFALSE;
1663    }
1664    // Get the pointer to AliAnalysisTaskSE::Class()
1665    TClass *badptr = (TClass*)gROOT->ProcessLine("AliAnalysisTaskSE::Class()");
1666    // Loop all tasks to check if their corresponding library was loaded
1667    TIter next(fTasks);
1668    TObject *obj;
1669    while ((obj=next())) {
1670       if (obj->IsA() == badptr) {
1671          Error("CheckTasks", "##################\n \
1672          Class for task %s NOT loaded. You probably forgot to load the library for this task (or compile it dynamically).\n###########################\n",obj->GetName());
1673          return kFALSE;
1674       }
1675    }
1676    return kTRUE;      
1677 }   
1678
1679 //______________________________________________________________________________
1680 void AliAnalysisManager::PrintStatus(Option_t *option) const
1681 {
1682 // Print task hierarchy.
1683    if (!fInitOK) {
1684       Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1685       return;
1686    }   
1687    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1688    if (getsysInfo)
1689       Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1690    TIter next(fTopTasks);
1691    AliAnalysisTask *task;
1692    while ((task=(AliAnalysisTask*)next()))
1693       task->PrintTask(option);
1694   
1695    if (!fAutoBranchHandling && !fRequestedBranches.IsNull()) 
1696       printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1697   
1698   TString sopt(option);
1699   sopt.ToUpper();
1700   
1701   if (sopt.Contains("ALL"))
1702   {
1703     if ( fOutputEventHandler )
1704     {
1705       cout << TString('_',78) << endl;
1706       cout << "OutputEventHandler:" << endl;
1707       fOutputEventHandler->Print("   ");
1708     }
1709   }
1710 }
1711
1712 //______________________________________________________________________________
1713 void AliAnalysisManager::ResetAnalysis()
1714 {
1715 // Reset all execution flags and clean containers.
1716    CleanContainers();
1717 }
1718
1719 //______________________________________________________________________________
1720 void AliAnalysisManager::RunLocalInit()
1721 {
1722 // Run LocalInit method for all tasks.
1723    TDirectory *cdir = gDirectory;
1724    if (IsTrainInitialized()) return;
1725    TIter nextTask(fTasks);
1726    AliAnalysisTask *task;
1727    while ((task=(AliAnalysisTask*)nextTask())) {
1728       gROOT->cd();
1729       task->LocalInit();
1730    }
1731    if (cdir) cdir->cd();
1732    TObject::SetBit(kTasksInitialized, kTRUE);
1733 }   
1734
1735 //______________________________________________________________________________
1736 void AliAnalysisManager::InputFileFromTree(TTree * const tree, TString &fname)
1737 {
1738 // Retrieves name of the file from tree
1739    fname = "";
1740    if (!tree) return;
1741    TFile *file = tree->GetCurrentFile();
1742    TString basename;
1743    if (!file) {
1744       TChain *chain = dynamic_cast<TChain*>(tree);
1745       if (!chain || !chain->GetNtrees()) return;
1746       basename = gSystem->BaseName(chain->GetListOfFiles()->First()->GetTitle());
1747    } else {   
1748       basename = gSystem->BaseName(file->GetName());
1749    }   
1750    Int_t index = basename.Index("#");
1751    fname = basename(index+1, basename.Length());
1752 }   
1753
1754 //______________________________________________________________________________
1755 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1756 {
1757 // Start analysis having a grid handler.
1758    if (!fGridHandler) {
1759       Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1760       Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1761       return -1;
1762    }
1763    TTree *tree = NULL;
1764    return StartAnalysis(type, tree, nentries, firstentry);
1765 }
1766
1767 //______________________________________________________________________________
1768 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1769 {
1770 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1771 // MIX. Process nentries starting from firstentry
1772    Long64_t retv = 0;
1773    // Backup current directory and make sure gDirectory points to gROOT
1774    TDirectory *cdir = gDirectory;
1775    gROOT->cd();
1776    if (!fInitOK) {
1777       Error("StartAnalysis","Analysis manager was not initialized !");
1778       if (cdir) cdir->cd();
1779       return -1;
1780    }
1781    if (!CheckTasks()) Fatal("StartAnalysis", "Not all needed libraries were loaded");
1782    if (fDebug > 1) {
1783       printf("StartAnalysis %s\n",GetName());
1784       AliLog::SetGlobalLogLevel(AliLog::kInfo);
1785    }
1786    fMaxEntries = nentries;
1787    fIsRemote = kFALSE;
1788    TString anaType = type;
1789    anaType.ToLower();
1790    fMode = kLocalAnalysis;
1791    if (anaType.Contains("file"))      fIsRemote = kTRUE;
1792    if (anaType.Contains("proof"))     fMode = kProofAnalysis;
1793    else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1794    else if (anaType.Contains("mix"))  fMode = kMixingAnalysis;
1795    if (fInputEventHandler) {
1796       TString fname;
1797       InputFileFromTree(tree, fname);
1798       if (fname.Length()) fInputEventHandler->SetInputFileName(fname);
1799    }
1800
1801    if (fMode == kGridAnalysis) {
1802       fIsRemote = kTRUE;
1803       if (!anaType.Contains("terminate")) {
1804          if (!fGridHandler) {
1805             Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1806             Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1807             if (cdir) cdir->cd();
1808             return -1;
1809          }
1810          // Write analysis manager in the analysis file
1811          cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1812          // run local task configuration
1813          RunLocalInit();
1814          if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1815             Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1816             if (cdir) cdir->cd();
1817             return -1;
1818          }   
1819
1820          // Terminate grid analysis
1821          if (fSelector && fSelector->GetStatus() == -1) {if (cdir) cdir->cd(); return -1;}
1822          if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {if (cdir) cdir->cd(); return 0;}
1823          cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1824          if (!fGridHandler->MergeOutputs()) {
1825             // Return if outputs could not be merged or if it alien handler
1826             // was configured for offline mode or local testing.
1827             if (cdir) cdir->cd();
1828             return 0;
1829          }
1830       }   
1831       cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1832       if (cdir) cdir->cd();
1833       ImportWrappers(NULL);
1834       Terminate();
1835       if (cdir) cdir->cd();
1836       return 0;
1837    }
1838    TString line;
1839    SetEventLoop(kFALSE);
1840    // Enable event loop mode if a tree was provided
1841    if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1842
1843    TChain *chain = 0;
1844    TString ttype = "TTree";
1845    if (tree && tree->IsA() == TChain::Class()) {
1846       chain = (TChain*)tree;
1847       if (!chain || !chain->GetListOfFiles()->First()) {
1848          Error("StartAnalysis", "Cannot process null or empty chain...");
1849          if (cdir) cdir->cd();
1850          return -1;
1851       }   
1852       ttype = "TChain";
1853    }   
1854
1855    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1856    if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1857    // Initialize locally all tasks (happens for all modes)
1858    TIter next(fTasks);
1859    AliAnalysisTask *task;
1860    RunLocalInit();
1861    
1862    switch (fMode) {
1863       case kLocalAnalysis:
1864          if (!tree && !fGridHandler) {
1865             TIter nextT(fTasks);
1866             // Call CreateOutputObjects for all tasks
1867             Int_t itask = 0;
1868             Bool_t dirStatus = TH1::AddDirectoryStatus();
1869             while ((task=(AliAnalysisTask*)nextT())) {
1870                TH1::AddDirectory(kFALSE);
1871                task->CreateOutputObjects();
1872                if (!task->CheckPostData()) {
1873                   Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
1874                         Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
1875                         ########### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ###########", task->GetName(), task->ClassName());
1876                }
1877                if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1878                gROOT->cd();
1879                itask++;
1880             }   
1881             TH1::AddDirectory(dirStatus);
1882             if (IsExternalLoop()) {
1883                Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1884                      \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1885                return 0;
1886             }
1887             ExecAnalysis();
1888             Terminate();
1889             return 0;
1890          } 
1891          fSelector = new AliAnalysisSelector(this);
1892          // Check if a plugin handler is used
1893          if (fGridHandler) {
1894             // Get the chain from the plugin
1895             TString dataType = "esdTree";
1896             if (fInputEventHandler) {
1897                dataType = fInputEventHandler->GetDataType();
1898                dataType.ToLower();
1899                dataType += "Tree";
1900             }   
1901             chain = fGridHandler->GetChainForTestMode(dataType);
1902             if (!chain) {
1903                Error("StartAnalysis", "No chain for test mode. Aborting.");
1904                return -1;
1905             }
1906             cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1907             retv = chain->Process(fSelector, "", nentries, firstentry);
1908             break;
1909          }
1910          // Run tree-based analysis via AliAnalysisSelector  
1911          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1912          retv = tree->Process(fSelector, "", nentries, firstentry);
1913          break;
1914       case kProofAnalysis:
1915          fIsRemote = kTRUE;
1916          // Check if the plugin is used
1917          if (fGridHandler) {
1918             return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1919          }
1920          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1921             Error("StartAnalysis", "No PROOF!!! Exiting.");
1922             if (cdir) cdir->cd();
1923             return -1;
1924          }   
1925          line = Form("gProof->AddInput((TObject*)%p);", this);
1926          gROOT->ProcessLine(line);
1927          if (chain) {
1928             chain->SetProof();
1929             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1930             retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1931          } else {
1932             Error("StartAnalysis", "No chain!!! Exiting.");
1933             if (cdir) cdir->cd();
1934             return -1;
1935          }      
1936          break;
1937       case kGridAnalysis:
1938          fIsRemote = kTRUE;
1939          if (!anaType.Contains("terminate")) {
1940             if (!fGridHandler) {
1941                Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1942                Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1943                if (cdir) cdir->cd();
1944                return -1;
1945             }
1946             // Write analysis manager in the analysis file
1947             cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1948             // Start the analysis via the handler
1949             if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1950                Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1951                if (cdir) cdir->cd();
1952                return -1;
1953             }   
1954
1955             // Terminate grid analysis
1956             if (fSelector && fSelector->GetStatus() == -1) {if (cdir) cdir->cd(); return -1;}
1957             if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {if (cdir) cdir->cd(); return 0;}
1958             cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1959             if (!fGridHandler->MergeOutputs()) {
1960                // Return if outputs could not be merged or if it alien handler
1961                // was configured for offline mode or local testing.
1962                if (cdir) cdir->cd();
1963                return 0;
1964             }
1965          }   
1966          cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1967          ImportWrappers(NULL);
1968          Terminate();
1969          if (cdir) cdir->cd();
1970          return 0;
1971       case kMixingAnalysis:   
1972          // Run event mixing analysis
1973          if (!fEventPool) {
1974             Error("StartAnalysis", "Cannot run event mixing without event pool");
1975             if (cdir) cdir->cd();
1976             return -1;
1977          }
1978          cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1979          fSelector = new AliAnalysisSelector(this);
1980          while ((chain=fEventPool->GetNextChain())) {
1981             next.Reset();
1982             // Call NotifyBinChange for all tasks
1983             while ((task=(AliAnalysisTask*)next()))
1984                if (!task->IsPostEventLoop()) task->NotifyBinChange();
1985             retv = chain->Process(fSelector);
1986             if (retv < 0) {
1987                Error("StartAnalysis", "Mixing analysis failed");
1988                if (cdir) cdir->cd();
1989                return retv;
1990             }   
1991          }
1992          PackOutput(fSelector->GetOutputList());
1993          Terminate();
1994    }
1995    if (cdir) cdir->cd();
1996    return retv;
1997 }   
1998
1999 //______________________________________________________________________________
2000 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
2001 {
2002 // Start analysis for this manager on a given dataset. Analysis task can be: 
2003 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
2004    if (!fInitOK) {
2005       Error("StartAnalysis","Analysis manager was not initialized !");
2006       return -1;
2007    }
2008    fIsRemote = kTRUE;
2009    if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
2010    TString anaType = type;
2011    anaType.ToLower();
2012    if (!anaType.Contains("proof")) {
2013       Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
2014       return -1;
2015    }   
2016    fMode = kProofAnalysis;
2017    TString line;
2018    TString proofProcessOpt;
2019    SetEventLoop(kTRUE);
2020    // Set the dataset flag
2021    TObject::SetBit(kUseDataSet);
2022    fTree = 0;
2023    if (fGridHandler) {
2024       // Start proof analysis using the grid handler
2025       if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
2026          Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
2027          return -1;
2028       }
2029       // Check if the plugin is in test mode
2030       if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
2031          dataset = "test_collection";
2032       } else {
2033          dataset = fGridHandler->GetProofDataSet();
2034       }
2035
2036       proofProcessOpt = fGridHandler->GetProofProcessOpt();
2037    }   
2038
2039    if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
2040       Error("StartAnalysis", "No PROOF!!! Exiting.");
2041       return -1;
2042    }   
2043
2044    // Initialize locally all tasks
2045    RunLocalInit();
2046       
2047    line = Form("gProof->AddInput((TObject*)%p);", this);
2048    gROOT->ProcessLine(line);
2049    Long_t retv;
2050    line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"%s\", %lld, %lld);",
2051                dataset,proofProcessOpt.Data(), nentries, firstentry);
2052    cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
2053    retv = (Long_t)gROOT->ProcessLine(line);
2054    return retv;
2055 }   
2056
2057 //______________________________________________________________________________
2058 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
2059 {
2060 // Opens according the option the file specified by cont->GetFileName() and changes
2061 // current directory to cont->GetFolderName(). If the file was already opened, it
2062 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
2063 // be optionally ignored.
2064   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2065   TString filename = cont->GetFileName();
2066   TFile *f = NULL;
2067   if (filename.IsNull()) {
2068     ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
2069     return NULL;
2070   }
2071   if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
2072       && !ignoreProof)
2073     f = mgr->OpenProofFile(cont,option);
2074   else {
2075     // Check first if the file is already opened
2076     f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2077     if (f) {
2078       // Check if option "UPDATE" was preserved 
2079       TString opt(option);
2080       opt.ToUpper();
2081       if ((opt=="UPDATE") && (opt!=f->GetOption())) 
2082         ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
2083     } else {
2084       f = TFile::Open(filename, option);
2085     }    
2086   }   
2087   if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
2088     cont->SetFile(f);
2089     // Cd to file
2090     f->cd();
2091     // Check for a folder request
2092     TString dir = cont->GetFolderName(); 
2093     if (!dir.IsNull()) {
2094       if (!f->GetDirectory(dir)) f->mkdir(dir);
2095       f->cd(dir);
2096     }
2097     return f;
2098   }
2099   ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
2100   cont->SetFile(NULL);
2101   return NULL;
2102 }    
2103  
2104 //______________________________________________________________________________
2105 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
2106 {
2107 // Opens a special output file used in PROOF.
2108   TString line;
2109   TString filename = cont->GetFileName();
2110   if (cont == fCommonOutput) {
2111      if (fOutputEventHandler) {
2112         if (strlen(extaod)) filename = extaod;
2113         filename = fOutputEventHandler->GetOutputFileName();
2114      }   
2115      else Fatal("OpenProofFile","No output container. Exiting.");
2116   }   
2117   TFile *f = NULL;
2118   if (fMode!=kProofAnalysis || !fSelector) {
2119     Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
2120     return NULL;
2121   } 
2122   if (fSpecialOutputLocation.Length()) {
2123     f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2124     if (f) {
2125       // Check if option "UPDATE" was preserved 
2126       TString opt(option);
2127       opt.ToUpper();
2128       if ((opt=="UPDATE") && (opt!=f->GetOption()))
2129         ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
2130     } else {
2131       f = new TFile(filename, option);
2132     }
2133     if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
2134       cont->SetFile(f);
2135       // Cd to file
2136       f->cd();
2137       // Check for a folder request
2138       TString dir = cont->GetFolderName(); 
2139       if (dir.Length()) {
2140         if (!f->GetDirectory(dir)) f->mkdir(dir);
2141         f->cd(dir);
2142       }      
2143       return f;
2144     }
2145     Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
2146     cont->SetFile(NULL);
2147     return NULL;       
2148   }
2149   // Check if there is already a proof output file in the output list
2150   TObject *pof = fSelector->GetOutputList()->FindObject(filename);
2151   if (pof) {
2152     // Get the actual file
2153     line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
2154     filename = (const char*)gROOT->ProcessLine(line);
2155     if (fDebug>1) {
2156       printf("File: %s already booked via TProofOutputFile\n", filename.Data());
2157     }  
2158     f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2159     if (!f) {
2160        Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
2161        return NULL;
2162     }   
2163     // Check if option "UPDATE" was preserved 
2164     TString opt(option);
2165     opt.ToUpper();
2166     if ((opt=="UPDATE") && (opt!=f->GetOption())) 
2167       Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
2168   } else {
2169     if (cont->IsRegisterDataset()) {
2170       TString dsetName = filename;
2171       dsetName.ReplaceAll(".root", cont->GetTitle());
2172       dsetName.ReplaceAll(":","_");
2173       if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
2174       line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
2175     } else {
2176       if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
2177       line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
2178     }
2179     if (fDebug > 1) printf("=== %s\n", line.Data());
2180     gROOT->ProcessLine(line);
2181     line = Form("pf->OpenFile(\"%s\");", option);
2182     gROOT->ProcessLine(line);
2183     f = gFile;
2184     if (fDebug > 1) {
2185       gROOT->ProcessLine("pf->Print()");
2186       printf(" == proof file name: %s", f->GetName());
2187     }   
2188     // Add to proof output list
2189     line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
2190     if (fDebug > 1) printf("=== %s\n", line.Data());
2191     gROOT->ProcessLine(line);
2192   }
2193   if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
2194     cont->SetFile(f);
2195     // Cd to file
2196     f->cd();
2197     // Check for a folder request
2198     TString dir = cont->GetFolderName(); 
2199     if (!dir.IsNull()) {
2200       if (!f->GetDirectory(dir)) f->mkdir(dir);
2201       f->cd(dir);
2202     }
2203     return f;
2204   }
2205   Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
2206   cont->SetFile(NULL);  
2207   return NULL;
2208 }   
2209
2210 //______________________________________________________________________________
2211 void AliAnalysisManager::ExecAnalysis(Option_t *option)
2212 {
2213 // Execute analysis.
2214    static Long64_t nentries = 0;
2215    static TTree *lastTree = 0;
2216    static TStopwatch *timer = new TStopwatch();
2217    // Only the first call to Process will trigger a true Notify. Other Notify
2218    // coming before is ignored.
2219    if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) {
2220       TObject::SetBit(AliAnalysisManager::kTrueNotify);
2221       Notify();
2222    }   
2223    if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
2224    else {
2225       if (fTree && (fTree != lastTree)) {
2226          nentries += fTree->GetEntries();
2227          lastTree = fTree;
2228       }   
2229       if (!fNcalls) timer->Start();
2230       if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, TMath::Min(fMaxEntries,nentries), timer, kFALSE);
2231    }
2232    fIOTimer->Start(kTRUE);
2233    gROOT->cd();
2234    TDirectory *cdir = gDirectory;
2235    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
2236    if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
2237    if (!fInitOK) {
2238       Error("ExecAnalysis", "Analysis manager was not initialized !");
2239       if (cdir) cdir->cd();
2240       return;
2241    }
2242    fNcalls++;
2243    AliAnalysisTask *task;
2244    // Check if the top tree is active.
2245    if (fTree) {
2246       if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2247          AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
2248       TIter next(fTasks);
2249    // De-activate all tasks
2250       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
2251       AliAnalysisDataContainer *cont = fCommonInput;
2252       if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
2253       if (!cont) {
2254               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
2255          if (cdir) cdir->cd();
2256          return;
2257       }   
2258       cont->SetData(fTree); // This will notify all consumers
2259       Long64_t entry = fTree->GetTree()->GetReadEntry();      
2260 //
2261 //    Call BeginEvent() for optional input/output and MC services 
2262       if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
2263       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
2264       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
2265       gROOT->cd();
2266       if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2267          AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2268 //
2269 //    Execute the tasks
2270 //      TIter next1(cont->GetConsumers());
2271       fIOTimer->Stop();
2272       fIOTime += fIOTimer->RealTime();
2273       fCPUTimer->Start(kTRUE);
2274       TIter next1(fTopTasks);
2275       Int_t itask = 0;
2276       while ((task=(AliAnalysisTask*)next1())) {
2277          if (fDebug >1) {
2278             cout << "    Executing task " << task->GetName() << endl;
2279          }
2280          if (fStatistics) fStatistics->StartTimer(GetTaskIndex(task), task->GetName(), task->ClassName());
2281          task->ExecuteTask(option);
2282          if (fStatistics) fStatistics->StopTimer();
2283          gROOT->cd();
2284          if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2285             AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
2286          itask++;   
2287       }
2288       fCPUTimer->Stop();
2289       fCPUTime += fCPUTimer->RealTime();
2290       fIOTimer->Start(kTRUE);
2291 //
2292 //    Call FinishEvent() for optional output and MC services 
2293       if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
2294       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
2295       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2296       // Gather system information if requested
2297       if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2298          AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
2299       if (cdir) cdir->cd();   
2300       fIOTimer->Stop();
2301       fIOTime += fIOTimer->RealTime();
2302       return;
2303    }   
2304    // The event loop is not controlled by TSelector   
2305 //
2306 //  Call BeginEvent() for optional input/output and MC services 
2307    fIOTimer->Start(kTRUE);
2308    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(-1);
2309    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(-1);
2310    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
2311    fIOTimer->Stop();
2312    fIOTime += fIOTimer->RealTime();
2313    gROOT->cd();
2314    if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2315       AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2316    fCPUTimer->Start(kTRUE);
2317    TIter next2(fTopTasks);
2318    while ((task=(AliAnalysisTask*)next2())) {
2319       task->SetActive(kTRUE);
2320       if (fDebug > 1) {
2321          cout << "    Executing task " << task->GetName() << endl;
2322       }   
2323       if (fStatistics) fStatistics->StartTimer(GetTaskIndex(task), task->GetName(), task->ClassName());
2324       task->ExecuteTask(option);
2325       if (fStatistics) fStatistics->StopTimer();
2326       gROOT->cd();
2327    }   
2328    fCPUTimer->Stop();
2329    fCPUTime += fCPUTimer->RealTime();
2330 //
2331 // Call FinishEvent() for optional output and MC services 
2332    fIOTimer->Start(kTRUE);
2333    if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
2334    if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
2335    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2336    if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2337       AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2338    if (cdir) cdir->cd();   
2339    fIOTimer->Stop();
2340    fIOTime += fIOTimer->RealTime();
2341 }
2342
2343 //______________________________________________________________________________
2344 Bool_t AliAnalysisManager::IsPipe(std::ostream &out)
2345 {
2346 // Check if the stdout is connected to a pipe (C.Holm)
2347   Bool_t ispipe = kFALSE;
2348   out.seekp(0, std::ios_base::cur);
2349   if (out.fail()) {
2350     out.clear();
2351     if (errno == ESPIPE) ispipe = kTRUE;
2352   }
2353   return ispipe;
2354 }
2355    
2356 //______________________________________________________________________________
2357 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2358 {
2359 // Set the input event handler and create a container for it.
2360    Changed();
2361    fInputEventHandler   = handler;
2362    if (!fCommonInput) fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2363 }
2364
2365 //______________________________________________________________________________
2366 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2367 {
2368 // Set the input event handler and create a container for it.
2369    Changed();
2370    fOutputEventHandler   = handler;
2371    if (!fCommonOutput) fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2372    fCommonOutput->SetSpecialOutput();
2373 }
2374
2375 //______________________________________________________________________________
2376 void AliAnalysisManager::SetDebugLevel(UInt_t level)
2377 {
2378 // Set verbosity of the analysis manager. If the progress bar is used, the call is ignored
2379    if (TObject::TestBit(kUseProgressBar)) {
2380       Info("SetDebugLevel","Ignored. Disable the progress bar first.");
2381       return;
2382    }
2383    fDebug = level;
2384 }
2385    
2386 //______________________________________________________________________________
2387 void AliAnalysisManager::SetUseProgressBar(Bool_t flag, Int_t freq)
2388 {
2389 // Enable a text mode progress bar. Resets debug level to 0.
2390    Info("SetUseProgressBar", "Progress bar enabled, updated every %d events.\n  ### NOTE: Debug level reset to 0 ###", freq);
2391    TObject::SetBit(kUseProgressBar,flag);
2392    fPBUpdateFreq = freq;
2393    fDebug = 0;
2394 }   
2395
2396 //______________________________________________________________________________
2397 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2398 {
2399 // This method is used externally to register output files which are not
2400 // connected to any output container, so that the manager can properly register,
2401 // retrieve or merge them when running in distributed mode. The file names are
2402 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2403    if (fExtraFiles.Contains(fname)) return;
2404    if (fExtraFiles.Length()) fExtraFiles += " ";
2405    fExtraFiles += fname;
2406 }
2407
2408 //______________________________________________________________________________
2409 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2410 {
2411 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2412    char fullPath[512];
2413    char chUrl[512];
2414    char tmp[1024];
2415    TObject *pof =  source->FindObject(filename);
2416    if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2417       Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2418       return kFALSE;
2419    }
2420    gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2421    gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2422    TString clientUrl(chUrl);
2423    TString fullPath_str(fullPath);
2424    if (clientUrl.Contains("localhost")){
2425       TObjArray* array = fullPath_str.Tokenize ( "//" );
2426       TObjString *strobj = ( TObjString *)array->At(1);
2427       TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2428       TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2429       fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2430       fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2431       if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2432       delete arrayPort;
2433       delete array;
2434    }
2435    else if (clientUrl.Contains("__lite__")) { 
2436      // Special case for ProofLite environement - get file info and copy. 
2437      gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2438      fullPath_str = Form("%s/%s", tmp, fullPath);
2439    }
2440    if (fDebug > 1) 
2441      Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2442    Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename); 
2443    if (!gotit)
2444       Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2445    return gotit;
2446 }
2447
2448 //______________________________________________________________________________
2449 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2450 {
2451 // Fill analysis type in the provided string.
2452    switch (fMode) {
2453       case kLocalAnalysis:
2454          type = "local";
2455          return;
2456       case kProofAnalysis:
2457          type = "proof";
2458          return;
2459       case kGridAnalysis:
2460          type = "grid";
2461          return;
2462       case kMixingAnalysis:
2463          type = "mix";
2464    }
2465 }
2466
2467 //______________________________________________________________________________
2468 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2469 {
2470 // Validate all output files.
2471    TIter next(fOutputs);
2472    AliAnalysisDataContainer *output;
2473    TDirectory *cdir = gDirectory;
2474    TString openedFiles;
2475    while ((output=(AliAnalysisDataContainer*)next())) {
2476       if (output->IsRegisterDataset()) continue;
2477       TString filename = output->GetFileName();
2478       if (filename == "default") {
2479          if (!fOutputEventHandler) continue;
2480          filename = fOutputEventHandler->GetOutputFileName();
2481          // Main AOD may not be there
2482          if (gSystem->AccessPathName(filename)) continue;
2483       }
2484       // Check if the file is closed
2485       if (openedFiles.Contains(filename)) continue;;
2486       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2487       if (file) {
2488          Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2489          // Clear file list to release object ownership to user.
2490 //         file->Clear();
2491          file->Close();
2492       }
2493       file = TFile::Open(filename);
2494       if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2495          Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2496          if (cdir) cdir->cd();
2497          return kFALSE;
2498       }
2499       file->Close();
2500       openedFiles += filename;
2501       openedFiles += " ";
2502    }
2503    if (cdir) cdir->cd();
2504    return kTRUE;
2505 }   
2506
2507 //______________________________________________________________________________
2508 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2509 {
2510 // Implements a nice text mode progress bar.
2511    static Long64_t icount = 0;
2512    static TString oname;
2513    static TString nname;
2514    static Long64_t ocurrent = 0;
2515    static Long64_t osize = 0;
2516    static Int_t oseconds = 0;
2517    static TStopwatch *owatch = 0;
2518    static Bool_t oneoftwo = kFALSE;
2519    static Int_t nrefresh = 0;
2520    static Int_t nchecks = 0;
2521    static char lastChar = 0;
2522    const char symbol[4] = {'-','\\','|','/'}; 
2523    
2524    if (!lastChar) lastChar = (IsPipe(std::cerr))?'\r':'\n';
2525    if (!refresh) {
2526       nrefresh = 0;
2527       if (!size) return;
2528       owatch = watch;
2529       oname = opname;
2530       ocurrent = TMath::Abs(current);
2531       osize = TMath::Abs(size);
2532       if (ocurrent > osize) ocurrent=osize;
2533    } else {
2534       nrefresh++;
2535       if (!osize) return;
2536    }     
2537    if ((current % fPBUpdateFreq) != 0) return;
2538    icount++;
2539    char progress[11] = "          ";
2540    Int_t ichar = icount%4;
2541    Double_t time = 0.;
2542    Int_t hours = 0;
2543    Int_t minutes = 0;
2544    Int_t seconds = 0;
2545    if (owatch && !last) {
2546       owatch->Stop();
2547       time = owatch->RealTime();
2548       seconds   = int(time) % 60;
2549       minutes   = (int(time) / 60) % 60;
2550       hours     = (int(time) / 60 / 60);
2551       if (refresh)  {
2552          if (oseconds==seconds) {
2553             owatch->Continue();
2554             return;
2555          }
2556          oneoftwo = !oneoftwo;   
2557       }
2558       oseconds = seconds;   
2559    }
2560    if (refresh && oneoftwo) {
2561       nname = oname;
2562       if (nchecks <= 0) nchecks = nrefresh+1;
2563       Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2564       oname = Form("     == %d%% ==", pctdone);
2565    }         
2566    Double_t percent = 100.0*ocurrent/osize;
2567    Int_t nchar = Int_t(percent/10);
2568    if (nchar>10) nchar=10;
2569    Int_t i;
2570    for (i=0; i<nchar; i++)  progress[i] = '=';
2571    progress[nchar] = symbol[ichar];
2572    for (i=nchar+1; i<10; i++) progress[i] = ' ';
2573    progress[10] = '\0';
2574    oname += "                    ";
2575    oname.Remove(20);
2576    if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2577    else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2578    else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2579    if (time>0.) {
2580      Int_t full   = Int_t(ocurrent > 0 ? 
2581                           time * (float(osize)/ocurrent) + .5 : 
2582                           99*3600+59*60+59); 
2583      Int_t remain = Int_t(full - time);
2584      Int_t rsec   = remain % 60;
2585      Int_t rmin   = (remain / 60) % 60;
2586      Int_t rhour  = (remain / 60 / 60);
2587      fprintf(stderr, "[%6.2f %%]   TIME %.2d:%.2d:%.2d  ETA %.2d:%.2d:%.2d%c",
2588              percent, hours, minutes, seconds, rhour, rmin, rsec, lastChar);
2589    }
2590    else fprintf(stderr, "[%6.2f %%]%c", percent, lastChar);
2591    if (refresh && oneoftwo) oname = nname;
2592    if (owatch) owatch->Continue();
2593    if (last) {
2594       icount = 0;
2595       owatch = 0;
2596       ocurrent = 0;
2597       osize = 0;
2598       oseconds = 0;
2599       oneoftwo = kFALSE;
2600       nrefresh = 0;
2601       fprintf(stderr, "\n");
2602    }   
2603 }
2604
2605 //______________________________________________________________________________
2606 void AliAnalysisManager::DoLoadBranch(const char *name) 
2607 {
2608   // Get tree and load branch if needed.
2609   static Long64_t crtEntry = -100;
2610
2611   if (fAutoBranchHandling || !fTree)
2612     return;
2613
2614   TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2615   if (!br) {
2616     br = fTree->GetBranch(name);
2617     if (!br) {
2618       Error("DoLoadBranch", "Could not find branch %s",name);
2619       return;
2620     }
2621     fTable.Add(br);
2622   }
2623   if (br->GetReadEntry()==fCurrentEntry) return;
2624   Long64_t readbytes = br->GetEntry(GetCurrentEntry());
2625   if (readbytes<0) {
2626     Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2627     if (crtEntry != fCurrentEntry) {
2628       CountEvent(1,0,1,0);
2629       crtEntry = fCurrentEntry;
2630     }  
2631   } else {
2632     if (crtEntry != fCurrentEntry) {
2633       CountEvent(1,1,0,0);
2634       crtEntry = fCurrentEntry;
2635     }
2636   }
2637 }
2638
2639 //______________________________________________________________________________
2640 void AliAnalysisManager::AddStatisticsTask(UInt_t offlineMask)
2641 {
2642 // Add the statistics task to the manager.
2643   if (fStatistics) {
2644      Info("AddStatisticsTask", "Already added");
2645      return;
2646   }
2647   TString line = Form("AliAnalysisTaskStat::AddToManager(%u);", offlineMask);
2648   gROOT->ProcessLine(line);
2649 }  
2650
2651 //______________________________________________________________________________
2652 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2653 {
2654 // Bookkeep current event;
2655    if (!fStatistics) return;
2656    fStatistics->AddInput(ninput);
2657    fStatistics->AddProcessed(nprocessed);
2658    fStatistics->AddFailed(nfailed);
2659    fStatistics->AddAccepted(naccepted);
2660 }   
2661
2662 //______________________________________________________________________________
2663 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2664 {
2665 // Add a line in the statistics message. If available, the statistics message is written
2666 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2667 // on the client.
2668    if (!strlen(line)) return;
2669    if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2670    fStatisticsMsg += line;
2671 }
2672
2673 //______________________________________________________________________________
2674 void AliAnalysisManager::WriteStatisticsMsg(Int_t)
2675 {
2676 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2677    static Bool_t done = kFALSE;
2678    if (done) return;
2679    done = kTRUE;
2680    if (!fStatistics) return;
2681    ofstream out;
2682    AddStatisticsMsg(Form("Number of input events:        %lld",fStatistics->GetNinput()));
2683    AddStatisticsMsg(Form("Number of processed events:    %lld",fStatistics->GetNprocessed()));      
2684    AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2685    AddStatisticsMsg(Form("Number of accepted events for mask %s: %lld", AliAnalysisStatistics::GetMaskAsString(fStatistics->GetOfflineMask()), fStatistics->GetNaccepted()));
2686    out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2687                  fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2688                  fStatistics->GetNaccepted()), ios::out);      
2689    out << fStatisticsMsg << endl;
2690    out.close();
2691 }
2692
2693 //______________________________________________________________________________
2694 const char* AliAnalysisManager::GetOADBPath()
2695 {
2696 // returns the path of the OADB
2697 // this static function just depends on environment variables
2698
2699    static TString oadbPath;
2700
2701    if (gSystem->Getenv("OADB_PATH"))
2702       oadbPath = gSystem->Getenv("OADB_PATH");
2703    else if (gSystem->Getenv("ALICE_ROOT"))
2704       oadbPath.Form("%s/OADB", gSystem->Getenv("ALICE_ROOT"));
2705    else
2706       ::Fatal("AliAnalysisManager::GetOADBPath", "Cannot figure out AODB path. Define ALICE_ROOT or OADB_PATH!");
2707       
2708    return oadbPath;
2709 }
2710
2711 //______________________________________________________________________________
2712 void AliAnalysisManager::SetGlobalStr(const char *key, const char *value)
2713 {
2714 // Define a custom string variable mapped to a global unique name. The variable
2715 // can be then retrieved by a given analysis macro via GetGlobalStr(key).
2716    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2717    if (!mgr) {
2718       ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2719       return;
2720    }   
2721    Bool_t valid = kFALSE;
2722    TString existing = AliAnalysisManager::GetGlobalStr(key, valid);
2723    if (valid) {
2724       ::Error("AliAnalysisManager::SetGlobalStr", "Global %s = %s already defined.", key, existing.Data());
2725       return;
2726    }
2727    mgr->GetGlobals()->Add(new TObjString(key), new TObjString(value));
2728 }
2729
2730 //______________________________________________________________________________
2731 const char *AliAnalysisManager::GetGlobalStr(const char *key, Bool_t &valid)
2732 {
2733 // Static method to retrieve a global variable defined via SetGlobalStr.
2734    valid = kFALSE;
2735    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2736    if (!mgr) return 0;
2737    TObject *value = mgr->GetGlobals()->GetValue(key);
2738    if (!value) return 0;
2739    valid = kTRUE;
2740    return value->GetName();
2741 }
2742
2743 //______________________________________________________________________________
2744 void AliAnalysisManager::SetGlobalInt(const char *key, Int_t value)
2745 {
2746 // Define a custom integer variable mapped to a global unique name. The variable
2747 // can be then retrieved by a given analysis macro via GetGlobalInt(key).
2748    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2749    if (!mgr) {
2750       ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2751       return;
2752    }   
2753    Bool_t valid = kFALSE;
2754    Int_t existing = AliAnalysisManager::GetGlobalInt(key, valid);
2755    if (valid) {
2756       ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %i already defined.", key, existing);
2757       return;
2758    }
2759    mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%i",value)));
2760 }
2761
2762 //______________________________________________________________________________
2763 Int_t AliAnalysisManager::GetGlobalInt(const char *key, Bool_t &valid)
2764 {
2765 // Static method to retrieve a global variable defined via SetGlobalInt.
2766    valid = kFALSE;
2767    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2768    if (!mgr) return 0;
2769    TObject *value = mgr->GetGlobals()->GetValue(key);
2770    if (!value) return 0;
2771    valid = kTRUE;
2772    TString s = value->GetName();
2773    return s.Atoi();
2774 }
2775
2776 //______________________________________________________________________________
2777 void AliAnalysisManager::SetGlobalDbl(const char *key, Double_t value)
2778 {
2779 // Define a custom double precision variable mapped to a global unique name. The variable
2780 // can be then retrieved by a given analysis macro via GetGlobalInt(key).
2781    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2782    if (!mgr) {
2783       ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2784       return;
2785    }   
2786    Bool_t valid = kFALSE;
2787    Double_t existing = AliAnalysisManager::GetGlobalDbl(key, valid);
2788    if (valid) {
2789       ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %g already defined.", key, existing);
2790       return;
2791    }
2792    mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%f.16",value)));
2793 }
2794
2795 //______________________________________________________________________________
2796 Double_t AliAnalysisManager::GetGlobalDbl(const char *key, Bool_t &valid)
2797 {
2798 // Static method to retrieve a global variable defined via SetGlobalDbl.
2799    valid = kFALSE;
2800    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2801    if (!mgr) return 0;
2802    TObject *value = mgr->GetGlobals()->GetValue(key);
2803    if (!value) return 0;
2804    valid = kTRUE;
2805    TString s = value->GetName();
2806    return s.Atof();
2807 }
2808
2809 //______________________________________________________________________________
2810 void AliAnalysisManager::AddClassDebug(const char *className, Int_t debugLevel)
2811 {
2812 // Sets Class debug level
2813
2814    if (!fDebugOptions) {
2815       fDebugOptions = new TObjArray();
2816       fDebugOptions->SetOwner(kTRUE);
2817    }
2818
2819    TNamed *debugOpt = (TNamed*)fDebugOptions->FindObject(className);
2820    if (!debugOpt) {
2821      AliInfo(TString::Format("Adding debug level %d for class %s",debugLevel,className).Data());
2822      fDebugOptions->Add(new TNamed(className,TString::Format("%d",debugLevel).Data()));
2823    } else {
2824       TString oldDebugStr = debugOpt->GetTitle();
2825       Int_t oldDebug = oldDebugStr.Atoi();
2826       if (debugLevel > oldDebug) {
2827          AliWarning(TString::Format("Overwriting debug level to %d class %s, because it is higher then previously set (%d).",debugLevel,className,oldDebug).Data());
2828          debugOpt->SetTitle(TString::Format("%d",debugLevel).Data());
2829       } else {
2830          AliWarning(TString::Format("Ignoring debug level to %d class %s, because it is smaller then previously set (%d).",debugLevel,className,oldDebug).Data());
2831       }
2832    }
2833 }
2834
2835 //______________________________________________________________________________
2836 void AliAnalysisManager::ApplyDebugOptions()
2837 {
2838 // Apply debug options
2839
2840    if (!fDebugOptions) return;
2841    
2842    TIter next(fDebugOptions);
2843    TNamed *debug;
2844    TString debugLevel;
2845    while ((debug=dynamic_cast<TNamed*>(next()))) {
2846       debugLevel = debug->GetTitle();
2847       AliInfo(TString::Format("Class=%s debulLevel=%d",debug->GetName(),debugLevel.Atoi()).Data());
2848       AliLog::SetClassDebugLevel(debug->GetName(), debugLevel.Atoi());
2849    }
2850 }
2851
2852 //______________________________________________________________________________
2853 Bool_t AliAnalysisManager::IsMacroLoaded(char filename)
2854 {
2855 // Check if a macro was loaded.
2856    return fgMacroNames.Contains(filename);
2857 }
2858    
2859 //______________________________________________________________________________
2860 Int_t AliAnalysisManager::LoadMacro(const char *filename, Int_t *error, Bool_t check)
2861 {
2862 // Redirection of gROOT->LoadMacro which makes sure the same macro is not loaded 
2863 // twice
2864    TString macroName = gSystem->BaseName(filename);
2865    // Strip appended +, ++, +g, +O
2866    Int_t index = macroName.Index("+");
2867    if (index>0) macroName.Remove(index);
2868    if (fgMacroNames.Contains(macroName)) {
2869       // Macro with the same name loaded already in this root session, do 
2870       // nothing
2871       error = 0;
2872       return 0;
2873    }
2874    Int_t ret = gROOT->LoadMacro(filename,error,check);
2875    // In case of error return the error code
2876    if (ret) return ret;
2877    // Append the macro name to the loaded macros list
2878    fgMacroNames += macroName;
2879    fgMacroNames += " ";
2880    return ret;
2881 }   
2882
2883 //______________________________________________________________________________
2884 void AliAnalysisManager::Lock()
2885 {
2886 // Security lock. This is to detect NORMAL user errors and not really to
2887 // protect against intentional hacks.
2888    if (fLocked) return;
2889    fLocked = kTRUE;
2890    if (fInputEventHandler)  fInputEventHandler->Lock();
2891    if (fOutputEventHandler) fOutputEventHandler->Lock();
2892    if (fMCtruthEventHandler) fMCtruthEventHandler->Lock();
2893    Info("Lock","====== ANALYSIS MANAGER LOCKED ======");
2894 }
2895
2896 //______________________________________________________________________________
2897 void AliAnalysisManager::UnLock()
2898 {
2899 // Verbose unlocking. Hackers will be punished ;-) ... 
2900    if (!fLocked) return;
2901    fLocked = kFALSE;
2902    if (fInputEventHandler)  fInputEventHandler->UnLock();
2903    if (fOutputEventHandler) fOutputEventHandler->UnLock();
2904    if (fMCtruthEventHandler) fMCtruthEventHandler->UnLock();
2905    Info("UnLock", "====== ANALYSIS MANAGER UNLOCKED ======");
2906 }
2907
2908 //______________________________________________________________________________
2909 void AliAnalysisManager::Changed()
2910 {
2911 // All critical setters pass through the Changed method that throws an exception 
2912 // in case the lock was set.
2913    if (fLocked) Fatal("Changed","Critical setter called in locked mode");
2914 }