]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
- Set of changes needed for merging aod files in CAF. Some fixes also for event mixing.
[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 //   AliAnalysysManager - 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 <Riostream.h>
29
30 #include <TClass.h>
31 #include <TFile.h>
32 #include <TMethodCall.h>
33 #include <TChain.h>
34 #include <TSystem.h>
35 #include <TROOT.h>
36 #include <TCanvas.h>
37
38 #include "AliAnalysisSelector.h"
39 #include "AliAnalysisTask.h"
40 #include "AliAnalysisDataContainer.h"
41 #include "AliAnalysisDataSlot.h"
42 #include "AliVEventHandler.h"
43 #include "AliVEventPool.h"
44 #include "AliSysInfo.h"
45 #include "AliAnalysisManager.h"
46
47 ClassImp(AliAnalysisManager)
48
49 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
50
51 //______________________________________________________________________________
52 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
53                    :TNamed(name,title),
54                     fTree(NULL),
55                     fInputEventHandler(NULL),
56                     fOutputEventHandler(NULL),
57                     fMCtruthEventHandler(NULL),
58                     fEventPool(NULL),
59                     fCurrentEntry(-1),
60                     fNSysInfo(0),
61                     fMode(kLocalAnalysis),
62                     fInitOK(kFALSE),
63                     fDebug(0),
64                     fSpecialOutputLocation(""), 
65                     fTasks(NULL),
66                     fTopTasks(NULL),
67                     fZombies(NULL),
68                     fContainers(NULL),
69                     fInputs(NULL),
70                     fOutputs(NULL),
71                     fSelector(NULL)
72 {
73 // Default constructor.
74    fgAnalysisManager = this;
75    fTasks      = new TObjArray();
76    fTopTasks   = new TObjArray();
77    fZombies    = new TObjArray();
78    fContainers = new TObjArray();
79    fInputs     = new TObjArray();
80    fOutputs    = new TObjArray();
81    SetEventLoop(kTRUE);
82 }
83
84 //______________________________________________________________________________
85 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
86                    :TNamed(other),
87                     fTree(NULL),
88                     fInputEventHandler(NULL),
89                     fOutputEventHandler(NULL),
90                     fMCtruthEventHandler(NULL),
91                     fEventPool(NULL),
92                     fCurrentEntry(-1),
93                     fNSysInfo(0),
94                     fMode(other.fMode),
95                     fInitOK(other.fInitOK),
96                     fDebug(other.fDebug),
97                     fSpecialOutputLocation(""), 
98                     fTasks(NULL),
99                     fTopTasks(NULL),
100                     fZombies(NULL),
101                     fContainers(NULL),
102                     fInputs(NULL),
103                     fOutputs(NULL),
104                     fSelector(NULL)
105 {
106 // Copy constructor.
107    fTasks      = new TObjArray(*other.fTasks);
108    fTopTasks   = new TObjArray(*other.fTopTasks);
109    fZombies    = new TObjArray(*other.fZombies);
110    fContainers = new TObjArray(*other.fContainers);
111    fInputs     = new TObjArray(*other.fInputs);
112    fOutputs    = new TObjArray(*other.fOutputs);
113    fgAnalysisManager = this;
114 }
115    
116 //______________________________________________________________________________
117 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
118 {
119 // Assignment
120    if (&other != this) {
121       TNamed::operator=(other);
122       fInputEventHandler   = other.fInputEventHandler;
123       fOutputEventHandler  = other.fOutputEventHandler;
124       fMCtruthEventHandler = other.fMCtruthEventHandler;
125       fEventPool           = other.fEventPool;
126       fTree       = NULL;
127       fCurrentEntry = -1;
128       fNSysInfo   = other.fNSysInfo;
129       fMode       = other.fMode;
130       fInitOK     = other.fInitOK;
131       fDebug      = other.fDebug;
132       fTasks      = new TObjArray(*other.fTasks);
133       fTopTasks   = new TObjArray(*other.fTopTasks);
134       fZombies    = new TObjArray(*other.fZombies);
135       fContainers = new TObjArray(*other.fContainers);
136       fInputs     = new TObjArray(*other.fInputs);
137       fOutputs    = new TObjArray(*other.fOutputs);
138       fSelector   = NULL;
139       fgAnalysisManager = this;
140    }
141    return *this;
142 }
143
144 //______________________________________________________________________________
145 AliAnalysisManager::~AliAnalysisManager()
146 {
147 // Destructor.
148    if (fTasks) {fTasks->Delete(); delete fTasks;}
149    if (fTopTasks) delete fTopTasks;
150    if (fZombies) delete fZombies;
151    if (fContainers) {fContainers->Delete(); delete fContainers;}
152    if (fInputs) delete fInputs;
153    if (fOutputs) delete fOutputs;
154    if (fgAnalysisManager==this) fgAnalysisManager = NULL;
155 }
156
157 //______________________________________________________________________________
158 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
159 {
160 // Read one entry of the tree or a whole branch.
161    if (fDebug > 0) printf("== AliAnalysisManager::GetEntry(%lld)\n", entry);
162    fCurrentEntry = entry;
163    return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
164 }
165    
166 //______________________________________________________________________________
167 void AliAnalysisManager::Init(TTree *tree)
168 {
169   // The Init() function is called when the selector needs to initialize
170   // a new tree or chain. Typically here the branch addresses of the tree
171   // will be set. It is normaly not necessary to make changes to the
172   // generated code, but the routine can be extended by the user if needed.
173   // Init() will be called many times when running with PROOF.
174    if (!tree) return;
175    if (fDebug > 0) {
176       printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
177    }
178
179    // Call InitTree of EventHandler
180    if (fOutputEventHandler) {
181       if (fMode == kProofAnalysis) {
182          fOutputEventHandler->Init(0x0, "proof");
183       } else {
184          fOutputEventHandler->Init(0x0, "local");
185       }
186    }
187
188    if (fInputEventHandler) {
189       if (fMode == kProofAnalysis) {
190          fInputEventHandler->Init(tree, "proof");
191       } else {
192          fInputEventHandler->Init(tree, "local");
193       }
194    } else {
195       // If no input event handler we need to get the tree once
196       // for the chain
197       if(!tree->GetTree()) tree->LoadTree(0);
198    }
199    
200
201    if (fMCtruthEventHandler) {
202       if (fMode == kProofAnalysis) {
203          fMCtruthEventHandler->Init(0x0, "proof");
204       } else {
205          fMCtruthEventHandler->Init(0x0, "local");
206       }
207    }
208
209    if (!fInitOK) InitAnalysis();
210    if (!fInitOK) return;
211    fTree = tree;
212    AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
213    if (!top) {
214       Error("Init","No top input container !");
215       return;
216    }
217    top->SetData(tree);
218    if (fDebug > 0) {
219       printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
220    }
221 }
222
223 //______________________________________________________________________________
224 void AliAnalysisManager::SlaveBegin(TTree *tree)
225 {
226   // The SlaveBegin() function is called after the Begin() function.
227   // When running with PROOF SlaveBegin() is called on each slave server.
228   // The tree argument is deprecated (on PROOF 0 is passed).
229    if (fDebug > 0) printf("->AliAnalysisManager::SlaveBegin()\n");
230    static Bool_t isCalled = kFALSE;
231    TDirectory *curdir = gDirectory;
232    // Call SlaveBegin only once in case of mixing
233    if (isCalled && fMode==kMixingAnalysis) return;
234    // Call Init of EventHandler
235    if (fOutputEventHandler) {
236       if (fMode == kProofAnalysis) {
237          TIter nextout(fOutputs);
238          AliAnalysisDataContainer *c_aod;
239          while ((c_aod=(AliAnalysisDataContainer*)nextout())) if (!strcmp(c_aod->GetFileName(),"default")) break;
240          if (c_aod && c_aod->IsSpecialOutput()) {
241             // Merging via files
242             if (fDebug > 1) printf("   Initializing special output file %s...\n", fOutputEventHandler->GetOutputFileName());
243             OpenProofFile(fOutputEventHandler->GetOutputFileName(), "RECREATE");
244             c_aod->SetFile(gFile);
245             fOutputEventHandler->Init("proofspecial");
246          } else {
247             // Merging in memory
248             fOutputEventHandler->Init("proof");
249          }   
250       } else {
251          fOutputEventHandler->Init("local");
252       }
253    }
254
255    if (fInputEventHandler) {
256       fInputEventHandler->SetInputTree(tree);
257       if (fMode == kProofAnalysis) {
258          fInputEventHandler->Init("proof");
259       } else {
260          fInputEventHandler->Init("local");
261       }
262    }
263
264    if (fMCtruthEventHandler) {
265       if (fMode == kProofAnalysis) {
266          fMCtruthEventHandler->Init("proof");
267       } else {
268          fMCtruthEventHandler->Init("local");
269       }
270    }
271    if (curdir) curdir->cd();
272    
273    TIter next(fTasks);
274    AliAnalysisTask *task;
275    // Call CreateOutputObjects for all tasks
276    while ((task=(AliAnalysisTask*)next())) {
277       curdir = gDirectory;
278       task->CreateOutputObjects();
279       if (curdir) curdir->cd();
280    }
281    isCalled = kTRUE;
282
283    if (fDebug > 0) printf("<-AliAnalysisManager::SlaveBegin()\n");
284 }
285
286 //______________________________________________________________________________
287 Bool_t AliAnalysisManager::Notify()
288 {
289    // The Notify() function is called when a new file is opened. This
290    // can be either for a new TTree in a TChain or when when a new TTree
291    // is started when using PROOF. It is normaly not necessary to make changes
292    // to the generated code, but the routine can be extended by the
293    // user if needed. The return value is currently not used.
294    if (!fTree) {
295       Error("Notify","No current tree.");
296       return kFALSE;
297    }   
298    TFile *curfile = fTree->GetCurrentFile();
299    if (!curfile) {
300       Error("Notify","No current file");
301       return kFALSE;
302    }   
303    
304    if (fDebug > 0) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
305    TIter next(fTasks);
306    AliAnalysisTask *task;
307    // Call Notify for all tasks
308    while ((task=(AliAnalysisTask*)next())) 
309       task->Notify();
310         
311    // Call Notify of the event handlers
312    if (fInputEventHandler) {
313        fInputEventHandler->Notify(curfile->GetName());
314    }
315
316    if (fOutputEventHandler) {
317        fOutputEventHandler->Notify(curfile->GetName());
318    }
319
320    if (fMCtruthEventHandler) {
321        fMCtruthEventHandler->Notify(curfile->GetName());
322    }
323    if (fDebug > 0) printf("<-AliAnalysisManager::Notify()\n");
324    return kTRUE;
325 }    
326
327 //______________________________________________________________________________
328 Bool_t AliAnalysisManager::Process(Long64_t entry)
329 {
330   // The Process() function is called for each entry in the tree (or possibly
331   // keyed object in the case of PROOF) to be processed. The entry argument
332   // specifies which entry in the currently loaded tree is to be processed.
333   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
334   // to read either all or the required parts of the data. When processing
335   // keyed objects with PROOF, the object is already loaded and is available
336   // via the fObject pointer.
337   //
338   // This function should contain the "body" of the analysis. It can contain
339   // simple or elaborate selection criteria, run algorithms on the data
340   // of the event and typically fill histograms.
341
342   // WARNING when a selector is used with a TChain, you must use
343   //  the pointer to the current TTree to call GetEntry(entry).
344   //  The entry is always the local entry number in the current tree.
345   //  Assuming that fChain is the pointer to the TChain being processed,
346   //  use fChain->GetTree()->GetEntry(entry).
347    if (fDebug > 0) printf("->AliAnalysisManager::Process(%lld)\n", entry);
348
349    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
350    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
351    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
352    
353    GetEntry(entry);
354    ExecAnalysis();
355    if (fDebug > 0) printf("<-AliAnalysisManager::Process()\n");
356    return kTRUE;
357 }
358
359 //______________________________________________________________________________
360 void AliAnalysisManager::PackOutput(TList *target)
361 {
362   // Pack all output data containers in the output list. Called at SlaveTerminate
363   // stage in PROOF case for each slave.
364    if (fDebug > 0) printf("->AliAnalysisManager::PackOutput()\n");
365    if (!target) {
366       Error("PackOutput", "No target. Aborting.");
367       return;
368    }
369    if (fInputEventHandler)   fInputEventHandler  ->Terminate();
370    if (fOutputEventHandler)  fOutputEventHandler ->Terminate();
371    if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
372
373    // Call FinishTaskOutput() for each event loop task (not called for 
374    // post-event loop tasks - use Terminate() fo those)
375    TIter nexttask(fTasks);
376    AliAnalysisTask *task;
377    while ((task=(AliAnalysisTask*)nexttask())) {
378       if (!task->IsPostEventLoop()) {
379          if (fDebug > 0) printf("->FinishTaskOutput: task %s\n", task->GetName());
380          task->FinishTaskOutput();
381          if (fDebug > 0) printf("<-FinishTaskOutput: task %s\n", task->GetName());
382       }
383    }      
384    
385    if (fMode == kProofAnalysis) {
386       TIter next(fOutputs);
387       AliAnalysisDataContainer *output;
388       Bool_t isManagedByHandler = kFALSE;
389       while ((output=(AliAnalysisDataContainer*)next())) {
390          // Do not consider outputs of post event loop tasks
391          if (output->GetProducer()->IsPostEventLoop()) continue;
392          const char *filename = output->GetFileName();
393          if (!(strcmp(filename, "default")) && fOutputEventHandler) {
394             isManagedByHandler = kTRUE;
395             filename = fOutputEventHandler->GetOutputFileName();
396          }
397          // Check if data was posted to this container. If not, issue an error.
398          if (!output->GetData() && !isManagedByHandler) {
399             Error("PackOutput", "No data for output container %s. Forgot to PostData ?\n", output->GetName());
400             continue;
401          }   
402          if (!output->IsSpecialOutput()) {
403             // Normal outputs
404             if (strlen(filename) && !isManagedByHandler) {
405                // File resident outputs
406                TFile *file = output->GetFile();
407                // Backup current folder
408                TDirectory *opwd = gDirectory;
409                // Create file if not existing and register to container.
410                if (file) file->cd();
411                else      file = new TFile(filename, "RECREATE"); 
412                if (file->IsZombie()) {
413                   Fatal("PackOutput", "Could not recreate file %s\n", filename);
414                   return;
415                }   
416                output->SetFile(file);
417                // Clear file list to release object ownership to user.
418                file->Clear();
419                // Save data to file, then close.
420                if (output->GetData()->InheritsFrom(TCollection::Class())) {
421                   // If data is a collection, we set the name of the collection 
422                   // as the one of the container and we save as a single key.
423                   TCollection *coll = (TCollection*)output->GetData();
424                   coll->SetName(output->GetName());
425                   coll->Write(output->GetName(), TObject::kSingleKey);
426                } else {
427                   output->GetData()->Write();
428                }      
429                if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
430                if (fDebug > 2) {
431                   printf("   file %s listing content:\n", filename);
432                   file->ls();
433                }   
434                file->Close();
435                // Restore current directory
436                if (opwd) opwd->cd();
437             } else {
438                // Memory-resident outputs   
439                if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
440             }   
441             AliAnalysisDataWrapper *wrap = 0;
442             if (isManagedByHandler) {
443                wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
444                wrap->SetName(output->GetName());
445             }   
446             else                    wrap =output->ExportData();
447             // Output wrappers must delete data after merging (AG 13/11/07)
448             wrap->SetDeleteData(kTRUE);
449             target->Add(wrap);
450          } else {
451          // Special outputs
452             TDirectory *opwd = gDirectory;
453             TFile *file = output->GetFile();
454             if (isManagedByHandler) {
455                // Terminate IO for files managed by the output handler
456                if (file) file->Write();
457                fOutputEventHandler->TerminateIO();
458                continue;
459             }   
460             
461             if (!file) {
462                AliAnalysisTask *producer = output->GetProducer();
463                Error("PackOutput", 
464                      "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
465                      output->GetFileName(), output->GetName(), producer->ClassName());
466                continue;
467             }   
468             file->cd();
469             // Release object ownership to users after writing data to file
470             if (output->GetData()->InheritsFrom(TCollection::Class())) {
471                // If data is a collection, we set the name of the collection 
472                // as the one of the container and we save as a single key.
473                TCollection *coll = (TCollection*)output->GetData();
474                coll->SetName(output->GetName());
475                coll->Write(output->GetName(), TObject::kSingleKey);
476             } else {
477                output->GetData()->Write();
478             }      
479             file->Clear();
480             if (fDebug > 2) {
481                printf("   file %s listing content:\n", output->GetFileName());
482                file->ls();
483             }   
484             file->Close();
485             // Restore current directory
486             if (opwd) opwd->cd();
487             // Check if a special output location was provided or the output files have to be merged
488             if (strlen(fSpecialOutputLocation.Data())) {
489                TString remote = fSpecialOutputLocation;
490                remote += "/";
491                Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
492                remote += Form("%s_%d_", gSystem->HostName(), gid);
493                remote += output->GetFileName();
494                TFile::Cp(output->GetFileName(), remote.Data());
495             } else {
496             // No special location specified-> use TProofOutputFile as merging utility
497             // The file at this output slot must be opened in CreateOutputObjects
498                if (fDebug > 1) printf("   File %s to be merged...\n", output->GetFileName());
499             }
500          }      
501       }
502    } 
503    if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
504 }
505
506 //______________________________________________________________________________
507 void AliAnalysisManager::ImportWrappers(TList *source)
508 {
509 // Import data in output containers from wrappers coming in source.
510    if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
511    TIter next(fOutputs);
512    AliAnalysisDataContainer *cont;
513    AliAnalysisDataWrapper   *wrap;
514    Int_t icont = 0;
515    while ((cont=(AliAnalysisDataContainer*)next())) {
516       wrap = 0;
517       if (cont->GetProducer()->IsPostEventLoop()) continue;
518       const char *filename = cont->GetFileName();
519       Bool_t isManagedByHandler = kFALSE;
520       if (!(strcmp(filename, "default")) && fOutputEventHandler) {
521          isManagedByHandler = kTRUE;
522          filename = fOutputEventHandler->GetOutputFileName();
523       }
524       if (cont->IsSpecialOutput()) {
525          if (strlen(fSpecialOutputLocation.Data()) && !isManagedByHandler) continue;
526          // Copy merged file from PROOF scratch space
527          char full_path[512];
528          TObject *pof =  source->FindObject(filename);
529          if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
530             Error("ImportWrappers", "TProofOutputFile object not found in output list for container %s", cont->GetName());
531             continue;
532          }
533          gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", ((TProofOutputFile*)0x%lx)->GetOutputFileName();)", full_path, pof));
534          if (fDebug > 1) 
535             printf("   Copying file %s from PROOF scratch space\n", full_path);
536          Bool_t gotit = TFile::Cp(full_path, filename); 
537          if (!gotit) {
538             Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
539          }
540          // Normally we should connect data from the copied file to the
541          // corresponding output container, but it is not obvious how to do this
542          // automatically if several objects in file...
543          TFile *f = new TFile(filename, "READ");
544          TObject *obj = f->Get(cont->GetName());
545          if (!obj) {
546             Error("ImportWrappers", "Could not find object %s in file %s", cont->GetName(), filename);
547             continue;
548          }
549          wrap = new AliAnalysisDataWrapper(obj);
550          wrap->SetDeleteData(kFALSE);
551       }   
552       if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
553       if (!wrap) {
554          Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
555          continue;
556       }
557       icont++;
558       if (fDebug > 1) {
559          printf("   Importing data for container %s", cont->GetName());
560          if (strlen(filename)) printf("    -> file %s\n", cont->GetFileName());
561          else printf("\n");
562       }   
563       cont->ImportData(wrap);
564    }         
565    if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
566 }
567
568 //______________________________________________________________________________
569 void AliAnalysisManager::UnpackOutput(TList *source)
570 {
571   // Called by AliAnalysisSelector::Terminate only on the client.
572    if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
573    if (!source) {
574       Error("UnpackOutput", "No target. Aborting.");
575       return;
576    }
577    if (fDebug > 1) printf("   Source list contains %d containers\n", source->GetSize());
578
579    if (fMode == kProofAnalysis) ImportWrappers(source);
580
581    TIter next(fOutputs);
582    AliAnalysisDataContainer *output;
583    while ((output=(AliAnalysisDataContainer*)next())) {
584       if (!output->GetData()) continue;
585       // Check if there are client tasks that run post event loop
586       if (output->HasConsumers()) {
587          // Disable event loop semaphore
588          output->SetPostEventLoop(kTRUE);
589          TObjArray *list = output->GetConsumers();
590          Int_t ncons = list->GetEntriesFast();
591          for (Int_t i=0; i<ncons; i++) {
592             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
593             task->CheckNotify(kTRUE);
594             // If task is active, execute it
595             if (task->IsPostEventLoop() && task->IsActive()) {
596                if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
597                task->ExecuteTask();
598             }   
599          }
600       }   
601    }
602    if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
603 }
604
605 //______________________________________________________________________________
606 void AliAnalysisManager::Terminate()
607 {
608   // The Terminate() function is the last function to be called during
609   // a query. It always runs on the client, it can be used to present
610   // the results graphically.
611    if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
612    AliAnalysisTask *task;
613    TIter next(fTasks);
614    // Call Terminate() for tasks
615    while ((task=(AliAnalysisTask*)next())) task->Terminate();
616    //
617    TIter next1(fOutputs);
618    AliAnalysisDataContainer *output;
619    while ((output=(AliAnalysisDataContainer*)next1())) {
620       // Close all files at output
621       // Special outputs have the files already closed and written.
622       if (output->IsSpecialOutput()) continue;
623       const char *filename = output->GetFileName();
624       if (!(strcmp(filename, "default"))) {
625          if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
626          TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
627          if (aodfile) {
628             if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
629             aodfile->Write();
630             continue;
631          }   
632       }      
633       if (!strlen(filename)) continue;
634       if (!output->GetData()) continue;
635       TFile *file = output->GetFile();
636       TDirectory *opwd = gDirectory;
637       if (file) {
638          file->cd();
639       } else {
640          file = new TFile(filename, "RECREATE");
641          if (file->IsZombie()) continue;
642          output->SetFile(file);
643       }   
644       if (fDebug > 1) printf("   writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
645       if (output->GetData()->InheritsFrom(TCollection::Class())) {
646       // If data is a collection, we set the name of the collection 
647       // as the one of the container and we save as a single key.
648          TCollection *coll = (TCollection*)output->GetData();
649          coll->SetName(output->GetName());
650          coll->Write(output->GetName(), TObject::kSingleKey);
651       } else {
652          output->GetData()->Write();
653       }      
654       file->Close();
655       if (opwd) opwd->cd();
656    }   
657
658    if (fInputEventHandler)   fInputEventHandler  ->TerminateIO();
659    if (fOutputEventHandler)  fOutputEventHandler ->TerminateIO();
660    if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
661
662    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
663    if (getsysInfo) {
664       TDirectory *cdir = gDirectory;
665       TFile f("syswatch.root", "RECREATE");
666       if (!f.IsZombie()) {
667          TTree *tree = AliSysInfo::MakeTree("syswatch.log");
668          tree->SetMarkerStyle(kCircle);
669          tree->SetMarkerColor(kBlue);
670          tree->SetMarkerSize(0.5);
671          if (!gROOT->IsBatch()) {
672             tree->SetAlias("event", "id0");
673             tree->SetAlias("memUSED", "pI.fMemVirtual");
674             tree->SetAlias("userCPU", "pI.fCpuUser");
675             TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
676             c->Divide(2,1,0.01,0.01);
677             c->cd(1);
678             tree->Draw("memUSED:event","","", 1234567890, 0);
679             c->cd(2);
680             tree->Draw("userCPU:event","","", 1234567890, 0);
681          }   
682          tree->Write();
683          f.Close();
684          delete tree;
685       }
686       if (cdir) cdir->cd();
687    }      
688    if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
689 }
690
691 //______________________________________________________________________________
692 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
693 {
694 // Adds a user task to the global list of tasks.
695    if (fTasks->FindObject(task)) {
696       Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
697       return;
698    }   
699    task->SetActive(kFALSE);
700    fTasks->Add(task);
701 }  
702
703 //______________________________________________________________________________
704 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
705 {
706 // Retreive task by name.
707    if (!fTasks) return NULL;
708    return (AliAnalysisTask*)fTasks->FindObject(name);
709 }
710
711 //______________________________________________________________________________
712 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
713                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
714 {
715 // Create a data container of a certain type. Types can be:
716 //   kExchangeContainer  = 0, used to exchange date between tasks
717 //   kInputContainer   = 1, used to store input data
718 //   kOutputContainer  = 2, used for posting results
719    if (fContainers->FindObject(name)) {
720       Error("CreateContainer","A container named %s already defined !\n",name);
721       return NULL;
722    }   
723    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
724    fContainers->Add(cont);
725    switch (type) {
726       case kInputContainer:
727          fInputs->Add(cont);
728          break;
729       case kOutputContainer:
730          fOutputs->Add(cont);
731          if (filename && strlen(filename)) {
732             cont->SetFileName(filename);
733             cont->SetDataOwned(kFALSE);  // data owned by the file
734          }   
735          break;
736       case kExchangeContainer:
737          break;   
738    }
739    return cont;
740 }
741          
742 //______________________________________________________________________________
743 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
744                                         AliAnalysisDataContainer *cont)
745 {
746 // Connect input of an existing task to a data container.
747    if (!fTasks->FindObject(task)) {
748       AddTask(task);
749       Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
750    } 
751    Bool_t connected = task->ConnectInput(islot, cont);
752    return connected;
753 }   
754
755 //______________________________________________________________________________
756 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
757                                         AliAnalysisDataContainer *cont)
758 {
759 // Connect output of an existing task to a data container.
760    if (!fTasks->FindObject(task)) {
761       AddTask(task);
762       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
763    } 
764    Bool_t connected = task->ConnectOutput(islot, cont);
765    return connected;
766 }   
767                                
768 //______________________________________________________________________________
769 void AliAnalysisManager::CleanContainers()
770 {
771 // Clean data from all containers that have already finished all client tasks.
772    TIter next(fContainers);
773    AliAnalysisDataContainer *cont;
774    while ((cont=(AliAnalysisDataContainer *)next())) {
775       if (cont->IsOwnedData() && 
776           cont->IsDataReady() && 
777           cont->ClientsExecuted()) cont->DeleteData();
778    }
779 }
780
781 //______________________________________________________________________________
782 Bool_t AliAnalysisManager::InitAnalysis()
783 {
784 // Initialization of analysis chain of tasks. Should be called after all tasks
785 // and data containers are properly connected
786    // Check for input/output containers
787    fInitOK = kFALSE;
788    // Check for top tasks (depending only on input data containers)
789    if (!fTasks->First()) {
790       Error("InitAnalysis", "Analysis has no tasks !");
791       return kFALSE;
792    }   
793    TIter next(fTasks);
794    AliAnalysisTask *task;
795    AliAnalysisDataContainer *cont;
796    Int_t ntop = 0;
797    Int_t nzombies = 0;
798    Bool_t iszombie = kFALSE;
799    Bool_t istop = kTRUE;
800    Int_t i;
801    while ((task=(AliAnalysisTask*)next())) {
802       istop = kTRUE;
803       iszombie = kFALSE;
804       Int_t ninputs = task->GetNinputs();
805       for (i=0; i<ninputs; i++) {
806          cont = task->GetInputSlot(i)->GetContainer();
807          if (!cont) {
808             if (!iszombie) {
809                task->SetZombie();
810                fZombies->Add(task);
811                nzombies++;
812                iszombie = kTRUE;
813             }   
814             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
815                   i, task->GetName()); 
816          }
817          if (iszombie) continue;
818          // Check if cont is an input container
819          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
820          // Connect to parent task
821       }
822       if (istop) {
823          ntop++;
824          fTopTasks->Add(task);
825       }
826    }
827    if (!ntop) {
828       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
829       return kFALSE;
830    }                        
831    // Check now if there are orphan tasks
832    for (i=0; i<ntop; i++) {
833       task = (AliAnalysisTask*)fTopTasks->At(i);
834       task->SetUsed();
835    }
836    Int_t norphans = 0;
837    next.Reset();
838    while ((task=(AliAnalysisTask*)next())) {
839       if (!task->IsUsed()) {
840          norphans++;
841          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
842       }   
843    }          
844    // Check the task hierarchy (no parent task should depend on data provided
845    // by a daughter task)
846    for (i=0; i<ntop; i++) {
847       task = (AliAnalysisTask*)fTopTasks->At(i);
848       if (task->CheckCircularDeps()) {
849          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
850          PrintStatus("dep");
851          return kFALSE;
852       }   
853    }
854    // Check that all containers feeding post-event loop tasks are in the outputs list
855    TIter nextcont(fContainers); // loop over all containers
856    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
857       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
858          if (cont->HasConsumers()) {
859          // Check if one of the consumers is post event loop
860             TIter nextconsumer(cont->GetConsumers());
861             while ((task=(AliAnalysisTask*)nextconsumer())) {
862                if (task->IsPostEventLoop()) {
863                   fOutputs->Add(cont);
864                   break;
865                }
866             }
867          }
868       }
869    }   
870    // Check if all special output containers have a file name provided
871    TIter nextout(fOutputs);
872    while ((cont=(AliAnalysisDataContainer*)nextout())) {
873       if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
874          Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
875          return kFALSE;
876       }
877    }      
878    fInitOK = kTRUE;
879    return kTRUE;
880 }   
881
882 //______________________________________________________________________________
883 void AliAnalysisManager::PrintStatus(Option_t *option) const
884 {
885 // Print task hierarchy.
886    if (!fInitOK) {
887       Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
888       return;
889    }   
890    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
891    if (getsysInfo)
892       Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
893    TIter next(fTopTasks);
894    AliAnalysisTask *task;
895    while ((task=(AliAnalysisTask*)next()))
896       task->PrintTask(option);
897 }
898
899 //______________________________________________________________________________
900 void AliAnalysisManager::ResetAnalysis()
901 {
902 // Reset all execution flags and clean containers.
903    CleanContainers();
904 }
905
906 //______________________________________________________________________________
907 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
908 {
909 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
910 // MIX. Process nentries starting from firstentry
911    if (!fInitOK) {
912       Error("StartAnalysis","Analysis manager was not initialized !");
913       return;
914    }
915    if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
916    TString anaType = type;
917    anaType.ToLower();
918    fMode = kLocalAnalysis;
919    if (anaType.Contains("proof"))     fMode = kProofAnalysis;
920    else if (anaType.Contains("grid")) fMode = kGridAnalysis;
921    else if (anaType.Contains("mix"))  fMode = kMixingAnalysis;
922
923    if (fMode == kGridAnalysis) {
924       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
925       fMode = kLocalAnalysis;
926    }
927    char line[256];
928    SetEventLoop(kFALSE);
929    // Enable event loop mode if a tree was provided
930    if (tree || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
931
932    TChain *chain = 0;
933    TString ttype = "TTree";
934    if (tree && tree->IsA() == TChain::Class()) {
935       chain = (TChain*)tree;
936       if (!chain || !chain->GetListOfFiles()->First()) {
937          Error("StartAnalysis", "Cannot process null or empty chain...");
938          return;
939       }   
940       ttype = "TChain";
941    }   
942
943    // Initialize locally all tasks (happens for all modes)
944    TIter next(fTasks);
945    AliAnalysisTask *task;
946    while ((task=(AliAnalysisTask*)next())) {
947       task->LocalInit();
948    }
949    
950    switch (fMode) {
951       case kLocalAnalysis:
952          if (!tree) {
953             TIter nextT(fTasks);
954             // Call CreateOutputObjects for all tasks
955             while ((task=(AliAnalysisTask*)nextT())) {
956                TDirectory *curdir = gDirectory;
957                task->CreateOutputObjects();
958                if (curdir) curdir->cd();
959             }   
960             ExecAnalysis();
961             Terminate();
962             return;
963          } 
964          // Run tree-based analysis via AliAnalysisSelector  
965          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
966          fSelector = new AliAnalysisSelector(this);
967          tree->Process(fSelector, "", nentries, firstentry);
968          break;
969       case kProofAnalysis:
970          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
971             printf("StartAnalysis: no PROOF!!!\n");
972             return;
973          }   
974          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
975          gROOT->ProcessLine(line);
976          if (chain) {
977             chain->SetProof();
978             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
979             chain->Process("AliAnalysisSelector", "", nentries, firstentry);
980          } else {
981             printf("StartAnalysis: no chain\n");
982             return;
983          }      
984          break;
985       case kGridAnalysis:
986          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
987          break;
988       case kMixingAnalysis:   
989          // Run event mixing analysis
990          if (!fEventPool) {
991             Error("StartAnalysis", "Cannot run event mixing without event pool");
992             return;
993          }
994          cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
995          fSelector = new AliAnalysisSelector(this);
996          while ((chain=fEventPool->GetNextChain())) {
997             next.Reset();
998             // Call NotifyBinChange for all tasks
999             while ((task=(AliAnalysisTask*)next()))
1000                if (!task->IsPostEventLoop()) task->NotifyBinChange();
1001             chain->Process(fSelector);
1002          }
1003          PackOutput(fSelector->GetOutputList());
1004          Terminate();
1005    }   
1006 }   
1007
1008 //______________________________________________________________________________
1009 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1010 {
1011 // Start analysis for this manager on a given dataset. Analysis task can be: 
1012 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1013    if (!fInitOK) {
1014       Error("StartAnalysis","Analysis manager was not initialized !");
1015       return;
1016    }
1017    if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
1018    TString anaType = type;
1019    anaType.ToLower();
1020    if (!anaType.Contains("proof")) {
1021       Error("Cannot process datasets in %s mode. Try PROOF.", type);
1022       return;
1023    }   
1024    fMode = kProofAnalysis;
1025    char line[256];
1026    SetEventLoop(kTRUE);
1027    // Set the dataset flag
1028    TObject::SetBit(kUseDataSet);
1029    fTree = 0;
1030
1031    // Initialize locally all tasks
1032    TIter next(fTasks);
1033    AliAnalysisTask *task;
1034    while ((task=(AliAnalysisTask*)next())) {
1035       task->LocalInit();
1036    }
1037    
1038    if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1039       printf("StartAnalysis: no PROOF!!!\n");
1040       return;
1041    }   
1042    sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1043    gROOT->ProcessLine(line);
1044    sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
1045    if (!gROOT->ProcessLine(line)) {
1046       Error("StartAnalysis", "Dataset %s not found", dataset);
1047       return;
1048    }   
1049    sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1050            dataset, nentries, firstentry);
1051    cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1052    gROOT->ProcessLine(line);
1053 }   
1054
1055 //______________________________________________________________________________
1056 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
1057 {
1058 // Opens a special output file used in PROOF.
1059    char line[256];
1060    if (fMode!=kProofAnalysis || !fSelector) {
1061       Error("OpenProofFile","Cannot open PROOF file %s",filename);
1062       return NULL;
1063    }   
1064    sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
1065    if (fDebug > 1) printf("=== %s\n", line);
1066    gROOT->ProcessLine(line);
1067    sprintf(line, "pf->OpenFile(\"%s\");", option);
1068    gROOT->ProcessLine(line);
1069    if (fDebug > 1) {
1070       gROOT->ProcessLine("pf->Print()");
1071       printf(" == proof file name: %s\n", gFile->GetName());
1072    }   
1073    sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1074    if (fDebug > 1) printf("=== %s\n", line);
1075    gROOT->ProcessLine(line);
1076    return gFile;
1077 }   
1078
1079 //______________________________________________________________________________
1080 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1081 {
1082 // Execute analysis.
1083    static Long64_t ncalls = 0;
1084    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1085    if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1086    ncalls++;
1087    if (!fInitOK) {
1088      Error("ExecAnalysis", "Analysis manager was not initialized !");
1089       return;
1090    }   
1091    AliAnalysisTask *task;
1092    // Check if the top tree is active.
1093    if (fTree) {
1094       TIter next(fTasks);
1095    // De-activate all tasks
1096       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1097       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1098       if (!cont) {
1099               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1100          return;
1101       }   
1102       cont->SetData(fTree); // This will notify all consumers
1103       Long64_t entry = fTree->GetTree()->GetReadEntry();
1104       
1105 //
1106 //    Call BeginEvent() for optional input/output and MC services 
1107       if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
1108       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
1109       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1110 //
1111 //    Execute the tasks
1112 //      TIter next1(cont->GetConsumers());
1113       TIter next1(fTopTasks);
1114       while ((task=(AliAnalysisTask*)next1())) {
1115          if (fDebug >1) {
1116             cout << "    Executing task " << task->GetName() << endl;
1117          }   
1118          
1119          task->ExecuteTask(option);
1120       }
1121 //
1122 //    Call FinishEvent() for optional output and MC services 
1123       if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
1124       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
1125       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1126       // Gather system information if requested
1127       if (getsysInfo && ((ncalls%fNSysInfo)==0)) 
1128          AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1129       return;
1130    }   
1131    // The event loop is not controlled by TSelector   
1132 //
1133 //  Call BeginEvent() for optional input/output and MC services 
1134    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(-1);
1135    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(-1);
1136    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1137    TIter next2(fTopTasks);
1138    while ((task=(AliAnalysisTask*)next2())) {
1139       task->SetActive(kTRUE);
1140       if (fDebug > 1) {
1141          cout << "    Executing task " << task->GetName() << endl;
1142       }   
1143       task->ExecuteTask(option);
1144    }   
1145 //
1146 // Call FinishEvent() for optional output and MC services 
1147    if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
1148    if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
1149    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1150 }
1151
1152 //______________________________________________________________________________
1153 void AliAnalysisManager::FinishAnalysis()
1154 {
1155 // Finish analysis.
1156 }