]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
Added job information to fileinfo.log: total duration, CPU time reported by ROOT...
[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       }
821    }
822               
823    if (cdir) cdir->cd();
824    if (fDebug > 1) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
825 }
826
827 //______________________________________________________________________________
828 void AliAnalysisManager::ImportWrappers(TList *source)
829 {
830 // Import data in output containers from wrappers coming in source.
831    if (fDebug > 1) printf("->AliAnalysisManager::ImportWrappers()\n");
832    fIOTimer->Start(kTRUE);
833    TIter next(fOutputs);
834    AliAnalysisDataContainer *cont;
835    AliAnalysisDataWrapper   *wrap;
836    Int_t icont = 0;
837    Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
838    TDirectory *cdir = gDirectory;
839    while ((cont=(AliAnalysisDataContainer*)next())) {
840       wrap = 0;
841       if (cont->GetProducer() && cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
842       if (cont->IsRegisterDataset()) continue;
843       const char *filename = cont->GetFileName();
844       Bool_t isManagedByHandler = kFALSE;
845       if (!(strcmp(filename, "default")) && fOutputEventHandler) {
846          isManagedByHandler = kTRUE;
847          filename = fOutputEventHandler->GetOutputFileName();
848       }
849       if (cont->IsSpecialOutput() || inGrid) {
850          if (strlen(fSpecialOutputLocation.Data())) continue;
851          // Copy merged file from PROOF scratch space. 
852          // In case of grid the files are already in the current directory.
853          if (!inGrid) {
854             if (isManagedByHandler && fExtraFiles.Length()) {
855                // Copy extra registered dAOD files.
856                TObjArray *arr = fExtraFiles.Tokenize(" ");
857                TObjString *os;
858                TIter nextfilename(arr);
859                while ((os=(TObjString*)nextfilename())) GetFileFromWrapper(os->GetString(), source);
860                delete arr;
861             }
862             if (!GetFileFromWrapper(filename, source)) continue;
863          }   
864          // Normally we should connect data from the copied file to the
865          // corresponding output container, but it is not obvious how to do this
866          // automatically if several objects in file...
867          TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
868          if (!f) f = TFile::Open(filename, "READ");
869          if (!f) {
870             Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
871             continue;
872          }   
873          f->cd();
874          TObject *obj = 0;
875          // Cd to the directory pointed by the container
876          TString folder = cont->GetFolderName();
877          if (!folder.IsNull()) f->cd(folder);
878          // Try to fetch first an object having the container name.
879          obj = gDirectory->Get(cont->GetName());
880          if (!obj) {
881             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", 
882                     cont->GetType()->GetName(), cont->GetName(), filename, cont->GetFolderName(), cont->GetName());
883             continue;
884          }  
885          wrap = new AliAnalysisDataWrapper(obj);
886          wrap->SetDeleteData(kFALSE);
887       }   
888       if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
889       if (!wrap) {
890          Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
891          continue;
892       }
893       icont++;
894       if (fDebug > 1) {
895          printf("   Importing data for container %s\n", cont->GetName());
896          if (strlen(filename)) printf("    -> file %s\n", filename);
897          else printf("\n");
898       }   
899       cont->ImportData(wrap);
900    }
901    if (cdir) cdir->cd();
902    fIOTimer->Stop();
903    fIOTime += fIOTimer->RealTime();
904    if (fDebug > 1) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
905 }
906
907 //______________________________________________________________________________
908 void AliAnalysisManager::UnpackOutput(TList *source)
909 {
910   // Called by AliAnalysisSelector::Terminate only on the client.
911    fIOTimer->Start(kTRUE);
912    if (fDebug > 1) printf("->AliAnalysisManager::UnpackOutput()\n");
913    if (!source) {
914       Error("UnpackOutput", "No target. Exiting.");
915       return;
916    }
917    if (fDebug > 1) printf("   Source list contains %d containers\n", source->GetSize());
918
919    if (fMode == kProofAnalysis) ImportWrappers(source);
920
921    TIter next(fOutputs);
922    AliAnalysisDataContainer *output;
923    while ((output=(AliAnalysisDataContainer*)next())) {
924       if (!output->GetData()) continue;
925       // Check if there are client tasks that run post event loop
926       if (output->HasConsumers()) {
927          // Disable event loop semaphore
928          output->SetPostEventLoop(kTRUE);
929          TObjArray *list = output->GetConsumers();
930          Int_t ncons = list->GetEntriesFast();
931          for (Int_t i=0; i<ncons; i++) {
932             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
933             task->CheckNotify(kTRUE);
934             // If task is active, execute it
935             if (task->IsPostEventLoop() && task->IsActive()) {
936                if (fDebug > 1) printf("== Executing post event loop task %s\n", task->GetName());
937                task->ExecuteTask();
938             }   
939          }
940       }   
941    }
942    fIOTimer->Stop();
943    fIOTime += fIOTimer->RealTime();
944    if (fDebug > 1) printf("<-AliAnalysisManager::UnpackOutput()\n");
945 }
946
947 //______________________________________________________________________________
948 void AliAnalysisManager::Terminate()
949 {
950   // The Terminate() function is the last function to be called during
951   // a query. It always runs on the client, it can be used to present
952   // the results graphically.
953    if (fDebug > 1) printf("->AliAnalysisManager::Terminate()\n");
954    fInitTimer->Start(kTRUE);
955    TDirectory *cdir = gDirectory;
956    gROOT->cd();
957    AliAnalysisTask *task;
958    AliAnalysisDataContainer *output;
959    TIter next(fTasks);
960    TStopwatch timer;
961    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
962    // Call Terminate() for tasks
963    Int_t itask = 0;
964    while (!IsSkipTerminate() && (task=(AliAnalysisTask*)next())) {
965       // Save all the canvases produced by the Terminate
966       TString pictname = Form("%s_%s", task->GetName(), task->ClassName());
967       task->Terminate();
968       gROOT->cd();
969       if (getsysInfo) 
970          AliSysInfo::AddStamp(Form("%s_TERMINATE",task->ClassName()),0, itask, 2);
971       itask++;   
972       if (TObject::TestBit(kSaveCanvases)) {
973          if (!gROOT->IsBatch()) {
974             if (fDebug>1) printf("Waiting 5 sec for %s::Terminate() to finish drawing ...\n", task->ClassName());
975             timer.Start(kTRUE);
976             while (timer.RealTime()<5) {
977                timer.Continue();
978                gSystem->ProcessEvents();
979             }
980          }
981          Int_t iend = gROOT->GetListOfCanvases()->GetEntries();
982          if (iend==0) continue;
983          TCanvas *canvas;
984          for (Int_t ipict=0; ipict<iend; ipict++) {
985             canvas = (TCanvas*)gROOT->GetListOfCanvases()->At(ipict);
986             if (!canvas) continue;         
987             canvas->SaveAs(Form("%s_%02d.gif", pictname.Data(),ipict));
988          } 
989          gROOT->GetListOfCanvases()->Delete(); 
990       }
991    }   
992    //
993    if (fInputEventHandler)   fInputEventHandler  ->TerminateIO();
994    if (fOutputEventHandler)  fOutputEventHandler ->TerminateIO();
995    if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
996    gROOT->cd();
997    TObjArray *allOutputs = new TObjArray();
998    Int_t icont;
999    for (icont=0; icont<fOutputs->GetEntriesFast(); icont++) allOutputs->Add(fOutputs->At(icont));
1000    if (!IsSkipTerminate())
1001       for (icont=0; icont<fParamCont->GetEntriesFast(); icont++) allOutputs->Add(fParamCont->At(icont));
1002    TIter next1(allOutputs);
1003    TString handlerFile = "";
1004    TString extraOutputs = "";
1005    if (fOutputEventHandler) {
1006       handlerFile = fOutputEventHandler->GetOutputFileName();
1007       extraOutputs = fOutputEventHandler->GetExtraOutputs();
1008    }
1009    icont = 0;
1010    TList filestmp;
1011    while ((output=(AliAnalysisDataContainer*)next1())) {
1012       // Special outputs or grid files have the files already closed and written.
1013       icont++;
1014       if (fMode == kGridAnalysis && icont<=fOutputs->GetEntriesFast()) continue;
1015       if (fMode == kProofAnalysis) {
1016         if (output->IsSpecialOutput() || output->IsRegisterDataset()) continue;
1017       }  
1018       const char *filename = output->GetFileName();
1019       TString openoption = "RECREATE";
1020       if (!(strcmp(filename, "default"))) continue;
1021       if (!strlen(filename)) continue;
1022       if (!output->GetData()) continue;
1023       TDirectory *opwd = gDirectory;
1024       TFile *file = output->GetFile();
1025       if (!file) file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1026       if (!file) {
1027               //if (handlerFile == filename && !gSystem->AccessPathName(filename)) openoption = "UPDATE";
1028          Bool_t firsttime = kTRUE;
1029          if (filestmp.FindObject(filename) || extraOutputs.Contains(filename)) {
1030             firsttime = kFALSE;
1031          } else {   
1032             filestmp.Add(new TNamed(filename,""));
1033          }   
1034          if (!gSystem->AccessPathName(filename) && !firsttime) openoption = "UPDATE";
1035               if (fDebug>1) printf("Opening file: %s  option=%s\n",filename, openoption.Data());
1036          file = new TFile(filename, openoption);
1037       } else {
1038          if (fDebug>1) printf("File <%s> already opened with option: <%s> \n", filename, file->GetOption());
1039          openoption = file->GetOption();
1040          if (openoption == "READ") {
1041             if (fDebug>1) printf("...reopening in UPDATE mode\n");
1042             file->ReOpen("UPDATE");            
1043          }
1044       }   
1045       if (file->IsZombie()) {
1046          Error("Terminate", "Cannot open output file %s", filename);
1047          continue;
1048       }   
1049       output->SetFile(file);
1050       file->cd();
1051       // Check for a folder request
1052       TString dir = output->GetFolderName();
1053       if (!dir.IsNull()) {
1054          if (!file->GetDirectory(dir)) file->mkdir(dir);
1055          file->cd(dir);
1056       }  
1057       if (fDebug > 1) printf("...writing container %s to file %s:%s\n", output->GetName(), file->GetName(), output->GetFolderName());
1058       if (output->GetData()->InheritsFrom(TCollection::Class())) {
1059       // If data is a collection, we set the name of the collection 
1060       // as the one of the container and we save as a single key.
1061          TCollection *coll = (TCollection*)output->GetData();
1062          coll->SetName(output->GetName());
1063          coll->Write(output->GetName(), TObject::kSingleKey);
1064       } else {
1065          if (output->GetData()->InheritsFrom(TTree::Class())) {
1066             TTree *tree = (TTree*)output->GetData();
1067             tree->SetDirectory(gDirectory);
1068             tree->AutoSave();
1069          } else {
1070             output->GetData()->Write();
1071          }   
1072       }      
1073       if (opwd) opwd->cd();
1074    }
1075    gROOT->cd();
1076    next1.Reset();
1077    TString copiedFiles;
1078    while ((output=(AliAnalysisDataContainer*)next1())) {
1079       // Close all files at output
1080       TDirectory *opwd = gDirectory;
1081       if (output->GetFile()) {
1082          // Clear file list to release object ownership to user.
1083 //         output->GetFile()->Clear();
1084          output->GetFile()->Close();
1085          // Copy merged outputs in alien if requested
1086          if (fSpecialOutputLocation.BeginsWith("alien://")) {
1087             if (copiedFiles.Contains(output->GetFile()->GetName())) {
1088                if (opwd) opwd->cd();
1089                output->SetFile(NULL);
1090                continue;
1091             } 
1092             Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data()); 
1093             gROOT->ProcessLine("if (!gGrid) TGrid::Connect(\"alien:\");");
1094             TFile::Cp(output->GetFile()->GetName(), 
1095                       Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
1096             copiedFiles += output->GetFile()->GetName();
1097          }             
1098          output->SetFile(NULL);
1099       }   
1100       if (opwd) opwd->cd();
1101    }   
1102    delete allOutputs;
1103    //Write statistics information on the client
1104    if (fStatistics) WriteStatisticsMsg(fNcalls);
1105    if (getsysInfo) {
1106       TDirectory *crtdir = gDirectory;
1107       TFile f("syswatch.root", "RECREATE");
1108       TH1 *hist;
1109       TString cut;
1110       if (!f.IsZombie()) {
1111          TTree *tree = AliSysInfo::MakeTree("syswatch.log");
1112          tree->SetName("syswatch");
1113          tree->SetMarkerStyle(kCircle);
1114          tree->SetMarkerColor(kBlue);
1115          tree->SetMarkerSize(0.5);
1116          if (!gROOT->IsBatch()) {
1117             tree->SetAlias("event", "id0");
1118             tree->SetAlias("task",  "id1");
1119             tree->SetAlias("stage", "id2");
1120             // Already defined aliases
1121             // tree->SetAlias("deltaT","stampSec-stampOldSec");
1122             // tree->SetAlias("T","stampSec-first");
1123             // tree->SetAlias("deltaVM","(pI.fMemVirtual-pIOld.fMemVirtual)");
1124             // tree->SetAlias("VM","pI.fMemVirtual");
1125             TCanvas *canvas = new TCanvas("SysInfo","SysInfo",10,10,1200,800);
1126             Int_t npads = 1 /*COO plot for all tasks*/ +
1127                           fTopTasks->GetEntries() /*Exec plot per task*/ +
1128                           1 /*Terminate plot for all tasks*/ +
1129                           1; /*vm plot*/
1130                           
1131             Int_t iopt = (Int_t)TMath::Sqrt((Double_t)npads);
1132             if (npads<iopt*(iopt+1))
1133                canvas->Divide(iopt, iopt+1, 0.01, 0.01);
1134             else
1135                canvas->Divide(iopt+1, iopt+1, 0.01, 0.01);
1136             Int_t ipad = 1;
1137             // draw the plot of deltaVM for Exec for each task
1138             for (itask=0; itask<fTopTasks->GetEntriesFast(); itask++) {
1139                task = (AliAnalysisTask*)fTopTasks->At(itask);
1140                canvas->cd(ipad++);
1141                cut = Form("task==%d && stage==1", itask);
1142                tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1143                hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");            
1144                if (hist) {
1145                   hist->SetTitle(Form("%s: Exec dVM[MB]/event", task->GetName()));
1146                   hist->GetYaxis()->SetTitle("deltaVM [MB]");
1147                }   
1148             }
1149             // Draw the plot of deltaVM for CreateOutputObjects for all tasks
1150             canvas->cd(ipad++);
1151             tree->SetMarkerStyle(kFullTriangleUp);
1152             tree->SetMarkerColor(kRed);
1153             tree->SetMarkerSize(0.8);
1154             cut = "task>=0 && task<1000 && stage==0";
1155             tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1156             hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");            
1157             if (hist) {
1158                hist->SetTitle("Memory in CreateOutputObjects()");
1159                hist->GetYaxis()->SetTitle("deltaVM [MB]");
1160                hist->GetXaxis()->SetTitle("task");
1161             }   
1162             // draw the plot of deltaVM for Terminate for all tasks
1163             canvas->cd(ipad++);
1164             tree->SetMarkerStyle(kOpenSquare);
1165             tree->SetMarkerColor(kMagenta);
1166             cut = "task>=0 && task<1000 && stage==2";
1167             tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1168             hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1169             if (hist) {
1170                hist->SetTitle("Memory in Terminate()");
1171                hist->GetYaxis()->SetTitle("deltaVM [MB]");
1172                hist->GetXaxis()->SetTitle("task");
1173             }   
1174             // Full VM profile
1175             canvas->cd(ipad++);
1176             tree->SetMarkerStyle(kFullCircle);
1177             tree->SetMarkerColor(kGreen);
1178             cut = Form("task==%d && stage==1",fTopTasks->GetEntriesFast()-1);            
1179             tree->Draw("VM:event",cut,"", 1234567890, 0);
1180             hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1181             if (hist) {
1182                hist->SetTitle("Virtual memory");
1183                hist->GetYaxis()->SetTitle("VM [MB]");
1184             }
1185             canvas->Modified();   
1186          }   
1187          tree->SetMarkerStyle(kCircle);
1188          tree->SetMarkerColor(kBlue);
1189          tree->SetMarkerSize(0.5);
1190          tree->Write();
1191          f.Close();
1192          delete tree;
1193       }
1194       if (crtdir) crtdir->cd();
1195    }
1196    // Validate the output files
1197    if (ValidateOutputFiles() && fIsRemote && fMode!=kProofAnalysis) {
1198       ofstream out;
1199       out.open("outputs_valid", ios::out);
1200       out.close();
1201    }
1202    if (cdir) cdir->cd();  
1203    fInitTimer->Stop();
1204    if (fDebug || IsCollectThroughput()) {
1205       printf("=Analysis %s= Terminate time:  %g[sec]\n", GetName(), fInitTimer->RealTime());
1206    }
1207    if (fDebug > 1) printf("<-AliAnalysisManager::Terminate()\n");
1208 }
1209 //______________________________________________________________________________
1210 void AliAnalysisManager::ProfileTask(Int_t itop, const char *option) const
1211 {
1212 // Profiles the task having the itop index in the list of top (first level) tasks.
1213    AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->At(itop);
1214    if (!task) {
1215       Error("ProfileTask", "There are only %d top tasks in the manager", fTopTasks->GetEntries());
1216       return;
1217    }
1218    ProfileTask(task->GetName(), option);
1219 }      
1220
1221 //______________________________________________________________________________
1222 void AliAnalysisManager::ProfileTask(const char *name, const char */*option*/) const
1223 {
1224 // Profile a managed task after the execution of the analysis in case NSysInfo
1225 // was used.
1226    if (gSystem->AccessPathName("syswatch.root")) {
1227       Error("ProfileTask", "No file syswatch.root found in the current directory");
1228       return;
1229    }
1230    if (gROOT->IsBatch()) return;
1231    AliAnalysisTask *task = (AliAnalysisTask*)fTopTasks->FindObject(name);
1232    if (!task) {
1233       Error("ProfileTask", "No top task named %s known by the manager.", name);
1234       return;
1235    }
1236    Int_t itop = fTopTasks->IndexOf(task);
1237    Int_t itask = fTasks->IndexOf(task);
1238    // Create canvas with 2 pads: first draw COO + Terminate, second Exec
1239    TDirectory *cdir = gDirectory;
1240    TFile f("syswatch.root");
1241    TTree *tree = (TTree*)f.Get("syswatch");
1242    if (!tree) {
1243       Error("ProfileTask", "No tree named <syswatch> found in file syswatch.root");
1244       return;
1245    }   
1246    if (fDebug > 1) printf("=== Profiling task %s (class %s)\n", name, task->ClassName());
1247    TCanvas *canvas = new TCanvas(Form("profile_%d",itop),Form("Profile of task %s (class %s)",name,task->ClassName()),10,10,800,600);
1248    canvas->Divide(2, 2, 0.01, 0.01);
1249    Int_t ipad = 1;
1250    TString cut;
1251    TH1 *hist;
1252    // VM profile for COO and Terminate methods
1253    canvas->cd(ipad++);
1254    cut = Form("task==%d && (stage==0 || stage==2)",itask);
1255    tree->Draw("deltaVM:sname",cut,"", 1234567890, 0);
1256    hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1257    if (hist) {
1258       hist->SetTitle("Alocated VM[MB] for COO and Terminate");
1259       hist->GetYaxis()->SetTitle("deltaVM [MB]");
1260       hist->GetXaxis()->SetTitle("method");
1261    }   
1262    // CPU profile per event
1263    canvas->cd(ipad++);
1264    cut = Form("task==%d && stage==1",itop);
1265    tree->Draw("deltaT:event",cut,"", 1234567890, 0);
1266    hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1267    if (hist) {
1268       hist->SetTitle("Execution time per event");
1269       hist->GetYaxis()->SetTitle("CPU/event [s]");
1270    }   
1271    // VM profile for Exec
1272    canvas->cd(ipad++);
1273    cut = Form("task==%d && stage==1",itop);
1274    tree->Draw("deltaVM:event",cut,"", 1234567890, 0);
1275    hist = (TH1*)gPad->GetListOfPrimitives()->FindObject("htemp");
1276    if (hist) {
1277       hist->SetTitle("Alocated VM[MB] per event");
1278       hist->GetYaxis()->SetTitle("deltaVM [MB]");
1279    }   
1280    canvas->Modified();
1281    delete tree;
1282    f.Close();
1283    if (cdir) cdir->cd();
1284 }     
1285
1286 //______________________________________________________________________________
1287 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
1288 {
1289 // Adds a user task to the global list of tasks.
1290    if (fInitOK) {
1291       Error("AddTask", "Cannot add task %s since InitAnalysis was already called", task->GetName());
1292       return;
1293    }   
1294       
1295    if (fTasks->FindObject(task)) {
1296       Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
1297       return;
1298    }   
1299    task->SetActive(kFALSE);
1300    fTasks->Add(task);
1301 }  
1302
1303 //______________________________________________________________________________
1304 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
1305 {
1306 // Retreive task by name.
1307    if (!fTasks) return NULL;
1308    return (AliAnalysisTask*)fTasks->FindObject(name);
1309 }
1310
1311 //______________________________________________________________________________
1312 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
1313                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
1314 {
1315 // Create a data container of a certain type. Types can be:
1316 //   kExchangeContainer  = 0, used to exchange data between tasks
1317 //   kInputContainer   = 1, used to store input data
1318 //   kOutputContainer  = 2, used for writing result to a file
1319 // filename: composed by file#folder (e.g. results.root#INCLUSIVE) - will write
1320 // the output object to a folder inside the output file
1321    if (fContainers->FindObject(name)) {
1322       Error("CreateContainer","A container named %s already defined !",name);
1323       return NULL;
1324    }   
1325    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
1326    fContainers->Add(cont);
1327    switch (type) {
1328       case kInputContainer:
1329          fInputs->Add(cont);
1330          break;
1331       case kOutputContainer:
1332          fOutputs->Add(cont);
1333          if (filename && strlen(filename)) {
1334             cont->SetFileName(filename);
1335             cont->SetDataOwned(kFALSE);  // data owned by the file
1336          }   
1337          break;
1338       case kParamContainer:
1339          fParamCont->Add(cont);
1340          if (filename && strlen(filename)) {
1341             cont->SetFileName(filename);
1342             cont->SetDataOwned(kFALSE);  // data owned by the file
1343          }   
1344          break;
1345       case kExchangeContainer:
1346          break;   
1347    }
1348    return cont;
1349 }
1350          
1351 //______________________________________________________________________________
1352 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
1353                                         AliAnalysisDataContainer *cont)
1354 {
1355 // Connect input of an existing task to a data container.
1356    if (!task) {
1357       Error("ConnectInput", "Task pointer is NULL");
1358       return kFALSE;
1359    }   
1360    if (!fTasks->FindObject(task)) {
1361       AddTask(task);
1362       Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
1363    } 
1364    Bool_t connected = task->ConnectInput(islot, cont);
1365    return connected;
1366 }   
1367
1368 //______________________________________________________________________________
1369 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
1370                                         AliAnalysisDataContainer *cont)
1371 {
1372 // Connect output of an existing task to a data container.
1373    if (!task) {
1374       Error("ConnectOutput", "Task pointer is NULL");
1375       return kFALSE;
1376    }   
1377    if (!fTasks->FindObject(task)) {
1378       AddTask(task);
1379       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
1380    } 
1381    Bool_t connected = task->ConnectOutput(islot, cont);
1382    return connected;
1383 }   
1384                                
1385 //______________________________________________________________________________
1386 void AliAnalysisManager::CleanContainers()
1387 {
1388 // Clean data from all containers that have already finished all client tasks.
1389    TIter next(fContainers);
1390    AliAnalysisDataContainer *cont;
1391    while ((cont=(AliAnalysisDataContainer *)next())) {
1392       if (cont->IsOwnedData() && 
1393           cont->IsDataReady() && 
1394           cont->ClientsExecuted()) cont->DeleteData();
1395    }
1396 }
1397
1398 //______________________________________________________________________________
1399 Bool_t AliAnalysisManager::InitAnalysis()
1400 {
1401 // Initialization of analysis chain of tasks. Should be called after all tasks
1402 // and data containers are properly connected
1403    // Reset flag and remove valid_outputs file if exists
1404    if (fInitOK) return kTRUE;
1405    if (!gSystem->AccessPathName("outputs_valid"))
1406       gSystem->Unlink("outputs_valid");
1407    // Check for top tasks (depending only on input data containers)
1408    if (!fTasks->First()) {
1409       Error("InitAnalysis", "Analysis has no tasks !");
1410       return kFALSE;
1411    }   
1412    TIter next(fTasks);
1413    AliAnalysisTask *task;
1414    AliAnalysisDataContainer *cont;
1415    Int_t ntop = 0;
1416    Int_t nzombies = 0;
1417    Bool_t iszombie = kFALSE;
1418    Bool_t istop = kTRUE;
1419    Int_t i;
1420    while ((task=(AliAnalysisTask*)next())) {
1421       istop = kTRUE;
1422       iszombie = kFALSE;
1423       Int_t ninputs = task->GetNinputs();
1424       for (i=0; i<ninputs; i++) {
1425          cont = task->GetInputSlot(i)->GetContainer();
1426          if (!cont) {
1427             if (!iszombie) {
1428                task->SetZombie();
1429                fZombies->Add(task);
1430                nzombies++;
1431                iszombie = kTRUE;
1432             }   
1433             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
1434                   i, task->GetName()); 
1435          }
1436          if (iszombie) continue;
1437          // Check if cont is an input container
1438          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
1439          // Connect to parent task
1440       }
1441       if (istop) {
1442          ntop++;
1443          fTopTasks->Add(task);
1444       }
1445    }
1446    if (!ntop) {
1447       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
1448       return kFALSE;
1449    }                        
1450    // Check now if there are orphan tasks
1451    for (i=0; i<ntop; i++) {
1452       task = (AliAnalysisTask*)fTopTasks->At(i);
1453       task->SetUsed();
1454    }
1455    Int_t norphans = 0;
1456    next.Reset();
1457    while ((task=(AliAnalysisTask*)next())) {
1458       if (!task->IsUsed()) {
1459          norphans++;
1460          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
1461       }   
1462    }          
1463    // Check the task hierarchy (no parent task should depend on data provided
1464    // by a daughter task)
1465    for (i=0; i<ntop; i++) {
1466       task = (AliAnalysisTask*)fTopTasks->At(i);
1467       if (task->CheckCircularDeps()) {
1468          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
1469          PrintStatus("dep");
1470          return kFALSE;
1471       }   
1472    }
1473    // Check that all containers feeding post-event loop tasks are in the outputs list
1474    TIter nextcont(fContainers); // loop over all containers
1475    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
1476       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
1477          if (cont->HasConsumers()) {
1478          // Check if one of the consumers is post event loop
1479             TIter nextconsumer(cont->GetConsumers());
1480             while ((task=(AliAnalysisTask*)nextconsumer())) {
1481                if (task->IsPostEventLoop()) {
1482                   fOutputs->Add(cont);
1483                   break;
1484                }
1485             }
1486          }
1487       }
1488    }   
1489    // Check if all special output containers have a file name provided
1490    TIter nextout(fOutputs);
1491    while ((cont=(AliAnalysisDataContainer*)nextout())) {
1492       if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
1493          Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
1494          return kFALSE;
1495       }
1496    }
1497    // Initialize requested branch list if needed
1498    if (!fAutoBranchHandling) {
1499       next.Reset();
1500       while ((task=(AliAnalysisTask*)next())) {
1501          if (!task->HasBranches()) {
1502             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\"",
1503                   task->GetName(), task->ClassName());
1504             return kFALSE;
1505          }
1506          if (!fInputEventHandler || !strlen(fInputEventHandler->GetDataType())) {
1507             Error("InitAnalysis", "Manual branch loading requested but no input handler defined or handler does not define data type.");
1508             return kFALSE;
1509          }
1510          TString taskbranches;
1511          task->GetBranches(fInputEventHandler->GetDataType(), taskbranches);
1512          if (taskbranches.IsNull()) {
1513             Error("InitAnalysis", "Manual branch loading requested but task %s of type %s does not define branches of type %s:",
1514                   task->GetName(), task->ClassName(), fInputEventHandler->GetDataType());
1515             return kFALSE;      
1516          }
1517          AddBranches(taskbranches);
1518       }         
1519    }
1520    fInitOK = kTRUE;
1521    return kTRUE;
1522 }   
1523
1524 //______________________________________________________________________________
1525 void AliAnalysisManager::AddBranches(const char *branches)
1526 {
1527 // Add branches to the existing fRequestedBranches.
1528    TString br(branches);
1529    TObjArray *arr = br.Tokenize(",");
1530    TIter next(arr);
1531    TObject *obj;
1532    while ((obj=next())) {
1533       if (!fRequestedBranches.Contains(obj->GetName())) {
1534          if (!fRequestedBranches.IsNull()) fRequestedBranches += ",";
1535          fRequestedBranches += obj->GetName();
1536       }
1537    }
1538   delete arr;
1539 }   
1540
1541 //______________________________________________________________________________
1542 void AliAnalysisManager::CheckBranches(Bool_t load)
1543 {
1544 // The method checks the input branches to be loaded during the analysis.
1545    if (fAutoBranchHandling || fRequestedBranches.IsNull() || !fTree) return;   
1546    TObjArray *arr = fRequestedBranches.Tokenize(",");
1547    TIter next(arr);
1548    TObject *obj;
1549    while ((obj=next())) {
1550       TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(obj->GetName()));
1551       if (!br) {
1552          br = fTree->GetBranch(obj->GetName());
1553          if (!br) {
1554             Error("CheckBranches", "Could not find branch %s",obj->GetName());
1555             continue;
1556          }
1557       }   
1558       fTable.Add(br);
1559       if (load && br->GetReadEntry()!=GetCurrentEntry()) {
1560          br->GetEntry(GetCurrentEntry());
1561       }      
1562    }
1563   delete arr;
1564 }
1565
1566 //______________________________________________________________________________
1567 Bool_t AliAnalysisManager::CheckTasks() const
1568 {
1569 // Check consistency of tasks.
1570    Int_t ntasks = fTasks->GetEntries();
1571    if (!ntasks) {
1572       Error("CheckTasks", "No tasks connected to the manager. This may be due to forgetting to compile the task or to load their library.");
1573       return kFALSE;
1574    }
1575    // Get the pointer to AliAnalysisTaskSE::Class()
1576    TClass *badptr = (TClass*)gROOT->ProcessLine("AliAnalysisTaskSE::Class()");
1577    // Loop all tasks to check if their corresponding library was loaded
1578    TIter next(fTasks);
1579    TObject *obj;
1580    while ((obj=next())) {
1581       if (obj->IsA() == badptr) {
1582          Error("CheckTasks", "##################\n \
1583          Class for task %s NOT loaded. You probably forgot to load the library for this task (or compile it dynamically).\n###########################\n",obj->GetName());
1584          return kFALSE;
1585       }
1586    }
1587    return kTRUE;      
1588 }   
1589
1590 //______________________________________________________________________________
1591 void AliAnalysisManager::PrintStatus(Option_t *option) const
1592 {
1593 // Print task hierarchy.
1594    if (!fInitOK) {
1595       Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
1596       return;
1597    }   
1598    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1599    if (getsysInfo)
1600       Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
1601    TIter next(fTopTasks);
1602    AliAnalysisTask *task;
1603    while ((task=(AliAnalysisTask*)next()))
1604       task->PrintTask(option);
1605   
1606    if (!fAutoBranchHandling && !fRequestedBranches.IsNull()) 
1607       printf("Requested input branches:\n%s\n", fRequestedBranches.Data());
1608   
1609   TString sopt(option);
1610   sopt.ToUpper();
1611   
1612   if (sopt.Contains("ALL"))
1613   {
1614     if ( fOutputEventHandler )
1615     {
1616       cout << TString('_',78) << endl;
1617       cout << "OutputEventHandler:" << endl;
1618       fOutputEventHandler->Print("   ");
1619     }
1620   }
1621 }
1622
1623 //______________________________________________________________________________
1624 void AliAnalysisManager::ResetAnalysis()
1625 {
1626 // Reset all execution flags and clean containers.
1627    CleanContainers();
1628 }
1629
1630 //______________________________________________________________________________
1631 void AliAnalysisManager::RunLocalInit()
1632 {
1633 // Run LocalInit method for all tasks.
1634    TDirectory *cdir = gDirectory;
1635    if (IsTrainInitialized()) return;
1636    TIter nextTask(fTasks);
1637    AliAnalysisTask *task;
1638    while ((task=(AliAnalysisTask*)nextTask())) {
1639       gROOT->cd();
1640       task->LocalInit();
1641    }
1642    if (cdir) cdir->cd();
1643    TObject::SetBit(kTasksInitialized, kTRUE);
1644 }   
1645
1646 //______________________________________________________________________________
1647 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
1648 {
1649 // Start analysis having a grid handler.
1650    if (!fGridHandler) {
1651       Error("StartAnalysis", "Cannot start analysis providing just the analysis type without a grid handler.");
1652       Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1653       return -1;
1654    }
1655    TTree *tree = NULL;
1656    return StartAnalysis(type, tree, nentries, firstentry);
1657 }
1658
1659 //______________________________________________________________________________
1660 Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree, Long64_t nentries, Long64_t firstentry)
1661 {
1662 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
1663 // MIX. Process nentries starting from firstentry
1664    Long64_t retv = 0;
1665    // Backup current directory and make sure gDirectory points to gROOT
1666    TDirectory *cdir = gDirectory;
1667    gROOT->cd();
1668    if (!fInitOK) {
1669       Error("StartAnalysis","Analysis manager was not initialized !");
1670       if (cdir) cdir->cd();
1671       return -1;
1672    }
1673    if (!CheckTasks()) Fatal("StartAnalysis", "Not all needed libraries were loaded");
1674    if (fDebug > 1) {
1675       printf("StartAnalysis %s\n",GetName());
1676       AliLog::SetGlobalLogLevel(AliLog::kInfo);
1677    }
1678    fMaxEntries = nentries;
1679    fIsRemote = kFALSE;
1680    TString anaType = type;
1681    anaType.ToLower();
1682    fMode = kLocalAnalysis;
1683    if (anaType.Contains("file"))      fIsRemote = kTRUE;
1684    if (anaType.Contains("proof"))     fMode = kProofAnalysis;
1685    else if (anaType.Contains("grid")) fMode = kGridAnalysis;
1686    else if (anaType.Contains("mix"))  fMode = kMixingAnalysis;
1687
1688    if (fMode == kGridAnalysis) {
1689       fIsRemote = kTRUE;
1690       if (!anaType.Contains("terminate")) {
1691          if (!fGridHandler) {
1692             Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1693             Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1694             if (cdir) cdir->cd();
1695             return -1;
1696          }
1697          // Write analysis manager in the analysis file
1698          cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1699          // run local task configuration
1700          RunLocalInit();
1701          if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1702             Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1703             if (cdir) cdir->cd();
1704             return -1;
1705          }   
1706
1707          // Terminate grid analysis
1708          if (fSelector && fSelector->GetStatus() == -1) {if (cdir) cdir->cd(); return -1;}
1709          if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {if (cdir) cdir->cd(); return 0;}
1710          cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1711          if (!fGridHandler->MergeOutputs()) {
1712             // Return if outputs could not be merged or if it alien handler
1713             // was configured for offline mode or local testing.
1714             if (cdir) cdir->cd();
1715             return 0;
1716          }
1717       }   
1718       cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1719       if (cdir) cdir->cd();
1720       ImportWrappers(NULL);
1721       Terminate();
1722       if (cdir) cdir->cd();
1723       return 0;
1724    }
1725    TString line;
1726    SetEventLoop(kFALSE);
1727    // Enable event loop mode if a tree was provided
1728    if (tree || fGridHandler || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1729
1730    TChain *chain = 0;
1731    TString ttype = "TTree";
1732    if (tree && tree->IsA() == TChain::Class()) {
1733       chain = (TChain*)tree;
1734       if (!chain || !chain->GetListOfFiles()->First()) {
1735          Error("StartAnalysis", "Cannot process null or empty chain...");
1736          if (cdir) cdir->cd();
1737          return -1;
1738       }   
1739       ttype = "TChain";
1740    }   
1741
1742    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1743    if (getsysInfo) AliSysInfo::AddStamp("Start", 0);
1744    // Initialize locally all tasks (happens for all modes)
1745    TIter next(fTasks);
1746    AliAnalysisTask *task;
1747    RunLocalInit();
1748    
1749    switch (fMode) {
1750       case kLocalAnalysis:
1751          if (!tree && !fGridHandler) {
1752             TIter nextT(fTasks);
1753             // Call CreateOutputObjects for all tasks
1754             Int_t itask = 0;
1755             Bool_t dirStatus = TH1::AddDirectoryStatus();
1756             while ((task=(AliAnalysisTask*)nextT())) {
1757                TH1::AddDirectory(kFALSE);
1758                task->CreateOutputObjects();
1759                if (!task->CheckPostData()) {
1760                   Error("SlaveBegin","####### IMPORTANT! ####### \n\n\n\
1761                         Task %s (%s) did not call PostData() for all its outputs in (User)CreateOutputObjects()\n\n\
1762                         ########### FIX YOUR CODE, THIS WILL PRODUCE A FATAL ERROR IN FUTURE! ###########", task->GetName(), task->ClassName());
1763                }
1764                if (getsysInfo) AliSysInfo::AddStamp(Form("%s_CREATEOUTOBJ",task->ClassName()), 0, itask, 0);
1765                gROOT->cd();
1766                itask++;
1767             }   
1768             TH1::AddDirectory(dirStatus);
1769             if (IsExternalLoop()) {
1770                Info("StartAnalysis", "Initialization done. Event loop is controlled externally.\
1771                      \nSetData for top container, call ExecAnalysis in a loop and then Terminate manually");
1772                return 0;
1773             }
1774             ExecAnalysis();
1775             Terminate();
1776             return 0;
1777          } 
1778          fSelector = new AliAnalysisSelector(this);
1779          // Check if a plugin handler is used
1780          if (fGridHandler) {
1781             // Get the chain from the plugin
1782             TString dataType = "esdTree";
1783             if (fInputEventHandler) {
1784                dataType = fInputEventHandler->GetDataType();
1785                dataType.ToLower();
1786                dataType += "Tree";
1787             }   
1788             chain = fGridHandler->GetChainForTestMode(dataType);
1789             if (!chain) {
1790                Error("StartAnalysis", "No chain for test mode. Aborting.");
1791                return -1;
1792             }
1793             cout << "===== RUNNING LOCAL ANALYSIS" << GetName() << " ON CHAIN " << chain->GetName() << endl;
1794             retv = chain->Process(fSelector, "", nentries, firstentry);
1795             break;
1796          }
1797          // Run tree-based analysis via AliAnalysisSelector  
1798          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1799          retv = tree->Process(fSelector, "", nentries, firstentry);
1800          break;
1801       case kProofAnalysis:
1802          fIsRemote = kTRUE;
1803          // Check if the plugin is used
1804          if (fGridHandler) {
1805             return StartAnalysis(type, fGridHandler->GetProofDataSet(), nentries, firstentry);
1806          }
1807          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1808             Error("StartAnalysis", "No PROOF!!! Exiting.");
1809             if (cdir) cdir->cd();
1810             return -1;
1811          }   
1812          line = Form("gProof->AddInput((TObject*)%p);", this);
1813          gROOT->ProcessLine(line);
1814          if (chain) {
1815             chain->SetProof();
1816             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1817             retv = chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1818          } else {
1819             Error("StartAnalysis", "No chain!!! Exiting.");
1820             if (cdir) cdir->cd();
1821             return -1;
1822          }      
1823          break;
1824       case kGridAnalysis:
1825          fIsRemote = kTRUE;
1826          if (!anaType.Contains("terminate")) {
1827             if (!fGridHandler) {
1828                Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
1829                Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
1830                if (cdir) cdir->cd();
1831                return -1;
1832             }
1833             // Write analysis manager in the analysis file
1834             cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1835             // Start the analysis via the handler
1836             if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1837                Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
1838                if (cdir) cdir->cd();
1839                return -1;
1840             }   
1841
1842             // Terminate grid analysis
1843             if (fSelector && fSelector->GetStatus() == -1) {if (cdir) cdir->cd(); return -1;}
1844             if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) {if (cdir) cdir->cd(); return 0;}
1845             cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1846             if (!fGridHandler->MergeOutputs()) {
1847                // Return if outputs could not be merged or if it alien handler
1848                // was configured for offline mode or local testing.
1849                if (cdir) cdir->cd();
1850                return 0;
1851             }
1852          }   
1853          cout << "===== TERMINATING GRID ANALYSIS JOB: " << GetName() << endl;
1854          ImportWrappers(NULL);
1855          Terminate();
1856          if (cdir) cdir->cd();
1857          return 0;
1858       case kMixingAnalysis:   
1859          // Run event mixing analysis
1860          if (!fEventPool) {
1861             Error("StartAnalysis", "Cannot run event mixing without event pool");
1862             if (cdir) cdir->cd();
1863             return -1;
1864          }
1865          cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1866          fSelector = new AliAnalysisSelector(this);
1867          while ((chain=fEventPool->GetNextChain())) {
1868             next.Reset();
1869             // Call NotifyBinChange for all tasks
1870             while ((task=(AliAnalysisTask*)next()))
1871                if (!task->IsPostEventLoop()) task->NotifyBinChange();
1872             retv = chain->Process(fSelector);
1873             if (retv < 0) {
1874                Error("StartAnalysis", "Mixing analysis failed");
1875                if (cdir) cdir->cd();
1876                return retv;
1877             }   
1878          }
1879          PackOutput(fSelector->GetOutputList());
1880          Terminate();
1881    }
1882    if (cdir) cdir->cd();
1883    return retv;
1884 }   
1885
1886 //______________________________________________________________________________
1887 Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1888 {
1889 // Start analysis for this manager on a given dataset. Analysis task can be: 
1890 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1891    if (!fInitOK) {
1892       Error("StartAnalysis","Analysis manager was not initialized !");
1893       return -1;
1894    }
1895    fIsRemote = kTRUE;
1896    if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
1897    TString anaType = type;
1898    anaType.ToLower();
1899    if (!anaType.Contains("proof")) {
1900       Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1901       return -1;
1902    }   
1903    fMode = kProofAnalysis;
1904    TString line;
1905    SetEventLoop(kTRUE);
1906    // Set the dataset flag
1907    TObject::SetBit(kUseDataSet);
1908    fTree = 0;
1909    if (fGridHandler) {
1910       // Start proof analysis using the grid handler
1911       if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
1912          Error("StartAnalysis", "The grid plugin could not start PROOF analysis");
1913          return -1;
1914       }
1915       // Check if the plugin is in test mode
1916       if (fGridHandler->GetRunMode() == AliAnalysisGrid::kTest) {
1917          dataset = "test_collection";
1918       } else {
1919          dataset = fGridHandler->GetProofDataSet();
1920       }
1921    }   
1922
1923    if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1924       Error("StartAnalysis", "No PROOF!!! Exiting.");
1925       return -1;
1926    }   
1927
1928    // Initialize locally all tasks
1929    RunLocalInit();
1930       
1931    line = Form("gProof->AddInput((TObject*)%p);", this);
1932    gROOT->ProcessLine(line);
1933    Long_t retv;
1934    line = Form("gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1935                dataset, nentries, firstentry);
1936    cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1937    retv = (Long_t)gROOT->ProcessLine(line);
1938    return retv;
1939 }   
1940
1941 //______________________________________________________________________________
1942 TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
1943 {
1944 // Opens according the option the file specified by cont->GetFileName() and changes
1945 // current directory to cont->GetFolderName(). If the file was already opened, it
1946 // checks if the option UPDATE was preserved. File open via TProofOutputFile can
1947 // be optionally ignored.
1948   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
1949   TString filename = cont->GetFileName();
1950   TFile *f = NULL;
1951   if (filename.IsNull()) {
1952     ::Error("AliAnalysisManager::OpenFile", "No file name specified for container %s", cont->GetName());
1953     return NULL;
1954   }
1955   if (mgr->GetAnalysisType()==AliAnalysisManager::kProofAnalysis && cont->IsSpecialOutput()
1956       && !ignoreProof)
1957     f = mgr->OpenProofFile(cont,option);
1958   else {
1959     // Check first if the file is already opened
1960     f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
1961     if (f) {
1962       // Check if option "UPDATE" was preserved 
1963       TString opt(option);
1964       opt.ToUpper();
1965       if ((opt=="UPDATE") && (opt!=f->GetOption())) 
1966         ::Info("AliAnalysisManager::OpenFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
1967     } else {
1968       f = TFile::Open(filename, option);
1969     }    
1970   }   
1971   if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
1972     cont->SetFile(f);
1973     // Cd to file
1974     f->cd();
1975     // Check for a folder request
1976     TString dir = cont->GetFolderName(); 
1977     if (!dir.IsNull()) {
1978       if (!f->GetDirectory(dir)) f->mkdir(dir);
1979       f->cd(dir);
1980     }
1981     return f;
1982   }
1983   ::Fatal("AliAnalysisManager::OpenFile", "File %s could not be opened", filename.Data());
1984   cont->SetFile(NULL);
1985   return NULL;
1986 }    
1987  
1988 //______________________________________________________________________________
1989 TFile *AliAnalysisManager::OpenProofFile(AliAnalysisDataContainer *cont, const char *option, const char *extaod)
1990 {
1991 // Opens a special output file used in PROOF.
1992   TString line;
1993   TString filename = cont->GetFileName();
1994   if (cont == fCommonOutput) {
1995      if (fOutputEventHandler) {
1996         if (strlen(extaod)) filename = extaod;
1997         filename = fOutputEventHandler->GetOutputFileName();
1998      }   
1999      else Fatal("OpenProofFile","No output container. Exiting.");
2000   }   
2001   TFile *f = NULL;
2002   if (fMode!=kProofAnalysis || !fSelector) {
2003     Fatal("OpenProofFile","Cannot open PROOF file %s: no PROOF or selector",filename.Data());
2004     return NULL;
2005   } 
2006   if (fSpecialOutputLocation.Length()) {
2007     f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2008     if (f) {
2009       // Check if option "UPDATE" was preserved 
2010       TString opt(option);
2011       opt.ToUpper();
2012       if ((opt=="UPDATE") && (opt!=f->GetOption()))
2013         ::Info("OpenProofFile", "File %s already opened in %s mode!", cont->GetFileName(), f->GetOption());
2014     } else {
2015       f = new TFile(filename, option);
2016     }
2017     if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
2018       cont->SetFile(f);
2019       // Cd to file
2020       f->cd();
2021       // Check for a folder request
2022       TString dir = cont->GetFolderName(); 
2023       if (dir.Length()) {
2024         if (!f->GetDirectory(dir)) f->mkdir(dir);
2025         f->cd(dir);
2026       }      
2027       return f;
2028     }
2029     Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
2030     cont->SetFile(NULL);
2031     return NULL;       
2032   }
2033   // Check if there is already a proof output file in the output list
2034   TObject *pof = fSelector->GetOutputList()->FindObject(filename);
2035   if (pof) {
2036     // Get the actual file
2037     line = Form("((TProofOutputFile*)%p)->GetFileName();", pof);
2038     filename = (const char*)gROOT->ProcessLine(line);
2039     if (fDebug>1) {
2040       printf("File: %s already booked via TProofOutputFile\n", filename.Data());
2041     }  
2042     f = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2043     if (!f) {
2044        Fatal("OpenProofFile", "Proof output file found but no file opened for %s", filename.Data());
2045        return NULL;
2046     }   
2047     // Check if option "UPDATE" was preserved 
2048     TString opt(option);
2049     opt.ToUpper();
2050     if ((opt=="UPDATE") && (opt!=f->GetOption())) 
2051       Fatal("OpenProofFile", "File %s already opened, but not in UPDATE mode!", cont->GetFileName());
2052   } else {
2053     if (cont->IsRegisterDataset()) {
2054       TString dsetName = filename;
2055       dsetName.ReplaceAll(".root", cont->GetTitle());
2056       dsetName.ReplaceAll(":","_");
2057       if (fDebug>1) printf("Booking dataset: %s\n", dsetName.Data());
2058       line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\", \"DROV\", \"%s\");", filename.Data(), dsetName.Data());
2059     } else {
2060       if (fDebug>1) printf("Booking TProofOutputFile: %s to be merged\n", filename.Data());
2061       line = Form("TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename.Data());
2062     }
2063     if (fDebug > 1) printf("=== %s\n", line.Data());
2064     gROOT->ProcessLine(line);
2065     line = Form("pf->OpenFile(\"%s\");", option);
2066     gROOT->ProcessLine(line);
2067     f = gFile;
2068     if (fDebug > 1) {
2069       gROOT->ProcessLine("pf->Print()");
2070       printf(" == proof file name: %s", f->GetName());
2071     }   
2072     // Add to proof output list
2073     line = Form("((TList*)%p)->Add(pf);",fSelector->GetOutputList());
2074     if (fDebug > 1) printf("=== %s\n", line.Data());
2075     gROOT->ProcessLine(line);
2076   }
2077   if (f && !f->IsZombie() && !f->TestBit(TFile::kRecovered)) {
2078     cont->SetFile(f);
2079     // Cd to file
2080     f->cd();
2081     // Check for a folder request
2082     TString dir = cont->GetFolderName(); 
2083     if (!dir.IsNull()) {
2084       if (!f->GetDirectory(dir)) f->mkdir(dir);
2085       f->cd(dir);
2086     }
2087     return f;
2088   }
2089   Fatal("OpenProofFile", "File %s could not be opened", cont->GetFileName());
2090   cont->SetFile(NULL);  
2091   return NULL;
2092 }   
2093
2094 //______________________________________________________________________________
2095 void AliAnalysisManager::ExecAnalysis(Option_t *option)
2096 {
2097 // Execute analysis.
2098    static Long64_t nentries = 0;
2099    static TTree *lastTree = 0;
2100    static TStopwatch *timer = new TStopwatch();
2101    // Only the first call to Process will trigger a true Notify. Other Notify
2102    // coming before is ignored.
2103    if (!TObject::TestBit(AliAnalysisManager::kTrueNotify)) {
2104       TObject::SetBit(AliAnalysisManager::kTrueNotify);
2105       Notify();
2106    }   
2107    if (fDebug > 0) printf("MGR: Processing event #%d\n", fNcalls);
2108    else {
2109       if (fTree && (fTree != lastTree)) {
2110          nentries += fTree->GetEntries();
2111          lastTree = fTree;
2112       }   
2113       if (!fNcalls) timer->Start();
2114       if (!fIsRemote && TObject::TestBit(kUseProgressBar)) ProgressBar("Processing event", fNcalls, TMath::Min(fMaxEntries,nentries), timer, kFALSE);
2115    }
2116    fIOTimer->Start(kTRUE);
2117    gROOT->cd();
2118    TDirectory *cdir = gDirectory;
2119    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
2120    if (getsysInfo && ((fNcalls%fNSysInfo)==0)) AliSysInfo::AddStamp("Exec_start", (Int_t)fNcalls);
2121    if (!fInitOK) {
2122       Error("ExecAnalysis", "Analysis manager was not initialized !");
2123       if (cdir) cdir->cd();
2124       return;
2125    }
2126    fNcalls++;
2127    AliAnalysisTask *task;
2128    // Check if the top tree is active.
2129    if (fTree) {
2130       if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2131          AliSysInfo::AddStamp("Handlers_BeginEventGroup",fNcalls, 1002, 0);
2132       TIter next(fTasks);
2133    // De-activate all tasks
2134       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
2135       AliAnalysisDataContainer *cont = fCommonInput;
2136       if (!cont) cont = (AliAnalysisDataContainer*)fInputs->At(0);
2137       if (!cont) {
2138               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
2139          if (cdir) cdir->cd();
2140          return;
2141       }   
2142       cont->SetData(fTree); // This will notify all consumers
2143       Long64_t entry = fTree->GetTree()->GetReadEntry();      
2144 //
2145 //    Call BeginEvent() for optional input/output and MC services 
2146       if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
2147       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
2148       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
2149       gROOT->cd();
2150       if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2151          AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2152 //
2153 //    Execute the tasks
2154 //      TIter next1(cont->GetConsumers());
2155       fIOTimer->Stop();
2156       fIOTime += fIOTimer->RealTime();
2157       fCPUTimer->Start(kTRUE);
2158       TIter next1(fTopTasks);
2159       Int_t itask = 0;
2160       while ((task=(AliAnalysisTask*)next1())) {
2161          if (fDebug >1) {
2162             cout << "    Executing task " << task->GetName() << endl;
2163          }       
2164          task->ExecuteTask(option);
2165          gROOT->cd();
2166          if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2167             AliSysInfo::AddStamp(task->ClassName(), fNcalls, itask, 1);
2168          itask++;   
2169       }
2170       fCPUTimer->Stop();
2171       fCPUTime += fCPUTimer->RealTime();
2172       fIOTimer->Start(kTRUE);
2173 //
2174 //    Call FinishEvent() for optional output and MC services 
2175       if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
2176       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
2177       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2178       // Gather system information if requested
2179       if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2180          AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1001, 1);
2181       if (cdir) cdir->cd();   
2182       fIOTimer->Stop();
2183       fIOTime += fIOTimer->RealTime();
2184       return;
2185    }   
2186    // The event loop is not controlled by TSelector   
2187 //
2188 //  Call BeginEvent() for optional input/output and MC services 
2189    fIOTimer->Start(kTRUE);
2190    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(-1);
2191    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(-1);
2192    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
2193    fIOTimer->Stop();
2194    fIOTime += fIOTimer->RealTime();
2195    gROOT->cd();
2196    if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2197       AliSysInfo::AddStamp("Handlers_BeginEvent",fNcalls, 1000, 0);
2198    fCPUTimer->Start(kTRUE);
2199    TIter next2(fTopTasks);
2200    while ((task=(AliAnalysisTask*)next2())) {
2201       task->SetActive(kTRUE);
2202       if (fDebug > 1) {
2203          cout << "    Executing task " << task->GetName() << endl;
2204       }   
2205       task->ExecuteTask(option);
2206       gROOT->cd();
2207    }   
2208    fCPUTimer->Stop();
2209    fCPUTime += fCPUTimer->RealTime();
2210 //
2211 // Call FinishEvent() for optional output and MC services 
2212    fIOTimer->Start(kTRUE);
2213    if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
2214    if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
2215    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
2216    if (getsysInfo && ((fNcalls%fNSysInfo)==0)) 
2217       AliSysInfo::AddStamp("Handlers_FinishEvent",fNcalls, 1000, 1);
2218    if (cdir) cdir->cd();   
2219    fIOTimer->Stop();
2220    fIOTime += fIOTimer->RealTime();
2221 }
2222
2223 //______________________________________________________________________________
2224 Bool_t AliAnalysisManager::IsPipe(std::ostream &out)
2225 {
2226 // Check if the stdout is connected to a pipe (C.Holm)
2227   Bool_t ispipe = kFALSE;
2228   out.seekp(0, std::ios_base::cur);
2229   if (out.fail()) {
2230     out.clear();
2231     if (errno == ESPIPE) ispipe = kTRUE;
2232   }
2233   return ispipe;
2234 }
2235    
2236 //______________________________________________________________________________
2237 void AliAnalysisManager::SetInputEventHandler(AliVEventHandler* const handler)
2238 {
2239 // Set the input event handler and create a container for it.
2240    Changed();
2241    fInputEventHandler   = handler;
2242    if (!fCommonInput) fCommonInput = CreateContainer("cAUTO_INPUT", TChain::Class(), AliAnalysisManager::kInputContainer);
2243 }
2244
2245 //______________________________________________________________________________
2246 void AliAnalysisManager::SetOutputEventHandler(AliVEventHandler* const handler)
2247 {
2248 // Set the input event handler and create a container for it.
2249    Changed();
2250    fOutputEventHandler   = handler;
2251    if (!fCommonOutput) fCommonOutput = CreateContainer("cAUTO_OUTPUT", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
2252    fCommonOutput->SetSpecialOutput();
2253 }
2254
2255 //______________________________________________________________________________
2256 void AliAnalysisManager::SetDebugLevel(UInt_t level)
2257 {
2258 // Set verbosity of the analysis manager. If the progress bar is used, the call is ignored
2259    if (TObject::TestBit(kUseProgressBar)) {
2260       Info("SetDebugLevel","Ignored. Disable the progress bar first.");
2261       return;
2262    }
2263    fDebug = level;
2264 }
2265    
2266 //______________________________________________________________________________
2267 void AliAnalysisManager::SetUseProgressBar(Bool_t flag, Int_t freq)
2268 {
2269 // Enable a text mode progress bar. Resets debug level to 0.
2270    Info("SetUseProgressBar", "Progress bar enabled, updated every %d events.\n  ### NOTE: Debug level reset to 0 ###", freq);
2271    TObject::SetBit(kUseProgressBar,flag);
2272    fPBUpdateFreq = freq;
2273    fDebug = 0;
2274 }   
2275
2276 //______________________________________________________________________________
2277 void AliAnalysisManager::RegisterExtraFile(const char *fname)
2278 {
2279 // This method is used externally to register output files which are not
2280 // connected to any output container, so that the manager can properly register,
2281 // retrieve or merge them when running in distributed mode. The file names are
2282 // separated by blancs. The method has to be called in MyAnalysisTask::LocalInit().
2283    if (fExtraFiles.Contains(fname)) return;
2284    if (fExtraFiles.Length()) fExtraFiles += " ";
2285    fExtraFiles += fname;
2286 }
2287
2288 //______________________________________________________________________________
2289 Bool_t AliAnalysisManager::GetFileFromWrapper(const char *filename, const TList *source)
2290 {
2291 // Copy a file from the location specified ina the wrapper with the same name from the source list.
2292    char fullPath[512];
2293    char chUrl[512];
2294    char tmp[1024];
2295    TObject *pof =  source->FindObject(filename);
2296    if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
2297       Error("GetFileFromWrapper", "TProofOutputFile object not found in output list for file %s", filename);
2298       return kFALSE;
2299    }
2300    gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", ((TProofOutputFile*)%p)->GetOutputFileName());", fullPath, pof));
2301    gROOT->ProcessLine(Form("sprintf((char*)%p, \"%%s\", gProof->GetUrl());",chUrl));
2302    TString clientUrl(chUrl);
2303    TString fullPath_str(fullPath);
2304    if (clientUrl.Contains("localhost")){
2305       TObjArray* array = fullPath_str.Tokenize ( "//" );
2306       TObjString *strobj = ( TObjString *)array->At(1);
2307       TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
2308       TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
2309       fullPath_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
2310       fullPath_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
2311       if (fDebug > 1) Info("GetFileFromWrapper","Using tunnel from %s to %s",fullPath_str.Data(),filename);
2312       delete arrayPort;
2313       delete array;
2314    }
2315    else if (clientUrl.Contains("__lite__")) { 
2316      // Special case for ProofLite environement - get file info and copy. 
2317      gROOT->ProcessLine(Form("sprintf((char*)%p,\"%%s\",((TProofOutputFile*)%p)->GetDir());", tmp, pof));
2318      fullPath_str = Form("%s/%s", tmp, fullPath);
2319    }
2320    if (fDebug > 1) 
2321      Info("GetFileFromWrapper","Copying file %s from PROOF scratch space to %s", fullPath_str.Data(),filename);
2322    Bool_t gotit = TFile::Cp(fullPath_str.Data(), filename); 
2323    if (!gotit)
2324       Error("GetFileFromWrapper", "Could not get file %s from proof scratch space", filename);
2325    return gotit;
2326 }
2327
2328 //______________________________________________________________________________
2329 void AliAnalysisManager::GetAnalysisTypeString(TString &type) const
2330 {
2331 // Fill analysis type in the provided string.
2332    switch (fMode) {
2333       case kLocalAnalysis:
2334          type = "local";
2335          return;
2336       case kProofAnalysis:
2337          type = "proof";
2338          return;
2339       case kGridAnalysis:
2340          type = "grid";
2341          return;
2342       case kMixingAnalysis:
2343          type = "mix";
2344    }
2345 }
2346
2347 //______________________________________________________________________________
2348 Bool_t AliAnalysisManager::ValidateOutputFiles() const
2349 {
2350 // Validate all output files.
2351    TIter next(fOutputs);
2352    AliAnalysisDataContainer *output;
2353    TDirectory *cdir = gDirectory;
2354    TString openedFiles;
2355    while ((output=(AliAnalysisDataContainer*)next())) {
2356       if (output->IsRegisterDataset()) continue;
2357       TString filename = output->GetFileName();
2358       if (filename == "default") {
2359          if (!fOutputEventHandler) continue;
2360          filename = fOutputEventHandler->GetOutputFileName();
2361          // Main AOD may not be there
2362          if (gSystem->AccessPathName(filename)) continue;
2363       }
2364       // Check if the file is closed
2365       if (openedFiles.Contains(filename)) continue;;
2366       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
2367       if (file) {
2368          Warning("ValidateOutputs", "File %s was not closed. Closing.", filename.Data());
2369          // Clear file list to release object ownership to user.
2370 //         file->Clear();
2371          file->Close();
2372       }
2373       file = TFile::Open(filename);
2374       if (!file || file->IsZombie() || file->TestBit(TFile::kRecovered)) {
2375          Error("ValidateOutputs", "Output file <%s> was not created or invalid", filename.Data());
2376          if (cdir) cdir->cd();
2377          return kFALSE;
2378       }
2379       file->Close();
2380       openedFiles += filename;
2381       openedFiles += " ";
2382    }
2383    if (cdir) cdir->cd();
2384    return kTRUE;
2385 }   
2386
2387 //______________________________________________________________________________
2388 void AliAnalysisManager::ProgressBar(const char *opname, Long64_t current, Long64_t size, TStopwatch * const watch, Bool_t last, Bool_t refresh)
2389 {
2390 // Implements a nice text mode progress bar.
2391    static Long64_t icount = 0;
2392    static TString oname;
2393    static TString nname;
2394    static Long64_t ocurrent = 0;
2395    static Long64_t osize = 0;
2396    static Int_t oseconds = 0;
2397    static TStopwatch *owatch = 0;
2398    static Bool_t oneoftwo = kFALSE;
2399    static Int_t nrefresh = 0;
2400    static Int_t nchecks = 0;
2401    static char lastChar = 0;
2402    const char symbol[4] = {'-','\\','|','/'}; 
2403    
2404    if (!lastChar) lastChar = (IsPipe(std::cerr))?'\r':'\n';
2405    if (!refresh) {
2406       nrefresh = 0;
2407       if (!size) return;
2408       owatch = watch;
2409       oname = opname;
2410       ocurrent = TMath::Abs(current);
2411       osize = TMath::Abs(size);
2412       if (ocurrent > osize) ocurrent=osize;
2413    } else {
2414       nrefresh++;
2415       if (!osize) return;
2416    }     
2417    if ((current % fPBUpdateFreq) != 0) return;
2418    icount++;
2419    char progress[11] = "          ";
2420    Int_t ichar = icount%4;
2421    Double_t time = 0.;
2422    Int_t hours = 0;
2423    Int_t minutes = 0;
2424    Int_t seconds = 0;
2425    if (owatch && !last) {
2426       owatch->Stop();
2427       time = owatch->RealTime();
2428       seconds   = int(time) % 60;
2429       minutes   = (int(time) / 60) % 60;
2430       hours     = (int(time) / 60 / 60);
2431       if (refresh)  {
2432          if (oseconds==seconds) {
2433             owatch->Continue();
2434             return;
2435          }
2436          oneoftwo = !oneoftwo;   
2437       }
2438       oseconds = seconds;   
2439    }
2440    if (refresh && oneoftwo) {
2441       nname = oname;
2442       if (nchecks <= 0) nchecks = nrefresh+1;
2443       Int_t pctdone = (Int_t)(100.*nrefresh/nchecks);
2444       oname = Form("     == %d%% ==", pctdone);
2445    }         
2446    Double_t percent = 100.0*ocurrent/osize;
2447    Int_t nchar = Int_t(percent/10);
2448    if (nchar>10) nchar=10;
2449    Int_t i;
2450    for (i=0; i<nchar; i++)  progress[i] = '=';
2451    progress[nchar] = symbol[ichar];
2452    for (i=nchar+1; i<10; i++) progress[i] = ' ';
2453    progress[10] = '\0';
2454    oname += "                    ";
2455    oname.Remove(20);
2456    if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
2457    else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
2458    else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
2459    if (time>0.) {
2460      Int_t full   = Int_t(ocurrent > 0 ? 
2461                           time * (float(osize)/ocurrent) + .5 : 
2462                           99*3600+59*60+59); 
2463      Int_t remain = Int_t(full - time);
2464      Int_t rsec   = remain % 60;
2465      Int_t rmin   = (remain / 60) % 60;
2466      Int_t rhour  = (remain / 60 / 60);
2467      fprintf(stderr, "[%6.2f %%]   TIME %.2d:%.2d:%.2d  ETA %.2d:%.2d:%.2d%c",
2468              percent, hours, minutes, seconds, rhour, rmin, rsec, lastChar);
2469    }
2470    else fprintf(stderr, "[%6.2f %%]%c", percent, lastChar);
2471    if (refresh && oneoftwo) oname = nname;
2472    if (owatch) owatch->Continue();
2473    if (last) {
2474       icount = 0;
2475       owatch = 0;
2476       ocurrent = 0;
2477       osize = 0;
2478       oseconds = 0;
2479       oneoftwo = kFALSE;
2480       nrefresh = 0;
2481       fprintf(stderr, "\n");
2482    }   
2483 }
2484
2485 //______________________________________________________________________________
2486 void AliAnalysisManager::DoLoadBranch(const char *name) 
2487 {
2488   // Get tree and load branch if needed.
2489   static Long64_t crtEntry = -100;
2490
2491   if (fAutoBranchHandling || !fTree)
2492     return;
2493
2494   TBranch *br = dynamic_cast<TBranch*>(fTable.FindObject(name));
2495   if (!br) {
2496     br = fTree->GetBranch(name);
2497     if (!br) {
2498       Error("DoLoadBranch", "Could not find branch %s",name);
2499       return;
2500     }
2501     fTable.Add(br);
2502   }
2503   if (br->GetReadEntry()==fCurrentEntry) return;
2504   Long64_t readbytes = br->GetEntry(GetCurrentEntry());
2505   if (readbytes<0) {
2506     Error("DoLoadBranch", "Could not load entry %lld from branch %s",GetCurrentEntry(), name);
2507     if (crtEntry != fCurrentEntry) {
2508       CountEvent(1,0,1,0);
2509       crtEntry = fCurrentEntry;
2510     }  
2511   } else {
2512     if (crtEntry != fCurrentEntry) {
2513       CountEvent(1,1,0,0);
2514       crtEntry = fCurrentEntry;
2515     }
2516   }
2517 }
2518
2519 //______________________________________________________________________________
2520 void AliAnalysisManager::AddStatisticsTask(UInt_t offlineMask)
2521 {
2522 // Add the statistics task to the manager.
2523   if (fStatistics) {
2524      Info("AddStatisticsTask", "Already added");
2525      return;
2526   }
2527   TString line = Form("AliAnalysisTaskStat::AddToManager(%u);", offlineMask);
2528   gROOT->ProcessLine(line);
2529 }  
2530
2531 //______________________________________________________________________________
2532 void AliAnalysisManager::CountEvent(Int_t ninput, Int_t nprocessed, Int_t nfailed, Int_t naccepted)
2533 {
2534 // Bookkeep current event;
2535    if (!fStatistics) return;
2536    fStatistics->AddInput(ninput);
2537    fStatistics->AddProcessed(nprocessed);
2538    fStatistics->AddFailed(nfailed);
2539    fStatistics->AddAccepted(naccepted);
2540 }   
2541
2542 //______________________________________________________________________________
2543 void AliAnalysisManager::AddStatisticsMsg(const char *line)
2544 {
2545 // Add a line in the statistics message. If available, the statistics message is written
2546 // at the end of the SlaveTerminate phase on workers AND at the end of Terminate
2547 // on the client.
2548    if (!strlen(line)) return;
2549    if (!fStatisticsMsg.IsNull()) fStatisticsMsg += "\n";
2550    fStatisticsMsg += line;
2551 }
2552
2553 //______________________________________________________________________________
2554 void AliAnalysisManager::WriteStatisticsMsg(Int_t)
2555 {
2556 // If fStatistics is present, write the file in the format ninput_nprocessed_nfailed_naccepted.stat
2557    static Bool_t done = kFALSE;
2558    if (done) return;
2559    done = kTRUE;
2560    if (!fStatistics) return;
2561    ofstream out;
2562    AddStatisticsMsg(Form("Number of input events:        %lld",fStatistics->GetNinput()));
2563    AddStatisticsMsg(Form("Number of processed events:    %lld",fStatistics->GetNprocessed()));      
2564    AddStatisticsMsg(Form("Number of failed events (I/O): %lld",fStatistics->GetNfailed()));
2565    AddStatisticsMsg(Form("Number of accepted events for mask %s: %lld", AliAnalysisStatistics::GetMaskAsString(fStatistics->GetOfflineMask()), fStatistics->GetNaccepted()));
2566    out.open(Form("%lld_%lld_%lld_%lld.stat",fStatistics->GetNinput(),
2567                  fStatistics->GetNprocessed(),fStatistics->GetNfailed(),
2568                  fStatistics->GetNaccepted()), ios::out);      
2569    out << fStatisticsMsg << endl;
2570    out.close();
2571 }
2572
2573 //______________________________________________________________________________
2574 const char* AliAnalysisManager::GetOADBPath()
2575 {
2576 // returns the path of the OADB
2577 // this static function just depends on environment variables
2578
2579    static TString oadbPath;
2580
2581    if (gSystem->Getenv("OADB_PATH"))
2582       oadbPath = gSystem->Getenv("OADB_PATH");
2583    else if (gSystem->Getenv("ALICE_ROOT"))
2584       oadbPath.Form("%s/OADB", gSystem->Getenv("ALICE_ROOT"));
2585    else
2586       ::Fatal("AliAnalysisManager::GetOADBPath", "Cannot figure out AODB path. Define ALICE_ROOT or OADB_PATH!");
2587       
2588    return oadbPath;
2589 }
2590
2591 //______________________________________________________________________________
2592 void AliAnalysisManager::SetGlobalStr(const char *key, const char *value)
2593 {
2594 // Define a custom string variable mapped to a global unique name. The variable
2595 // can be then retrieved by a given analysis macro via GetGlobalStr(key).
2596    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2597    if (!mgr) {
2598       ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2599       return;
2600    }   
2601    Bool_t valid = kFALSE;
2602    TString existing = AliAnalysisManager::GetGlobalStr(key, valid);
2603    if (valid) {
2604       ::Error("AliAnalysisManager::SetGlobalStr", "Global %s = %s already defined.", key, existing.Data());
2605       return;
2606    }
2607    mgr->GetGlobals()->Add(new TObjString(key), new TObjString(value));
2608 }
2609
2610 //______________________________________________________________________________
2611 const char *AliAnalysisManager::GetGlobalStr(const char *key, Bool_t &valid)
2612 {
2613 // Static method to retrieve a global variable defined via SetGlobalStr.
2614    valid = kFALSE;
2615    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2616    if (!mgr) return 0;
2617    TObject *value = mgr->GetGlobals()->GetValue(key);
2618    if (!value) return 0;
2619    valid = kTRUE;
2620    return value->GetName();
2621 }
2622
2623 //______________________________________________________________________________
2624 void AliAnalysisManager::SetGlobalInt(const char *key, Int_t value)
2625 {
2626 // Define a custom integer variable mapped to a global unique name. The variable
2627 // can be then retrieved by a given analysis macro via GetGlobalInt(key).
2628    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2629    if (!mgr) {
2630       ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2631       return;
2632    }   
2633    Bool_t valid = kFALSE;
2634    Int_t existing = AliAnalysisManager::GetGlobalInt(key, valid);
2635    if (valid) {
2636       ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %i already defined.", key, existing);
2637       return;
2638    }
2639    mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%i",value)));
2640 }
2641
2642 //______________________________________________________________________________
2643 Int_t AliAnalysisManager::GetGlobalInt(const char *key, Bool_t &valid)
2644 {
2645 // Static method to retrieve a global variable defined via SetGlobalInt.
2646    valid = kFALSE;
2647    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2648    if (!mgr) return 0;
2649    TObject *value = mgr->GetGlobals()->GetValue(key);
2650    if (!value) return 0;
2651    valid = kTRUE;
2652    TString s = value->GetName();
2653    return s.Atoi();
2654 }
2655
2656 //______________________________________________________________________________
2657 void AliAnalysisManager::SetGlobalDbl(const char *key, Double_t value)
2658 {
2659 // Define a custom double precision variable mapped to a global unique name. The variable
2660 // can be then retrieved by a given analysis macro via GetGlobalInt(key).
2661    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2662    if (!mgr) {
2663       ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
2664       return;
2665    }   
2666    Bool_t valid = kFALSE;
2667    Double_t existing = AliAnalysisManager::GetGlobalDbl(key, valid);
2668    if (valid) {
2669       ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %g already defined.", key, existing);
2670       return;
2671    }
2672    mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%f.16",value)));
2673 }
2674
2675 //______________________________________________________________________________
2676 Double_t AliAnalysisManager::GetGlobalDbl(const char *key, Bool_t &valid)
2677 {
2678 // Static method to retrieve a global variable defined via SetGlobalDbl.
2679    valid = kFALSE;
2680    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
2681    if (!mgr) return 0;
2682    TObject *value = mgr->GetGlobals()->GetValue(key);
2683    if (!value) return 0;
2684    valid = kTRUE;
2685    TString s = value->GetName();
2686    return s.Atof();
2687 }
2688
2689 //______________________________________________________________________________
2690 void AliAnalysisManager::AddClassDebug(const char *className, Int_t debugLevel)
2691 {
2692 // Sets Class debug level
2693
2694    if (!fDebugOptions) {
2695       fDebugOptions = new TObjArray();
2696       fDebugOptions->SetOwner(kTRUE);
2697    }
2698
2699    // substracting DebugOffset, beacuse of AliLog::SetClassDebugLevel()
2700    debugLevel -= AliLog::kDebug-1;
2701
2702    TNamed *debugOpt = (TNamed*)fDebugOptions->FindObject(className);
2703    if (!debugOpt) {
2704      AliInfo(TString::Format("Adding debug level %d for class %s",debugLevel+AliLog::kDebug-1,className).Data());
2705      fDebugOptions->Add(new TNamed(className,TString::Format("%d",debugLevel).Data()));
2706    } else {
2707       TString oldDebugStr = debugOpt->GetTitle();
2708       Int_t oldDebug = oldDebugStr.Atoi();
2709       if (debugLevel > oldDebug) {
2710          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());
2711          debugOpt->SetTitle(TString::Format("%d",debugLevel).Data());
2712       } else {
2713          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());
2714       }
2715    }
2716 }
2717
2718 //______________________________________________________________________________
2719 void AliAnalysisManager::ApplyDebugOptions()
2720 {
2721 // Apply debug options
2722
2723    if (!fDebugOptions) return;
2724    
2725    TIter next(fDebugOptions);
2726    TNamed *debug;
2727    TString debugLevel;
2728    while ((debug=dynamic_cast<TNamed*>(next()))) {
2729       debugLevel = debug->GetTitle();
2730       AliInfo(TString::Format("ApplyDebugOptions : Class=%s debulLevel=%d",debug->GetName(),debugLevel.Atoi()+AliLog::kDebug-1).Data());
2731       AliLog::SetClassDebugLevel(debug->GetName(), debugLevel.Atoi());
2732    }
2733 }
2734
2735 //______________________________________________________________________________
2736 void AliAnalysisManager::Lock()
2737 {
2738 // Security lock. This is to detect NORMAL user errors and not really to
2739 // protect against intentional hacks.
2740    if (fLocked) return;
2741    fLocked = kTRUE;
2742    if (fInputEventHandler)  fInputEventHandler->Lock();
2743    if (fOutputEventHandler) fOutputEventHandler->Lock();
2744    if (fMCtruthEventHandler) fMCtruthEventHandler->Lock();
2745    Info("Lock","====== ANALYSIS MANAGER LOCKED ======");
2746 }
2747
2748 //______________________________________________________________________________
2749 void AliAnalysisManager::UnLock()
2750 {
2751 // Verbose unlocking. Hackers will be punished ;-) ... 
2752    if (!fLocked) return;
2753    fLocked = kFALSE;
2754    if (fInputEventHandler)  fInputEventHandler->UnLock();
2755    if (fOutputEventHandler) fOutputEventHandler->UnLock();
2756    if (fMCtruthEventHandler) fMCtruthEventHandler->UnLock();
2757    Info("UnLock", "====== ANALYSIS MANAGER UNLOCKED ======");
2758 }
2759
2760 //______________________________________________________________________________
2761 void AliAnalysisManager::Changed()
2762 {
2763 // All critical setters pass through the Changed method that throws an exception 
2764 // in case the lock was set.
2765    if (fLocked) Fatal("Changed","Critical setter called in locked mode");
2766 }