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