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