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