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