1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 // Author: Andrei Gheata, 31/05/2006
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 //==============================================================================
26 //==============================================================================
28 #include <Riostream.h>
32 #include <TMethodCall.h>
38 #include "AliAnalysisTask.h"
39 #include "AliAnalysisDataContainer.h"
40 #include "AliAnalysisDataSlot.h"
41 #include "AliVEventHandler.h"
42 #include "AliSysInfo.h"
43 #include "AliAnalysisManager.h"
45 ClassImp(AliAnalysisManager)
47 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
49 //______________________________________________________________________________
50 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
53 fInputEventHandler(NULL),
54 fOutputEventHandler(NULL),
55 fMCtruthEventHandler(NULL),
58 fMode(kLocalAnalysis),
61 fSpecialOutputLocation(""),
69 // Default constructor.
70 fgAnalysisManager = this;
71 fTasks = new TObjArray();
72 fTopTasks = new TObjArray();
73 fZombies = new TObjArray();
74 fContainers = new TObjArray();
75 fInputs = new TObjArray();
76 fOutputs = new TObjArray();
80 //______________________________________________________________________________
81 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
84 fInputEventHandler(NULL),
85 fOutputEventHandler(NULL),
86 fMCtruthEventHandler(NULL),
90 fInitOK(other.fInitOK),
92 fSpecialOutputLocation(""),
101 fTasks = new TObjArray(*other.fTasks);
102 fTopTasks = new TObjArray(*other.fTopTasks);
103 fZombies = new TObjArray(*other.fZombies);
104 fContainers = new TObjArray(*other.fContainers);
105 fInputs = new TObjArray(*other.fInputs);
106 fOutputs = new TObjArray(*other.fOutputs);
107 fgAnalysisManager = this;
110 //______________________________________________________________________________
111 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
114 if (&other != this) {
115 TNamed::operator=(other);
116 fInputEventHandler = other.fInputEventHandler;
117 fOutputEventHandler = other.fOutputEventHandler;
118 fMCtruthEventHandler = other.fMCtruthEventHandler;
121 fNSysInfo = other.fNSysInfo;
123 fInitOK = other.fInitOK;
124 fDebug = other.fDebug;
125 fTasks = new TObjArray(*other.fTasks);
126 fTopTasks = new TObjArray(*other.fTopTasks);
127 fZombies = new TObjArray(*other.fZombies);
128 fContainers = new TObjArray(*other.fContainers);
129 fInputs = new TObjArray(*other.fInputs);
130 fOutputs = new TObjArray(*other.fOutputs);
131 fgAnalysisManager = this;
136 //______________________________________________________________________________
137 AliAnalysisManager::~AliAnalysisManager()
140 if (fTasks) {fTasks->Delete(); delete fTasks;}
141 if (fTopTasks) delete fTopTasks;
142 if (fZombies) delete fZombies;
143 if (fContainers) {fContainers->Delete(); delete fContainers;}
144 if (fInputs) delete fInputs;
145 if (fOutputs) delete fOutputs;
146 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
149 //______________________________________________________________________________
150 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
152 // Read one entry of the tree or a whole branch.
154 cout << "== AliAnalysisManager::GetEntry()" << endl;
156 fCurrentEntry = entry;
157 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
160 //______________________________________________________________________________
161 void AliAnalysisManager::Init(TTree *tree)
163 // The Init() function is called when the selector needs to initialize
164 // a new tree or chain. Typically here the branch addresses of the tree
165 // will be set. It is normaly not necessary to make changes to the
166 // generated code, but the routine can be extended by the user if needed.
167 // Init() will be called many times when running with PROOF.
170 printf("->AliAnalysisManager::InitTree(%s)\n", tree->GetName());
173 // Call InitTree of EventHandler
174 if (fOutputEventHandler) {
175 if (fMode == kProofAnalysis) {
176 fOutputEventHandler->Init(0x0, "proof");
178 fOutputEventHandler->Init(0x0, "local");
182 if (fInputEventHandler) {
183 if (fMode == kProofAnalysis) {
184 fInputEventHandler->Init(tree, "proof");
186 fInputEventHandler->Init(tree, "local");
189 // If no input event handler we need to get the tree once
191 if(!tree->GetTree()) tree->LoadTree(0);
195 if (fMCtruthEventHandler) {
196 if (fMode == kProofAnalysis) {
197 fMCtruthEventHandler->Init(0x0, "proof");
199 fMCtruthEventHandler->Init(0x0, "local");
203 if (!fInitOK) InitAnalysis();
204 if (!fInitOK) return;
206 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
208 cout<<"Error: No top input container !" <<endl;
213 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
217 //______________________________________________________________________________
218 void AliAnalysisManager::SlaveBegin(TTree *tree)
220 // The SlaveBegin() function is called after the Begin() function.
221 // When running with PROOF SlaveBegin() is called on each slave server.
222 // The tree argument is deprecated (on PROOF 0 is passed).
224 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
227 // Call Init of EventHandler
228 if (fOutputEventHandler) {
229 if (fMode == kProofAnalysis) {
230 fOutputEventHandler->Init("proof");
232 fOutputEventHandler->Init("local");
236 if (fInputEventHandler) {
237 fInputEventHandler->SetInputTree(tree);
238 if (fMode == kProofAnalysis) {
239 fInputEventHandler->Init("proof");
241 fInputEventHandler->Init("local");
245 if (fMCtruthEventHandler) {
246 if (fMode == kProofAnalysis) {
247 fMCtruthEventHandler->Init("proof");
249 fMCtruthEventHandler->Init("local");
254 AliAnalysisTask *task;
255 // Call CreateOutputObjects for all tasks
256 while ((task=(AliAnalysisTask*)next())) {
257 TDirectory *curdir = gDirectory;
258 task->CreateOutputObjects();
259 if (curdir) curdir->cd();
263 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
267 //______________________________________________________________________________
268 Bool_t AliAnalysisManager::Notify()
270 // The Notify() function is called when a new file is opened. This
271 // can be either for a new TTree in a TChain or when when a new TTree
272 // is started when using PROOF. It is normaly not necessary to make changes
273 // to the generated code, but the routine can be extended by the
274 // user if needed. The return value is currently not used.
276 TFile *curfile = fTree->GetCurrentFile();
277 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
279 AliAnalysisTask *task;
280 // Call Notify for all tasks
281 while ((task=(AliAnalysisTask*)next()))
284 // Call Notify of the event handlers
285 if (fInputEventHandler) {
286 fInputEventHandler->Notify(curfile->GetName());
289 if (fOutputEventHandler) {
290 fOutputEventHandler->Notify(curfile->GetName());
293 if (fMCtruthEventHandler) {
294 fMCtruthEventHandler->Notify(curfile->GetName());
301 //______________________________________________________________________________
302 Bool_t AliAnalysisManager::Process(Long64_t entry)
304 // The Process() function is called for each entry in the tree (or possibly
305 // keyed object in the case of PROOF) to be processed. The entry argument
306 // specifies which entry in the currently loaded tree is to be processed.
307 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
308 // to read either all or the required parts of the data. When processing
309 // keyed objects with PROOF, the object is already loaded and is available
310 // via the fObject pointer.
312 // This function should contain the "body" of the analysis. It can contain
313 // simple or elaborate selection criteria, run algorithms on the data
314 // of the event and typically fill histograms.
316 // WARNING when a selector is used with a TChain, you must use
317 // the pointer to the current TTree to call GetEntry(entry).
318 // The entry is always the local entry number in the current tree.
319 // Assuming that fChain is the pointer to the TChain being processed,
320 // use fChain->GetTree()->GetEntry(entry).
322 cout << "->AliAnalysisManager::Process(" << entry << ")" << endl;
324 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
325 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
326 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
331 cout << "<-AliAnalysisManager::Process()" << endl;
336 //______________________________________________________________________________
337 void AliAnalysisManager::PackOutput(TList *target)
339 // Pack all output data containers in the output list. Called at SlaveTerminate
340 // stage in PROOF case for each slave.
342 cout << "->AliAnalysisManager::PackOutput()" << endl;
345 Error("PackOutput", "No target. Aborting.");
348 if (fInputEventHandler) fInputEventHandler ->Terminate();
349 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
350 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
352 if (fMode == kProofAnalysis) {
353 TIter next(fOutputs);
354 AliAnalysisDataContainer *output;
355 while ((output=(AliAnalysisDataContainer*)next())) {
356 if (output->GetData() && !output->IsSpecialOutput()) {
357 if (output->GetProducer()->IsPostEventLoop()) continue;
359 const char *filename = output->GetFileName();
360 if (!(strcmp(filename, "default"))) {
361 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
363 if (strlen(filename)) {
364 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
365 TDirectory *opwd = gDirectory;
366 if (file) file->cd();
367 else file = new TFile(filename, "RECREATE");
368 if (file->IsZombie()) continue;
369 // Clear file list to release object ownership to user.
370 // Save data to file, then close.
372 output->GetData()->Write();
374 // Set null directory to histograms and trees.
376 if (output->GetData()->IsA())
377 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
378 if (callEnv.IsValid()) {
379 callEnv.SetParam(Long_t(0));
380 callEnv.Execute(output->GetData());
382 // Restore current directory
383 if (opwd) opwd->cd();
385 AliAnalysisDataWrapper *wrap = output->ExportData();
386 // Output wrappers must delete data after merging (AG 13/11/07)
387 wrap->SetDeleteData(kTRUE);
388 if (fDebug > 1) printf(" Packing container %s...\n", output->GetName());
391 // Special outputs files are closed and copied on the remote location
392 if (output->IsSpecialOutput() && strlen(output->GetFileName())) {
393 TDirectory *opwd = gDirectory;
394 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(output->GetFileName());
397 if (output->GetData()) output->GetData()->Write();
399 if (opwd) opwd->cd();
400 if (strlen(fSpecialOutputLocation.Data())) {
401 TString remote = fSpecialOutputLocation;
403 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
404 remote += Form("%s_%d_", gSystem->HostName(), gid);
405 remote += output->GetFileName();
406 TFile::Cp(output->GetFileName(), remote.Data());
408 // No special location specified-> use TProofFile as merging utility
410 sprintf(line, "((TList*)0x%lx)->Add(new TProofFile(\"%s\"));",
411 (ULong_t)target, output->GetFileName());
412 gROOT->ProcessLine(line);
418 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
422 //______________________________________________________________________________
423 void AliAnalysisManager::ImportWrappers(TList *source)
425 // Import data in output containers from wrappers coming in source.
427 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
429 TIter next(fOutputs);
430 AliAnalysisDataContainer *cont;
431 AliAnalysisDataWrapper *wrap;
433 while ((cont=(AliAnalysisDataContainer*)next())) {
434 if (cont->GetProducer()->IsPostEventLoop() ||
435 cont->IsSpecialOutput()) continue;
436 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
437 if (!wrap && fDebug>1) {
438 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
442 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
443 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
444 cont->ImportData(wrap);
447 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
451 //______________________________________________________________________________
452 void AliAnalysisManager::UnpackOutput(TList *source)
454 // Called by AliAnalysisSelector::Terminate only on the client.
456 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
459 Error("UnpackOutput", "No target. Aborting.");
463 printf(" Source list contains %d containers\n", source->GetSize());
466 if (fMode == kProofAnalysis) ImportWrappers(source);
468 TIter next(fOutputs);
469 AliAnalysisDataContainer *output;
470 while ((output=(AliAnalysisDataContainer*)next())) {
471 if (!output->GetData()) continue;
472 // Check if there are client tasks that run post event loop
473 if (output->HasConsumers()) {
474 // Disable event loop semaphore
475 output->SetPostEventLoop(kTRUE);
476 TObjArray *list = output->GetConsumers();
477 Int_t ncons = list->GetEntriesFast();
478 for (Int_t i=0; i<ncons; i++) {
479 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
480 task->CheckNotify(kTRUE);
481 // If task is active, execute it
482 if (task->IsPostEventLoop() && task->IsActive()) {
484 cout << "== Executing post event loop task " << task->GetName() << endl;
490 // Check if the output need to be written to a file.
491 const char *filename = output->GetFileName();
492 if (!(strcmp(filename, "default"))) {
493 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
496 if (!filename || !strlen(filename)) continue;
497 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
498 TDirectory *opwd = gDirectory;
499 if (file) file->cd();
500 else file = new TFile(filename, "RECREATE");
501 if (file->IsZombie()) continue;
502 // Reparent data to this file
504 if (output->GetData()->IsA())
505 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
506 if (callEnv.IsValid()) {
507 callEnv.SetParam((Long_t) file);
508 callEnv.Execute(output->GetData());
510 output->GetData()->Write();
511 if (opwd) opwd->cd();
514 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
518 //______________________________________________________________________________
519 void AliAnalysisManager::Terminate()
521 // The Terminate() function is the last function to be called during
522 // a query. It always runs on the client, it can be used to present
523 // the results graphically.
525 cout << "->AliAnalysisManager::Terminate()" << endl;
527 AliAnalysisTask *task;
529 // Call Terminate() for tasks
530 while ((task=(AliAnalysisTask*)next())) task->Terminate();
532 cout << "<-AliAnalysisManager::Terminate()" << endl;
535 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
536 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
537 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
538 TIter next1(fOutputs);
539 AliAnalysisDataContainer *output;
540 while ((output=(AliAnalysisDataContainer*)next1())) {
541 // Close all files at output
542 const char *filename = output->GetFileName();
543 if (!(strcmp(filename, "default"))) {
544 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
547 if (!filename || !strlen(filename)) continue;
548 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
549 if (!file || file->IsZombie()) continue;
553 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
555 TDirectory *cdir = gDirectory;
556 TFile f("syswatch.root", "RECREATE");
558 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
559 tree->SetMarkerStyle(kCircle);
560 tree->SetMarkerColor(kBlue);
561 tree->SetMarkerSize(0.5);
562 if (!gROOT->IsBatch()) {
563 tree->SetAlias("event", "id0");
564 tree->SetAlias("memUSED", "mi.fMemUsed");
565 tree->SetAlias("userCPU", "pI.fCpuUser");
566 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
567 c->Divide(2,1,0.01,0.01);
569 tree->Draw("memUSED:event","","", 1234567890, 0);
571 tree->Draw("userCPU:event","","", 1234567890, 0);
577 if (cdir) cdir->cd();
581 //______________________________________________________________________________
582 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
584 // Adds a user task to the global list of tasks.
585 task->SetActive(kFALSE);
589 //______________________________________________________________________________
590 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
592 // Retreive task by name.
593 if (!fTasks) return NULL;
594 return (AliAnalysisTask*)fTasks->FindObject(name);
597 //______________________________________________________________________________
598 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
599 TClass *datatype, EAliAnalysisContType type, const char *filename)
601 // Create a data container of a certain type. Types can be:
602 // kExchangeContainer = 0, used to exchange date between tasks
603 // kInputContainer = 1, used to store input data
604 // kOutputContainer = 2, used for posting results
605 if (fContainers->FindObject(name)) {
606 Error("CreateContainer","A container named %s already defined !\n",name);
609 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
610 fContainers->Add(cont);
612 case kInputContainer:
615 case kOutputContainer:
617 if (filename && strlen(filename)) {
618 cont->SetFileName(filename);
619 cont->SetDataOwned(kFALSE); // data owned by the file
622 case kExchangeContainer:
628 //______________________________________________________________________________
629 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
630 AliAnalysisDataContainer *cont)
632 // Connect input of an existing task to a data container.
633 if (!fTasks->FindObject(task)) {
635 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
637 Bool_t connected = task->ConnectInput(islot, cont);
641 //______________________________________________________________________________
642 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
643 AliAnalysisDataContainer *cont)
645 // Connect output of an existing task to a data container.
646 if (!fTasks->FindObject(task)) {
648 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
650 Bool_t connected = task->ConnectOutput(islot, cont);
654 //______________________________________________________________________________
655 void AliAnalysisManager::CleanContainers()
657 // Clean data from all containers that have already finished all client tasks.
658 TIter next(fContainers);
659 AliAnalysisDataContainer *cont;
660 while ((cont=(AliAnalysisDataContainer *)next())) {
661 if (cont->IsOwnedData() &&
662 cont->IsDataReady() &&
663 cont->ClientsExecuted()) cont->DeleteData();
667 //______________________________________________________________________________
668 Bool_t AliAnalysisManager::InitAnalysis()
670 // Initialization of analysis chain of tasks. Should be called after all tasks
671 // and data containers are properly connected
672 // Check for input/output containers
674 // Check for top tasks (depending only on input data containers)
675 if (!fTasks->First()) {
676 Error("InitAnalysis", "Analysis has no tasks !");
680 AliAnalysisTask *task;
681 AliAnalysisDataContainer *cont;
684 Bool_t iszombie = kFALSE;
685 Bool_t istop = kTRUE;
687 while ((task=(AliAnalysisTask*)next())) {
690 Int_t ninputs = task->GetNinputs();
691 for (i=0; i<ninputs; i++) {
692 cont = task->GetInputSlot(i)->GetContainer();
700 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
703 if (iszombie) continue;
704 // Check if cont is an input container
705 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
706 // Connect to parent task
710 fTopTasks->Add(task);
714 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
717 // Check now if there are orphan tasks
718 for (i=0; i<ntop; i++) {
719 task = (AliAnalysisTask*)fTopTasks->At(i);
724 while ((task=(AliAnalysisTask*)next())) {
725 if (!task->IsUsed()) {
727 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
730 // Check the task hierarchy (no parent task should depend on data provided
731 // by a daughter task)
732 for (i=0; i<ntop; i++) {
733 task = (AliAnalysisTask*)fTopTasks->At(i);
734 if (task->CheckCircularDeps()) {
735 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
740 // Check that all containers feeding post-event loop tasks are in the outputs list
741 TIter nextcont(fContainers); // loop over all containers
742 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
743 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
744 if (cont->HasConsumers()) {
745 // Check if one of the consumers is post event loop
746 TIter nextconsumer(cont->GetConsumers());
747 while ((task=(AliAnalysisTask*)nextconsumer())) {
748 if (task->IsPostEventLoop()) {
760 //______________________________________________________________________________
761 void AliAnalysisManager::PrintStatus(Option_t *option) const
763 // Print task hierarchy.
765 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
768 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
770 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
771 TIter next(fTopTasks);
772 AliAnalysisTask *task;
773 while ((task=(AliAnalysisTask*)next()))
774 task->PrintTask(option);
777 //______________________________________________________________________________
778 void AliAnalysisManager::ResetAnalysis()
780 // Reset all execution flags and clean containers.
784 //______________________________________________________________________________
785 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
787 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
788 // Process nentries starting from firstentry
790 Error("StartAnalysis","Analysis manager was not initialized !");
794 cout << "StartAnalysis: " << GetName() << endl;
796 TString anaType = type;
798 fMode = kLocalAnalysis;
800 if (anaType.Contains("proof")) fMode = kProofAnalysis;
801 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
803 if (fMode == kGridAnalysis) {
804 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
805 fMode = kLocalAnalysis;
808 SetEventLoop(kFALSE);
809 // Disable all branches if requested and set event loop mode
811 if (TestBit(kDisableBranches)) {
812 printf("Disabling all branches...\n");
813 // tree->SetBranchStatus("*",0); // not yet working
819 TString ttype = "TTree";
820 if (tree->IsA() == TChain::Class()) {
821 chain = (TChain*)tree;
822 if (!chain || !chain->GetListOfFiles()->First()) {
823 Error("StartAnalysis", "Cannot process null or empty chain...");
829 // Initialize locally all tasks
831 AliAnalysisTask *task;
832 while ((task=(AliAnalysisTask*)next())) {
840 AliAnalysisTask *task;
841 // Call CreateOutputObjects for all tasks
842 while ((task=(AliAnalysisTask*)next())) {
843 TDirectory *curdir = gDirectory;
844 task->CreateOutputObjects();
845 if (curdir) curdir->cd();
851 // Run tree-based analysis via AliAnalysisSelector
852 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
853 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
854 gROOT->ProcessLine(line);
855 sprintf(line, "((%s*)0x%lx)->Process(selector, \"\",%lld, %lld);",ttype.Data(),(ULong_t)tree, nentries, firstentry);
856 gROOT->ProcessLine(line);
859 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
860 printf("StartAnalysis: no PROOF!!!\n");
863 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
864 gROOT->ProcessLine(line);
867 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
868 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
870 printf("StartAnalysis: no chain\n");
875 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
879 //______________________________________________________________________________
880 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
882 // Start analysis for this manager on a given dataset. Analysis task can be:
883 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
885 Error("StartAnalysis","Analysis manager was not initialized !");
889 cout << "StartAnalysis: " << GetName() << endl;
891 TString anaType = type;
893 if (!anaType.Contains("proof")) {
894 Error("Cannot process datasets in %s mode. Try PROOF.", type);
897 fMode = kProofAnalysis;
900 // Set the dataset flag
901 TObject::SetBit(kUseDataSet);
904 // Initialize locally all tasks
906 AliAnalysisTask *task;
907 while ((task=(AliAnalysisTask*)next())) {
911 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
912 printf("StartAnalysis: no PROOF!!!\n");
915 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
916 gROOT->ProcessLine(line);
917 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
918 if (!gROOT->ProcessLine(line)) {
919 Error("StartAnalysis", "Dataset %s not found", dataset);
922 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
923 dataset, nentries, firstentry);
924 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
925 gROOT->ProcessLine(line);
928 //______________________________________________________________________________
929 void AliAnalysisManager::ExecAnalysis(Option_t *option)
932 static Long64_t ncalls = 0;
933 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
934 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
937 Error("ExecAnalysis", "Analysis manager was not initialized !");
940 AliAnalysisTask *task;
941 // Check if the top tree is active.
944 // De-activate all tasks
945 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
946 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
948 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
951 cont->SetData(fTree); // This will notify all consumers
952 Long64_t entry = fTree->GetTree()->GetReadEntry();
955 // Call BeginEvent() for optional input/output and MC services
956 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
957 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
958 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
961 // TIter next1(cont->GetConsumers());
962 TIter next1(fTopTasks);
963 while ((task=(AliAnalysisTask*)next1())) {
965 cout << " Executing task " << task->GetName() << endl;
968 task->ExecuteTask(option);
971 // Call FinishEvent() for optional output and MC services
972 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
973 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
974 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
975 // Gather system information if requested
976 if (getsysInfo && ((ncalls%fNSysInfo)==0))
977 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
980 // The event loop is not controlled by TSelector
982 // Call BeginEvent() for optional input/output and MC services
983 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
984 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
985 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
986 TIter next2(fTopTasks);
987 while ((task=(AliAnalysisTask*)next2())) {
988 task->SetActive(kTRUE);
990 cout << " Executing task " << task->GetName() << endl;
992 task->ExecuteTask(option);
995 // Call FinishEvent() for optional output and MC services
996 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
997 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
998 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1001 //______________________________________________________________________________
1002 void AliAnalysisManager::FinishAnalysis()