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