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