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),
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();
79 //______________________________________________________________________________
80 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
83 fInputEventHandler(NULL),
84 fOutputEventHandler(NULL),
85 fMCtruthEventHandler(NULL),
89 fInitOK(other.fInitOK),
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;
108 //______________________________________________________________________________
109 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
112 if (&other != this) {
113 TNamed::operator=(other);
114 fInputEventHandler = other.fInputEventHandler;
115 fOutputEventHandler = other.fOutputEventHandler;
116 fMCtruthEventHandler = other.fMCtruthEventHandler;
119 fNSysInfo = other.fNSysInfo;
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;
134 //______________________________________________________________________________
135 AliAnalysisManager::~AliAnalysisManager()
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;
147 //______________________________________________________________________________
148 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
150 // Read one entry of the tree or a whole branch.
152 cout << "== AliAnalysisManager::GetEntry()" << endl;
154 fCurrentEntry = entry;
155 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
158 //______________________________________________________________________________
159 void AliAnalysisManager::Init(TTree *tree)
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.
168 printf("->AliAnalysisManager::InitTree(%s)\n", tree->GetName());
171 // Call InitTree of EventHandler
172 if (fOutputEventHandler) {
173 if (fMode == kProofAnalysis) {
174 fOutputEventHandler->Init(0x0, "proof");
176 fOutputEventHandler->Init(0x0, "local");
180 if (fInputEventHandler) {
181 if (fMode == kProofAnalysis) {
182 fInputEventHandler->Init(tree, "proof");
184 fInputEventHandler->Init(tree, "local");
187 // If no input event handler we need to get the tree once
189 if(!tree->GetTree()) tree->LoadTree(0);
193 if (fMCtruthEventHandler) {
194 if (fMode == kProofAnalysis) {
195 fMCtruthEventHandler->Init(0x0, "proof");
197 fMCtruthEventHandler->Init(0x0, "local");
201 if (!fInitOK) InitAnalysis();
202 if (!fInitOK) return;
204 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
206 cout<<"Error: No top input container !" <<endl;
211 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
215 //______________________________________________________________________________
216 void AliAnalysisManager::SlaveBegin(TTree *tree)
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).
222 cout << "->AliAnalysisManager::SlaveBegin()" << endl;
225 // Call Init of EventHandler
226 if (fOutputEventHandler) {
227 if (fMode == kProofAnalysis) {
228 fOutputEventHandler->Init("proof");
230 fOutputEventHandler->Init("local");
234 if (fInputEventHandler) {
235 fInputEventHandler->SetInputTree(tree);
236 if (fMode == kProofAnalysis) {
237 fInputEventHandler->Init("proof");
239 fInputEventHandler->Init("local");
243 if (fMCtruthEventHandler) {
244 if (fMode == kProofAnalysis) {
245 fMCtruthEventHandler->Init("proof");
247 fMCtruthEventHandler->Init("local");
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();
261 cout << "<-AliAnalysisManager::SlaveBegin()" << endl;
265 //______________________________________________________________________________
266 Bool_t AliAnalysisManager::Notify()
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.
274 TFile *curfile = fTree->GetCurrentFile();
275 if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
277 AliAnalysisTask *task;
278 // Call Notify for all tasks
279 while ((task=(AliAnalysisTask*)next()))
282 // Call Notify of the event handlers
283 if (fInputEventHandler) {
284 fInputEventHandler->Notify(curfile->GetName());
287 if (fOutputEventHandler) {
288 fOutputEventHandler->Notify(curfile->GetName());
291 if (fMCtruthEventHandler) {
292 fMCtruthEventHandler->Notify(curfile->GetName());
299 //______________________________________________________________________________
300 Bool_t AliAnalysisManager::Process(Long64_t entry)
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.
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.
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).
320 cout << "->AliAnalysisManager::Process(" << entry << ")" << endl;
322 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
323 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
324 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
329 cout << "<-AliAnalysisManager::Process()" << endl;
334 //______________________________________________________________________________
335 void AliAnalysisManager::PackOutput(TList *target)
337 // Pack all output data containers in the output list. Called at SlaveTerminate
338 // stage in PROOF case for each slave.
340 cout << "->AliAnalysisManager::PackOutput()" << endl;
343 Error("PackOutput", "No target. Aborting.");
346 if (fInputEventHandler) fInputEventHandler ->Terminate();
347 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
348 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
350 if (fMode == kProofAnalysis) {
351 TIter next(fOutputs);
352 AliAnalysisDataContainer *output;
353 while ((output=(AliAnalysisDataContainer*)next())) {
354 if (output->GetData()) {
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());
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());
367 if (strlen(fSpecialOutputLocation.Data())) {
368 TString remote = fSpecialOutputLocation;
370 remote += gSystem->HostName();
371 remote += output->GetFileName();
372 TFile::Cp(output->GetFileName(), remote.Data());
378 printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
382 //______________________________________________________________________________
383 void AliAnalysisManager::ImportWrappers(TList *source)
385 // Import data in output containers from wrappers coming in source.
387 cout << "->AliAnalysisManager::ImportWrappers()" << endl;
389 TIter next(fOutputs);
390 AliAnalysisDataContainer *cont;
391 AliAnalysisDataWrapper *wrap;
393 while ((cont=(AliAnalysisDataContainer*)next())) {
394 if (cont->GetProducer()->IsPostEventLoop()) continue;
395 wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
396 if (!wrap && fDebug>1) {
397 printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName());
401 if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName());
402 if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName());
403 cont->ImportData(wrap);
406 cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl;
410 //______________________________________________________________________________
411 void AliAnalysisManager::UnpackOutput(TList *source)
413 // Called by AliAnalysisSelector::Terminate. Output containers should
414 // be in source in the same order as in fOutputs.
416 cout << "->AliAnalysisManager::UnpackOutput()" << endl;
419 Error("UnpackOutput", "No target. Aborting.");
423 printf(" Source list contains %d containers\n", source->GetSize());
426 if (fMode == kProofAnalysis) ImportWrappers(source);
428 TIter next(fOutputs);
429 AliAnalysisDataContainer *output;
430 while ((output=(AliAnalysisDataContainer*)next())) {
431 if (!output->GetData()) continue;
432 // Check if there are client tasks that run post event loop
433 if (output->HasConsumers()) {
434 // Disable event loop semaphore
435 output->SetPostEventLoop(kTRUE);
436 TObjArray *list = output->GetConsumers();
437 Int_t ncons = list->GetEntriesFast();
438 for (Int_t i=0; i<ncons; i++) {
439 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
440 task->CheckNotify(kTRUE);
441 // If task is active, execute it
442 if (task->IsPostEventLoop() && task->IsActive()) {
444 cout << "== Executing post event loop task " << task->GetName() << endl;
450 // Check if the output need to be written to a file.
451 const char *filename = output->GetFileName();
452 if (!(strcmp(filename, "default"))) {
453 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
456 if (!filename || !strlen(filename)) continue;
457 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
458 TDirectory *opwd = gDirectory;
459 if (file) file->cd();
460 else file = new TFile(filename, "RECREATE");
461 if (file->IsZombie()) continue;
462 // Reparent data to this file
464 if (output->GetData()->IsA())
465 callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*");
466 if (callEnv.IsValid()) {
467 callEnv.SetParam((Long_t) file);
468 callEnv.Execute(output->GetData());
470 output->GetData()->Write();
471 if (opwd) opwd->cd();
474 cout << "<-AliAnalysisManager::UnpackOutput()" << endl;
478 //______________________________________________________________________________
479 void AliAnalysisManager::Terminate()
481 // The Terminate() function is the last function to be called during
482 // a query. It always runs on the client, it can be used to present
483 // the results graphically.
485 cout << "->AliAnalysisManager::Terminate()" << endl;
487 AliAnalysisTask *task;
489 // Call Terminate() for tasks
490 while ((task=(AliAnalysisTask*)next())) task->Terminate();
492 cout << "<-AliAnalysisManager::Terminate()" << endl;
495 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
496 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
497 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
498 TIter next1(fOutputs);
499 AliAnalysisDataContainer *output;
500 while ((output=(AliAnalysisDataContainer*)next1())) {
501 // Close all files at output
502 const char *filename = output->GetFileName();
503 if (!(strcmp(filename, "default"))) {
504 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
507 if (!filename || !strlen(filename)) continue;
508 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
509 if (!file || file->IsZombie()) continue;
513 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
515 TDirectory *cdir = gDirectory;
516 TFile f("syswatch.root", "RECREATE");
518 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
519 tree->SetMarkerStyle(kCircle);
520 tree->SetMarkerColor(kBlue);
521 tree->SetMarkerSize(0.5);
522 if (!gROOT->IsBatch()) {
523 tree->SetAlias("event", "id0");
524 tree->SetAlias("memUSED", "mi.fMemUsed");
525 tree->SetAlias("userCPU", "pI.fCpuUser");
526 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
527 c->Divide(2,1,0.01,0.01);
529 tree->Draw("memUSED:event","","", 1234567890, 0);
531 tree->Draw("userCPU:event","","", 1234567890, 0);
537 if (cdir) cdir->cd();
541 //______________________________________________________________________________
542 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
544 // Adds a user task to the global list of tasks.
545 task->SetActive(kFALSE);
549 //______________________________________________________________________________
550 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
552 // Retreive task by name.
553 if (!fTasks) return NULL;
554 return (AliAnalysisTask*)fTasks->FindObject(name);
557 //______________________________________________________________________________
558 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
559 TClass *datatype, EAliAnalysisContType type, const char *filename)
561 // Create a data container of a certain type. Types can be:
562 // kExchangeContainer = 0, used to exchange date between tasks
563 // kInputContainer = 1, used to store input data
564 // kOutputContainer = 2, used for posting results
565 if (fContainers->FindObject(name)) {
566 Error("CreateContainer","A container named %s already defined !\n",name);
569 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
570 fContainers->Add(cont);
572 case kInputContainer:
575 case kOutputContainer:
577 if (filename && strlen(filename)) {
578 cont->SetFileName(filename);
579 cont->SetDataOwned(kFALSE); // data owned by the file
582 case kExchangeContainer:
588 //______________________________________________________________________________
589 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
590 AliAnalysisDataContainer *cont)
592 // Connect input of an existing task to a data container.
593 if (!fTasks->FindObject(task)) {
595 Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName());
597 Bool_t connected = task->ConnectInput(islot, cont);
601 //______________________________________________________________________________
602 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
603 AliAnalysisDataContainer *cont)
605 // Connect output of an existing task to a data container.
606 if (!fTasks->FindObject(task)) {
608 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
610 Bool_t connected = task->ConnectOutput(islot, cont);
614 //______________________________________________________________________________
615 void AliAnalysisManager::CleanContainers()
617 // Clean data from all containers that have already finished all client tasks.
618 TIter next(fContainers);
619 AliAnalysisDataContainer *cont;
620 while ((cont=(AliAnalysisDataContainer *)next())) {
621 if (cont->IsOwnedData() &&
622 cont->IsDataReady() &&
623 cont->ClientsExecuted()) cont->DeleteData();
627 //______________________________________________________________________________
628 Bool_t AliAnalysisManager::InitAnalysis()
630 // Initialization of analysis chain of tasks. Should be called after all tasks
631 // and data containers are properly connected
632 // Check for input/output containers
634 // Check for top tasks (depending only on input data containers)
635 if (!fTasks->First()) {
636 Error("InitAnalysis", "Analysis has no tasks !");
640 AliAnalysisTask *task;
641 AliAnalysisDataContainer *cont;
644 Bool_t iszombie = kFALSE;
645 Bool_t istop = kTRUE;
647 while ((task=(AliAnalysisTask*)next())) {
650 Int_t ninputs = task->GetNinputs();
651 for (i=0; i<ninputs; i++) {
652 cont = task->GetInputSlot(i)->GetContainer();
660 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
663 if (iszombie) continue;
664 // Check if cont is an input container
665 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
666 // Connect to parent task
670 fTopTasks->Add(task);
674 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
677 // Check now if there are orphan tasks
678 for (i=0; i<ntop; i++) {
679 task = (AliAnalysisTask*)fTopTasks->At(i);
684 while ((task=(AliAnalysisTask*)next())) {
685 if (!task->IsUsed()) {
687 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
690 // Check the task hierarchy (no parent task should depend on data provided
691 // by a daughter task)
692 for (i=0; i<ntop; i++) {
693 task = (AliAnalysisTask*)fTopTasks->At(i);
694 if (task->CheckCircularDeps()) {
695 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
700 // Check that all containers feeding post-event loop tasks are in the outputs list
701 TIter nextcont(fContainers); // loop over all containers
702 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
703 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
704 if (cont->HasConsumers()) {
705 // Check if one of the consumers is post event loop
706 TIter nextconsumer(cont->GetConsumers());
707 while ((task=(AliAnalysisTask*)nextconsumer())) {
708 if (task->IsPostEventLoop()) {
720 //______________________________________________________________________________
721 void AliAnalysisManager::PrintStatus(Option_t *option) const
723 // Print task hierarchy.
725 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
728 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
730 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
731 TIter next(fTopTasks);
732 AliAnalysisTask *task;
733 while ((task=(AliAnalysisTask*)next()))
734 task->PrintTask(option);
737 //______________________________________________________________________________
738 void AliAnalysisManager::ResetAnalysis()
740 // Reset all execution flags and clean containers.
744 //______________________________________________________________________________
745 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
747 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID.
748 // Process nentries starting from firstentry
750 Error("StartAnalysis","Analysis manager was not initialized !");
754 cout << "StartAnalysis: " << GetName() << endl;
756 TString anaType = type;
758 fMode = kLocalAnalysis;
760 if (anaType.Contains("proof")) fMode = kProofAnalysis;
761 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
763 if (fMode == kGridAnalysis) {
764 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
765 fMode = kLocalAnalysis;
768 SetEventLoop(kFALSE);
769 // Disable all branches if requested and set event loop mode
771 if (TestBit(kDisableBranches)) {
772 printf("Disabling all branches...\n");
773 // tree->SetBranchStatus("*",0); // not yet working
779 TString ttype = "TTree";
780 if (tree->IsA() == TChain::Class()) {
781 chain = (TChain*)tree;
785 // Initialize locally all tasks
787 AliAnalysisTask *task;
788 while ((task=(AliAnalysisTask*)next())) {
796 AliAnalysisTask *task;
797 // Call CreateOutputObjects for all tasks
798 while ((task=(AliAnalysisTask*)next())) {
799 TDirectory *curdir = gDirectory;
800 task->CreateOutputObjects();
801 if (curdir) curdir->cd();
807 // Run tree-based analysis via AliAnalysisSelector
808 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
809 sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this);
810 gROOT->ProcessLine(line);
811 sprintf(line, "((%s*)0x%lx)->Process(selector, \"\",%lld, %lld);",ttype.Data(),(ULong_t)tree, nentries, firstentry);
812 gROOT->ProcessLine(line);
815 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
816 printf("StartAnalysis: no PROOF!!!\n");
819 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
820 gROOT->ProcessLine(line);
823 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
824 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
826 printf("StartAnalysis: no chain\n");
831 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
835 //______________________________________________________________________________
836 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
838 // Start analysis for this manager on a given dataset. Analysis task can be:
839 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
841 Error("StartAnalysis","Analysis manager was not initialized !");
845 cout << "StartAnalysis: " << GetName() << endl;
847 TString anaType = type;
849 if (!anaType.Contains("proof")) {
850 Error("Cannot process datasets in %s mode. Try PROOF.", type);
853 fMode = kProofAnalysis;
856 // Set the dataset flag
857 TObject::SetBit(kUseDataSet);
860 // Initialize locally all tasks
862 AliAnalysisTask *task;
863 while ((task=(AliAnalysisTask*)next())) {
867 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
868 printf("StartAnalysis: no PROOF!!!\n");
871 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
872 gROOT->ProcessLine(line);
873 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
874 if (!gROOT->ProcessLine(line)) {
875 Error("StartAnalysis", "Dataset %s not found", dataset);
878 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
879 dataset, nentries, firstentry);
880 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
881 gROOT->ProcessLine(line);
884 //______________________________________________________________________________
885 void AliAnalysisManager::ExecAnalysis(Option_t *option)
888 static Long64_t ncalls = 0;
889 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
890 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
893 Error("ExecAnalysis", "Analysis manager was not initialized !");
896 AliAnalysisTask *task;
897 // Check if the top tree is active.
900 // De-activate all tasks
901 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
902 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
904 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
907 cont->SetData(fTree); // This will notify all consumers
908 Long64_t entry = fTree->GetTree()->GetReadEntry();
911 // Call BeginEvent() for optional input/output and MC services
912 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
913 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
914 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
917 TIter next1(cont->GetConsumers());
918 while ((task=(AliAnalysisTask*)next1())) {
920 cout << " Executing task " << task->GetName() << endl;
923 task->ExecuteTask(option);
926 // Call FinishEvent() for optional output and MC services
927 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
928 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
929 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
930 // Gather system information if requested
931 if (getsysInfo && ((ncalls%fNSysInfo)==0))
932 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
935 // The event loop is not controlled by TSelector
937 // Call BeginEvent() for optional input/output and MC services
938 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
939 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
940 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
941 TIter next2(fTopTasks);
942 while ((task=(AliAnalysisTask*)next2())) {
943 task->SetActive(kTRUE);
945 cout << " Executing task " << task->GetName() << endl;
947 task->ExecuteTask(option);
950 // Call FinishEvent() for optional output and MC services
951 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
952 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
953 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
956 //______________________________________________________________________________
957 void AliAnalysisManager::FinishAnalysis()