]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.cxx
Updating the AliAnalysisFemtoTask to work with
[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             TDirectory *opwd = gDirectory;
365             TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(output->GetFileName());
366             if (!file) continue;
367             file->cd();
368             if (output->GetData()) output->GetData()->Write();
369             file->Close();
370             if (opwd) opwd->cd();
371             if (strlen(fSpecialOutputLocation.Data())) {
372                TString remote = fSpecialOutputLocation;
373                remote += "/";
374                Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
375                remote += Form("%s_%d_", gSystem->HostName(), gid);
376                remote += output->GetFileName();
377                TFile::Cp(output->GetFileName(), remote.Data());
378             }
379          }      
380       }
381    } 
382    if (fDebug > 1) {
383       printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
384    }
385 }
386
387 //______________________________________________________________________________
388 void AliAnalysisManager::ImportWrappers(TList *source)
389 {
390 // Import data in output containers from wrappers coming in source.
391    if (fDebug > 1) {
392       cout << "->AliAnalysisManager::ImportWrappers()" << endl;
393    }   
394    TIter next(fOutputs);
395    AliAnalysisDataContainer *cont;
396    AliAnalysisDataWrapper   *wrap;
397    Int_t icont = 0;
398    while ((cont=(AliAnalysisDataContainer*)next())) {
399       if (cont->GetProducer()->IsPostEventLoop() ||
400           cont->IsSpecialOutput()) continue;
401       wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
402       if (!wrap && fDebug>1) {
403          printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
404          continue;
405       }
406       icont++;
407       if (fDebug > 1) printf("   Importing data for container %s\n", wrap->GetName());
408       if (cont->GetFileName()) printf("    -> %s\n", cont->GetFileName());
409       cont->ImportData(wrap);
410    }         
411    if (fDebug > 1) {
412       cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
413    }   
414 }
415
416 //______________________________________________________________________________
417 void AliAnalysisManager::UnpackOutput(TList *source)
418 {
419   // Called by AliAnalysisSelector::Terminate. Output containers should
420   // be in source in the same order as in fOutputs.
421    if (fDebug > 1) {
422       cout << "->AliAnalysisManager::UnpackOutput()" << endl;
423    }   
424    if (!source) {
425       Error("UnpackOutput", "No target. Aborting.");
426       return;
427    }
428    if (fDebug > 1) {
429       printf("   Source list contains %d containers\n", source->GetSize());
430    }   
431
432    if (fMode == kProofAnalysis) ImportWrappers(source);
433
434    TIter next(fOutputs);
435    AliAnalysisDataContainer *output;
436    while ((output=(AliAnalysisDataContainer*)next())) {
437       if (!output->GetData()) continue;
438       // Check if there are client tasks that run post event loop
439       if (output->HasConsumers()) {
440          // Disable event loop semaphore
441          output->SetPostEventLoop(kTRUE);
442          TObjArray *list = output->GetConsumers();
443          Int_t ncons = list->GetEntriesFast();
444          for (Int_t i=0; i<ncons; i++) {
445             AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
446             task->CheckNotify(kTRUE);
447             // If task is active, execute it
448             if (task->IsPostEventLoop() && task->IsActive()) {
449                if (fDebug > 1) {
450                   cout << "== Executing post event loop task " << task->GetName() << endl;
451                }                  
452                task->ExecuteTask();
453             }   
454          }
455       }   
456       // Check if the output need to be written to a file.
457       const char *filename = output->GetFileName();
458       if (!(strcmp(filename, "default"))) {
459           if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
460       }
461       
462       if (!filename || !strlen(filename)) continue;
463       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
464       TDirectory *opwd = gDirectory;
465       if (file) file->cd();
466       else      file = new TFile(filename, "RECREATE");
467       if (file->IsZombie()) continue;
468       // Reparent data to this file
469       TMethodCall callEnv;
470       if (output->GetData()->IsA())
471          callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
472       if (callEnv.IsValid()) {
473          callEnv.SetParam((Long_t) file);
474          callEnv.Execute(output->GetData());
475       }
476       output->GetData()->Write();
477       if (opwd) opwd->cd();
478    }
479    if (fDebug > 1) {
480       cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
481    }   
482 }
483
484 //______________________________________________________________________________
485 void AliAnalysisManager::Terminate()
486 {
487   // The Terminate() function is the last function to be called during
488   // a query. It always runs on the client, it can be used to present
489   // the results graphically.
490    if (fDebug > 1) {
491       cout << "->AliAnalysisManager::Terminate()" << endl;
492    }   
493    AliAnalysisTask *task;
494    TIter next(fTasks);
495    // Call Terminate() for tasks
496    while ((task=(AliAnalysisTask*)next())) task->Terminate();
497    if (fDebug > 1) {
498       cout << "<-AliAnalysisManager::Terminate()" << endl;
499    }   
500    //
501    if (fInputEventHandler)   fInputEventHandler  ->TerminateIO();
502    if (fOutputEventHandler)  fOutputEventHandler ->TerminateIO();
503    if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
504    TIter next1(fOutputs);
505    AliAnalysisDataContainer *output;
506    while ((output=(AliAnalysisDataContainer*)next1())) {
507       // Close all files at output
508       const char *filename = output->GetFileName();
509       if (!(strcmp(filename, "default"))) {
510          if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
511       }
512       
513       if (!filename || !strlen(filename)) continue;
514       TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
515       if (!file || file->IsZombie()) continue;
516       file->Close();
517    }   
518
519    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
520    if (getsysInfo) {
521       TDirectory *cdir = gDirectory;
522       TFile f("syswatch.root", "RECREATE");
523       if (!f.IsZombie()) {
524          TTree *tree = AliSysInfo::MakeTree("syswatch.log");
525          tree->SetMarkerStyle(kCircle);
526          tree->SetMarkerColor(kBlue);
527          tree->SetMarkerSize(0.5);
528          if (!gROOT->IsBatch()) {
529             tree->SetAlias("event", "id0");
530             tree->SetAlias("memUSED", "mi.fMemUsed");
531             tree->SetAlias("userCPU", "pI.fCpuUser");
532             TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
533             c->Divide(2,1,0.01,0.01);
534             c->cd(1);
535             tree->Draw("memUSED:event","","", 1234567890, 0);
536             c->cd(2);
537             tree->Draw("userCPU:event","","", 1234567890, 0);
538          }   
539          tree->Write();
540          f.Close();
541          delete tree;
542       }
543       if (cdir) cdir->cd();
544    }      
545 }
546
547 //______________________________________________________________________________
548 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
549 {
550 // Adds a user task to the global list of tasks.
551    task->SetActive(kFALSE);
552    fTasks->Add(task);
553 }  
554
555 //______________________________________________________________________________
556 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
557 {
558 // Retreive task by name.
559    if (!fTasks) return NULL;
560    return (AliAnalysisTask*)fTasks->FindObject(name);
561 }
562
563 //______________________________________________________________________________
564 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, 
565                                 TClass *datatype, EAliAnalysisContType type, const char *filename)
566 {
567 // Create a data container of a certain type. Types can be:
568 //   kExchangeContainer  = 0, used to exchange date between tasks
569 //   kInputContainer   = 1, used to store input data
570 //   kOutputContainer  = 2, used for posting results
571    if (fContainers->FindObject(name)) {
572       Error("CreateContainer","A container named %s already defined !\n",name);
573       return NULL;
574    }   
575    AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
576    fContainers->Add(cont);
577    switch (type) {
578       case kInputContainer:
579          fInputs->Add(cont);
580          break;
581       case kOutputContainer:
582          fOutputs->Add(cont);
583          if (filename && strlen(filename)) {
584             cont->SetFileName(filename);
585             cont->SetDataOwned(kFALSE);  // data owned by the file
586          }   
587          break;
588       case kExchangeContainer:
589          break;   
590    }
591    return cont;
592 }
593          
594 //______________________________________________________________________________
595 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
596                                         AliAnalysisDataContainer *cont)
597 {
598 // Connect input of an existing task to a data container.
599    if (!fTasks->FindObject(task)) {
600       AddTask(task);
601       Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
602    } 
603    Bool_t connected = task->ConnectInput(islot, cont);
604    return connected;
605 }   
606
607 //______________________________________________________________________________
608 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
609                                         AliAnalysisDataContainer *cont)
610 {
611 // Connect output of an existing task to a data container.
612    if (!fTasks->FindObject(task)) {
613       AddTask(task);
614       Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
615    } 
616    Bool_t connected = task->ConnectOutput(islot, cont);
617    return connected;
618 }   
619                                
620 //______________________________________________________________________________
621 void AliAnalysisManager::CleanContainers()
622 {
623 // Clean data from all containers that have already finished all client tasks.
624    TIter next(fContainers);
625    AliAnalysisDataContainer *cont;
626    while ((cont=(AliAnalysisDataContainer *)next())) {
627       if (cont->IsOwnedData() && 
628           cont->IsDataReady() && 
629           cont->ClientsExecuted()) cont->DeleteData();
630    }
631 }
632
633 //______________________________________________________________________________
634 Bool_t AliAnalysisManager::InitAnalysis()
635 {
636 // Initialization of analysis chain of tasks. Should be called after all tasks
637 // and data containers are properly connected
638    // Check for input/output containers
639    fInitOK = kFALSE;
640    // Check for top tasks (depending only on input data containers)
641    if (!fTasks->First()) {
642       Error("InitAnalysis", "Analysis has no tasks !");
643       return kFALSE;
644    }   
645    TIter next(fTasks);
646    AliAnalysisTask *task;
647    AliAnalysisDataContainer *cont;
648    Int_t ntop = 0;
649    Int_t nzombies = 0;
650    Bool_t iszombie = kFALSE;
651    Bool_t istop = kTRUE;
652    Int_t i;
653    while ((task=(AliAnalysisTask*)next())) {
654       istop = kTRUE;
655       iszombie = kFALSE;
656       Int_t ninputs = task->GetNinputs();
657       for (i=0; i<ninputs; i++) {
658          cont = task->GetInputSlot(i)->GetContainer();
659          if (!cont) {
660             if (!iszombie) {
661                task->SetZombie();
662                fZombies->Add(task);
663                nzombies++;
664                iszombie = kTRUE;
665             }   
666             Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", 
667                   i, task->GetName()); 
668          }
669          if (iszombie) continue;
670          // Check if cont is an input container
671          if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
672          // Connect to parent task
673       }
674       if (istop) {
675          ntop++;
676          fTopTasks->Add(task);
677       }
678    }
679    if (!ntop) {
680       Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
681       return kFALSE;
682    }                        
683    // Check now if there are orphan tasks
684    for (i=0; i<ntop; i++) {
685       task = (AliAnalysisTask*)fTopTasks->At(i);
686       task->SetUsed();
687    }
688    Int_t norphans = 0;
689    next.Reset();
690    while ((task=(AliAnalysisTask*)next())) {
691       if (!task->IsUsed()) {
692          norphans++;
693          Warning("InitAnalysis", "Task %s is orphan", task->GetName());
694       }   
695    }          
696    // Check the task hierarchy (no parent task should depend on data provided
697    // by a daughter task)
698    for (i=0; i<ntop; i++) {
699       task = (AliAnalysisTask*)fTopTasks->At(i);
700       if (task->CheckCircularDeps()) {
701          Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
702          PrintStatus("dep");
703          return kFALSE;
704       }   
705    }
706    // Check that all containers feeding post-event loop tasks are in the outputs list
707    TIter nextcont(fContainers); // loop over all containers
708    while ((cont=(AliAnalysisDataContainer*)nextcont())) {
709       if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
710          if (cont->HasConsumers()) {
711          // Check if one of the consumers is post event loop
712             TIter nextconsumer(cont->GetConsumers());
713             while ((task=(AliAnalysisTask*)nextconsumer())) {
714                if (task->IsPostEventLoop()) {
715                   fOutputs->Add(cont);
716                   break;
717                }
718             }
719          }
720       }
721    }   
722    fInitOK = kTRUE;
723    return kTRUE;
724 }   
725
726 //______________________________________________________________________________
727 void AliAnalysisManager::PrintStatus(Option_t *option) const
728 {
729 // Print task hierarchy.
730    if (!fInitOK) {
731       Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
732       return;
733    }   
734    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
735    if (getsysInfo)
736       Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
737    TIter next(fTopTasks);
738    AliAnalysisTask *task;
739    while ((task=(AliAnalysisTask*)next()))
740       task->PrintTask(option);
741 }
742
743 //______________________________________________________________________________
744 void AliAnalysisManager::ResetAnalysis()
745 {
746 // Reset all execution flags and clean containers.
747    CleanContainers();
748 }
749
750 //______________________________________________________________________________
751 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
752 {
753 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
754 // Process nentries starting from firstentry
755    if (!fInitOK) {
756       Error("StartAnalysis","Analysis manager was not initialized !");
757       return;
758    }
759    if (fDebug>1) {
760       cout << "StartAnalysis: " << GetName() << endl;   
761    }   
762    TString anaType = type;
763    anaType.ToLower();
764    fMode = kLocalAnalysis;
765    if (tree) {
766       if (anaType.Contains("proof"))     fMode = kProofAnalysis;
767       else if (anaType.Contains("grid")) fMode = kGridAnalysis;
768    }
769    if (fMode == kGridAnalysis) {
770       Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
771       fMode = kLocalAnalysis;
772    }
773    char line[256];
774    SetEventLoop(kFALSE);
775    // Disable all branches if requested and set event loop mode
776    if (tree) {
777       if (TestBit(kDisableBranches)) {
778          printf("Disabling all branches...\n");
779 //         tree->SetBranchStatus("*",0); // not yet working
780       }   
781       SetEventLoop(kTRUE);
782    }   
783
784    TChain *chain = 0;
785    TString ttype = "TTree";
786    if (tree->IsA() == TChain::Class()) {
787       chain = (TChain*)tree;
788       if (!chain || !chain->GetListOfFiles()->First()) {
789          Error("StartAnalysis", "Cannot process null or empty chain...");
790          return;
791       }   
792       ttype = "TChain";
793    }   
794
795    // Initialize locally all tasks
796    TIter next(fTasks);
797    AliAnalysisTask *task;
798    while ((task=(AliAnalysisTask*)next())) {
799       task->LocalInit();
800    }
801    
802    switch (fMode) {
803       case kLocalAnalysis:
804          if (!tree) {
805             TIter next(fTasks);
806             AliAnalysisTask *task;
807             // Call CreateOutputObjects for all tasks
808             while ((task=(AliAnalysisTask*)next())) {
809                TDirectory *curdir = gDirectory;
810                task->CreateOutputObjects();
811                if (curdir) curdir->cd();
812             }   
813             ExecAnalysis();
814             Terminate();
815             return;
816          } 
817          // Run tree-based analysis via AliAnalysisSelector  
818          cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
819          sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
820          gROOT->ProcessLine(line);
821          sprintf(line, "((%s*)0x%lx)->Process(selector, \"\",%lld, %lld);",ttype.Data(),(ULong_t)tree, nentries, firstentry);
822          gROOT->ProcessLine(line);
823          break;
824       case kProofAnalysis:
825          if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
826             printf("StartAnalysis: no PROOF!!!\n");
827             return;
828          }   
829          sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
830          gROOT->ProcessLine(line);
831          if (chain) {
832             chain->SetProof();
833             cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
834             chain->Process("AliAnalysisSelector", "", nentries, firstentry);
835          } else {
836             printf("StartAnalysis: no chain\n");
837             return;
838          }      
839          break;
840       case kGridAnalysis:
841          Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
842    }   
843 }   
844
845 //______________________________________________________________________________
846 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
847 {
848 // Start analysis for this manager on a given dataset. Analysis task can be: 
849 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
850    if (!fInitOK) {
851       Error("StartAnalysis","Analysis manager was not initialized !");
852       return;
853    }
854    if (fDebug>1) {
855       cout << "StartAnalysis: " << GetName() << endl;   
856    }   
857    TString anaType = type;
858    anaType.ToLower();
859    if (!anaType.Contains("proof")) {
860       Error("Cannot process datasets in %s mode. Try PROOF.", type);
861       return;
862    }   
863    fMode = kProofAnalysis;
864    char line[256];
865    SetEventLoop(kTRUE);
866    // Set the dataset flag
867    TObject::SetBit(kUseDataSet);
868    fTree = 0;
869
870    // Initialize locally all tasks
871    TIter next(fTasks);
872    AliAnalysisTask *task;
873    while ((task=(AliAnalysisTask*)next())) {
874       task->LocalInit();
875    }
876    
877    if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
878       printf("StartAnalysis: no PROOF!!!\n");
879       return;
880    }   
881    sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
882    gROOT->ProcessLine(line);
883    sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
884    if (!gROOT->ProcessLine(line)) {
885       Error("StartAnalysis", "Dataset %s not found", dataset);
886       return;
887    }   
888    sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
889            dataset, nentries, firstentry);
890    cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
891    gROOT->ProcessLine(line);
892 }   
893
894 //______________________________________________________________________________
895 void AliAnalysisManager::ExecAnalysis(Option_t *option)
896 {
897 // Execute analysis.
898    static Long64_t ncalls = 0;
899    Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
900    if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
901    ncalls++;
902    if (!fInitOK) {
903      Error("ExecAnalysis", "Analysis manager was not initialized !");
904       return;
905    }   
906    AliAnalysisTask *task;
907    // Check if the top tree is active.
908    if (fTree) {
909       TIter next(fTasks);
910    // De-activate all tasks
911       while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
912       AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
913       if (!cont) {
914               Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
915          return;
916       }   
917       cont->SetData(fTree); // This will notify all consumers
918       Long64_t entry = fTree->GetTree()->GetReadEntry();
919       
920 //
921 //    Call BeginEvent() for optional input/output and MC services 
922       if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(entry);
923       if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(entry);
924       if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
925 //
926 //    Execute the tasks
927 //      TIter next1(cont->GetConsumers());
928       TIter next1(fTopTasks);
929       while ((task=(AliAnalysisTask*)next1())) {
930          if (fDebug >1) {
931             cout << "    Executing task " << task->GetName() << endl;
932          }   
933          
934          task->ExecuteTask(option);
935       }
936 //
937 //    Call FinishEvent() for optional output and MC services 
938       if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
939       if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
940       if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
941       // Gather system information if requested
942       if (getsysInfo && ((ncalls%fNSysInfo)==0)) 
943          AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
944       return;
945    }   
946    // The event loop is not controlled by TSelector   
947 //
948 //  Call BeginEvent() for optional input/output and MC services 
949    if (fInputEventHandler)   fInputEventHandler  ->BeginEvent(-1);
950    if (fOutputEventHandler)  fOutputEventHandler ->BeginEvent(-1);
951    if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
952    TIter next2(fTopTasks);
953    while ((task=(AliAnalysisTask*)next2())) {
954       task->SetActive(kTRUE);
955       if (fDebug > 1) {
956          cout << "    Executing task " << task->GetName() << endl;
957       }   
958       task->ExecuteTask(option);
959    }   
960 //
961 // Call FinishEvent() for optional output and MC services 
962    if (fInputEventHandler)   fInputEventHandler  ->FinishEvent();
963    if (fOutputEventHandler)  fOutputEventHandler ->FinishEvent();
964    if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
965 }
966
967 //______________________________________________________________________________
968 void AliAnalysisManager::FinishAnalysis()
969 {
970 // Finish analysis.
971 }