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