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 "AliAnalysisSelector.h"
39 #include "AliAnalysisTask.h"
40 #include "AliAnalysisDataContainer.h"
41 #include "AliAnalysisDataSlot.h"
42 #include "AliVEventHandler.h"
43 #include "AliVEventPool.h"
44 #include "AliSysInfo.h"
45 #include "AliAnalysisManager.h"
47 ClassImp(AliAnalysisManager)
49 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
51 //______________________________________________________________________________
52 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
55 fInputEventHandler(NULL),
56 fOutputEventHandler(NULL),
57 fMCtruthEventHandler(NULL),
61 fMode(kLocalAnalysis),
64 fSpecialOutputLocation(""),
73 // Default constructor.
74 fgAnalysisManager = this;
75 fTasks = new TObjArray();
76 fTopTasks = new TObjArray();
77 fZombies = new TObjArray();
78 fContainers = new TObjArray();
79 fInputs = new TObjArray();
80 fOutputs = new TObjArray();
84 //______________________________________________________________________________
85 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
88 fInputEventHandler(NULL),
89 fOutputEventHandler(NULL),
90 fMCtruthEventHandler(NULL),
95 fInitOK(other.fInitOK),
97 fSpecialOutputLocation(""),
107 fTasks = new TObjArray(*other.fTasks);
108 fTopTasks = new TObjArray(*other.fTopTasks);
109 fZombies = new TObjArray(*other.fZombies);
110 fContainers = new TObjArray(*other.fContainers);
111 fInputs = new TObjArray(*other.fInputs);
112 fOutputs = new TObjArray(*other.fOutputs);
113 fgAnalysisManager = this;
116 //______________________________________________________________________________
117 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
120 if (&other != this) {
121 TNamed::operator=(other);
122 fInputEventHandler = other.fInputEventHandler;
123 fOutputEventHandler = other.fOutputEventHandler;
124 fMCtruthEventHandler = other.fMCtruthEventHandler;
125 fEventPool = other.fEventPool;
128 fNSysInfo = other.fNSysInfo;
130 fInitOK = other.fInitOK;
131 fDebug = other.fDebug;
132 fTasks = new TObjArray(*other.fTasks);
133 fTopTasks = new TObjArray(*other.fTopTasks);
134 fZombies = new TObjArray(*other.fZombies);
135 fContainers = new TObjArray(*other.fContainers);
136 fInputs = new TObjArray(*other.fInputs);
137 fOutputs = new TObjArray(*other.fOutputs);
139 fgAnalysisManager = this;
144 //______________________________________________________________________________
145 AliAnalysisManager::~AliAnalysisManager()
148 if (fTasks) {fTasks->Delete(); delete fTasks;}
149 if (fTopTasks) delete fTopTasks;
150 if (fZombies) delete fZombies;
151 if (fContainers) {fContainers->Delete(); delete fContainers;}
152 if (fInputs) delete fInputs;
153 if (fOutputs) delete fOutputs;
154 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
157 //______________________________________________________________________________
158 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
160 // Read one entry of the tree or a whole branch.
161 if (fDebug > 0) printf("== AliAnalysisManager::GetEntry(%lld)\n", entry);
162 fCurrentEntry = entry;
163 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
166 //______________________________________________________________________________
167 void AliAnalysisManager::Init(TTree *tree)
169 // The Init() function is called when the selector needs to initialize
170 // a new tree or chain. Typically here the branch addresses of the tree
171 // will be set. It is normaly not necessary to make changes to the
172 // generated code, but the routine can be extended by the user if needed.
173 // Init() will be called many times when running with PROOF.
176 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
179 // Call InitTree of EventHandler
180 if (fOutputEventHandler) {
181 if (fMode == kProofAnalysis) {
182 fOutputEventHandler->Init(0x0, "proof");
184 fOutputEventHandler->Init(0x0, "local");
188 if (fInputEventHandler) {
189 if (fMode == kProofAnalysis) {
190 fInputEventHandler->Init(tree, "proof");
192 fInputEventHandler->Init(tree, "local");
195 // If no input event handler we need to get the tree once
197 if(!tree->GetTree()) tree->LoadTree(0);
201 if (fMCtruthEventHandler) {
202 if (fMode == kProofAnalysis) {
203 fMCtruthEventHandler->Init(0x0, "proof");
205 fMCtruthEventHandler->Init(0x0, "local");
209 if (!fInitOK) InitAnalysis();
210 if (!fInitOK) return;
212 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
214 Error("Init","No top input container !");
219 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
223 //______________________________________________________________________________
224 void AliAnalysisManager::SlaveBegin(TTree *tree)
226 // The SlaveBegin() function is called after the Begin() function.
227 // When running with PROOF SlaveBegin() is called on each slave server.
228 // The tree argument is deprecated (on PROOF 0 is passed).
229 if (fDebug > 0) printf("->AliAnalysisManager::SlaveBegin()\n");
230 static Bool_t isCalled = kFALSE;
231 TDirectory *curdir = gDirectory;
232 // Call SlaveBegin only once in case of mixing
233 if (isCalled && fMode==kMixingAnalysis) return;
234 // Call Init of EventHandler
235 if (fOutputEventHandler) {
236 if (fMode == kProofAnalysis) {
237 TIter nextout(fOutputs);
238 AliAnalysisDataContainer *c_aod;
239 while ((c_aod=(AliAnalysisDataContainer*)nextout())) if (!strcmp(c_aod->GetFileName(),"default")) break;
240 if (c_aod && c_aod->IsSpecialOutput()) {
242 if (fDebug > 1) printf(" Initializing special output file %s...\n", fOutputEventHandler->GetOutputFileName());
243 OpenProofFile(fOutputEventHandler->GetOutputFileName(), "RECREATE");
244 c_aod->SetFile(gFile);
245 fOutputEventHandler->Init("proofspecial");
248 fOutputEventHandler->Init("proof");
251 fOutputEventHandler->Init("local");
255 if (fInputEventHandler) {
256 fInputEventHandler->SetInputTree(tree);
257 if (fMode == kProofAnalysis) {
258 fInputEventHandler->Init("proof");
260 fInputEventHandler->Init("local");
264 if (fMCtruthEventHandler) {
265 if (fMode == kProofAnalysis) {
266 fMCtruthEventHandler->Init("proof");
268 fMCtruthEventHandler->Init("local");
271 if (curdir) curdir->cd();
274 AliAnalysisTask *task;
275 // Call CreateOutputObjects for all tasks
276 while ((task=(AliAnalysisTask*)next())) {
278 task->CreateOutputObjects();
279 if (curdir) curdir->cd();
283 if (fDebug > 0) printf("<-AliAnalysisManager::SlaveBegin()\n");
286 //______________________________________________________________________________
287 Bool_t AliAnalysisManager::Notify()
289 // The Notify() function is called when a new file is opened. This
290 // can be either for a new TTree in a TChain or when when a new TTree
291 // is started when using PROOF. It is normaly not necessary to make changes
292 // to the generated code, but the routine can be extended by the
293 // user if needed. The return value is currently not used.
295 Error("Notify","No current tree.");
298 TFile *curfile = fTree->GetCurrentFile();
300 Error("Notify","No current file");
304 if (fDebug > 0) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
306 AliAnalysisTask *task;
307 // Call Notify for all tasks
308 while ((task=(AliAnalysisTask*)next()))
311 // Call Notify of the event handlers
312 if (fInputEventHandler) {
313 fInputEventHandler->Notify(curfile->GetName());
316 if (fOutputEventHandler) {
317 fOutputEventHandler->Notify(curfile->GetName());
320 if (fMCtruthEventHandler) {
321 fMCtruthEventHandler->Notify(curfile->GetName());
323 if (fDebug > 0) printf("<-AliAnalysisManager::Notify()\n");
327 //______________________________________________________________________________
328 Bool_t AliAnalysisManager::Process(Long64_t entry)
330 // The Process() function is called for each entry in the tree (or possibly
331 // keyed object in the case of PROOF) to be processed. The entry argument
332 // specifies which entry in the currently loaded tree is to be processed.
333 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
334 // to read either all or the required parts of the data. When processing
335 // keyed objects with PROOF, the object is already loaded and is available
336 // via the fObject pointer.
338 // This function should contain the "body" of the analysis. It can contain
339 // simple or elaborate selection criteria, run algorithms on the data
340 // of the event and typically fill histograms.
342 // WARNING when a selector is used with a TChain, you must use
343 // the pointer to the current TTree to call GetEntry(entry).
344 // The entry is always the local entry number in the current tree.
345 // Assuming that fChain is the pointer to the TChain being processed,
346 // use fChain->GetTree()->GetEntry(entry).
347 if (fDebug > 0) printf("->AliAnalysisManager::Process(%lld)\n", entry);
349 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
350 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
351 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
355 if (fDebug > 0) printf("<-AliAnalysisManager::Process()\n");
359 //______________________________________________________________________________
360 void AliAnalysisManager::PackOutput(TList *target)
362 // Pack all output data containers in the output list. Called at SlaveTerminate
363 // stage in PROOF case for each slave.
364 if (fDebug > 0) printf("->AliAnalysisManager::PackOutput()\n");
366 Error("PackOutput", "No target. Aborting.");
369 if (fInputEventHandler) fInputEventHandler ->Terminate();
370 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
371 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
373 // Call FinishTaskOutput() for each event loop task (not called for
374 // post-event loop tasks - use Terminate() fo those)
375 TIter nexttask(fTasks);
376 AliAnalysisTask *task;
377 while ((task=(AliAnalysisTask*)nexttask())) {
378 if (!task->IsPostEventLoop()) {
379 if (fDebug > 0) printf("->FinishTaskOutput: task %s\n", task->GetName());
380 task->FinishTaskOutput();
381 if (fDebug > 0) printf("<-FinishTaskOutput: task %s\n", task->GetName());
385 if (fMode == kProofAnalysis) {
386 TIter next(fOutputs);
387 AliAnalysisDataContainer *output;
388 Bool_t isManagedByHandler = kFALSE;
389 while ((output=(AliAnalysisDataContainer*)next())) {
390 // Do not consider outputs of post event loop tasks
391 isManagedByHandler = kFALSE;
392 if (output->GetProducer()->IsPostEventLoop()) continue;
393 const char *filename = output->GetFileName();
394 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
395 isManagedByHandler = kTRUE;
396 filename = fOutputEventHandler->GetOutputFileName();
398 // Check if data was posted to this container. If not, issue an error.
399 if (!output->GetData() && !isManagedByHandler) {
400 Error("PackOutput", "No data for output container %s. Forgot to PostData ?\n", output->GetName());
403 if (!output->IsSpecialOutput()) {
405 if (strlen(filename) && !isManagedByHandler) {
406 // File resident outputs
407 TFile *file = output->GetFile();
408 // Backup current folder
409 TDirectory *opwd = gDirectory;
410 // Create file if not existing and register to container.
411 if (file) file->cd();
412 else file = new TFile(filename, "RECREATE");
413 if (file->IsZombie()) {
414 Fatal("PackOutput", "Could not recreate file %s\n", filename);
417 output->SetFile(file);
418 // Clear file list to release object ownership to user.
420 // Save data to file, then close.
421 if (output->GetData()->InheritsFrom(TCollection::Class())) {
422 // If data is a collection, we set the name of the collection
423 // as the one of the container and we save as a single key.
424 TCollection *coll = (TCollection*)output->GetData();
425 coll->SetName(output->GetName());
426 coll->Write(output->GetName(), TObject::kSingleKey);
428 if (output->GetData()->InheritsFrom(TTree::Class())) {
429 TTree *tree = (TTree*)output->GetData();
430 tree->SetDirectory(file);
433 output->GetData()->Write();
436 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
438 printf(" file %s listing content:\n", filename);
442 // Restore current directory
443 if (opwd) opwd->cd();
445 // Memory-resident outputs
446 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
448 AliAnalysisDataWrapper *wrap = 0;
449 if (isManagedByHandler) {
450 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
451 wrap->SetName(output->GetName());
453 else wrap =output->ExportData();
454 // Output wrappers must NOT delete data after merging - the user owns them
455 wrap->SetDeleteData(kFALSE);
459 TDirectory *opwd = gDirectory;
460 TFile *file = output->GetFile();
461 if (fDebug > 1 && file) printf("PackOutput %s: file merge, special output\n", output->GetName());
462 if (isManagedByHandler) {
463 // Terminate IO for files managed by the output handler
464 if (file) file->Write();
465 if (file && fDebug > 2) {
466 printf(" handled file %s listing content:\n", file->GetName());
469 fOutputEventHandler->TerminateIO();
474 AliAnalysisTask *producer = output->GetProducer();
476 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
477 output->GetFileName(), output->GetName(), producer->ClassName());
481 // Release object ownership to users after writing data to file
482 if (output->GetData()->InheritsFrom(TCollection::Class())) {
483 // If data is a collection, we set the name of the collection
484 // as the one of the container and we save as a single key.
485 TCollection *coll = (TCollection*)output->GetData();
486 coll->SetName(output->GetName());
487 coll->Write(output->GetName(), TObject::kSingleKey);
489 if (output->GetData()->InheritsFrom(TTree::Class())) {
490 TTree *tree = (TTree*)output->GetData();
491 tree->SetDirectory(file);
494 output->GetData()->Write();
499 printf(" file %s listing content:\n", output->GetFileName());
503 // Restore current directory
504 if (opwd) opwd->cd();
505 // Check if a special output location was provided or the output files have to be merged
506 if (strlen(fSpecialOutputLocation.Data())) {
507 TString remote = fSpecialOutputLocation;
509 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
510 remote += Form("%s_%d_", gSystem->HostName(), gid);
511 remote += output->GetFileName();
512 TFile::Cp(output->GetFileName(), remote.Data());
514 // No special location specified-> use TProofOutputFile as merging utility
515 // The file at this output slot must be opened in CreateOutputObjects
516 if (fDebug > 1) printf(" File %s to be merged...\n", output->GetFileName());
521 if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
524 //______________________________________________________________________________
525 void AliAnalysisManager::ImportWrappers(TList *source)
527 // Import data in output containers from wrappers coming in source.
528 if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
529 TIter next(fOutputs);
530 AliAnalysisDataContainer *cont;
531 AliAnalysisDataWrapper *wrap;
533 while ((cont=(AliAnalysisDataContainer*)next())) {
535 if (cont->GetProducer()->IsPostEventLoop()) continue;
536 const char *filename = cont->GetFileName();
537 Bool_t isManagedByHandler = kFALSE;
538 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
539 isManagedByHandler = kTRUE;
540 filename = fOutputEventHandler->GetOutputFileName();
542 if (cont->IsSpecialOutput()) {
543 if (strlen(fSpecialOutputLocation.Data()) && !isManagedByHandler) continue;
544 // Copy merged file from PROOF scratch space
546 TObject *pof = source->FindObject(filename);
547 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
548 Error("ImportWrappers", "TProofOutputFile object not found in output list for container %s", cont->GetName());
551 gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", ((TProofOutputFile*)0x%lx)->GetOutputFileName();)", full_path, pof));
553 printf(" Copying file %s from PROOF scratch space\n", full_path);
554 Bool_t gotit = TFile::Cp(full_path, filename);
556 Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
558 // Normally we should connect data from the copied file to the
559 // corresponding output container, but it is not obvious how to do this
560 // automatically if several objects in file...
561 TFile *f = new TFile(filename, "READ");
563 if (!isManagedByHandler) obj = f->Get(cont->GetName());
564 if (!obj && !isManagedByHandler) {
565 Error("ImportWrappers", "Could not find object %s in file %s", cont->GetName(), filename);
568 wrap = new AliAnalysisDataWrapper(obj);
569 wrap->SetDeleteData(kFALSE);
571 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
573 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
578 printf(" Importing data for container %s", cont->GetName());
579 if (strlen(filename)) printf(" -> file %s\n", cont->GetFileName());
582 cont->ImportData(wrap);
584 if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
587 //______________________________________________________________________________
588 void AliAnalysisManager::UnpackOutput(TList *source)
590 // Called by AliAnalysisSelector::Terminate only on the client.
591 if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
593 Error("UnpackOutput", "No target. Aborting.");
596 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
598 if (fMode == kProofAnalysis) ImportWrappers(source);
600 TIter next(fOutputs);
601 AliAnalysisDataContainer *output;
602 while ((output=(AliAnalysisDataContainer*)next())) {
603 if (!output->GetData()) continue;
604 // Check if there are client tasks that run post event loop
605 if (output->HasConsumers()) {
606 // Disable event loop semaphore
607 output->SetPostEventLoop(kTRUE);
608 TObjArray *list = output->GetConsumers();
609 Int_t ncons = list->GetEntriesFast();
610 for (Int_t i=0; i<ncons; i++) {
611 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
612 task->CheckNotify(kTRUE);
613 // If task is active, execute it
614 if (task->IsPostEventLoop() && task->IsActive()) {
615 if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
621 if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
624 //______________________________________________________________________________
625 void AliAnalysisManager::Terminate()
627 // The Terminate() function is the last function to be called during
628 // a query. It always runs on the client, it can be used to present
629 // the results graphically.
630 if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
631 AliAnalysisTask *task;
633 // Call Terminate() for tasks
634 while ((task=(AliAnalysisTask*)next())) task->Terminate();
636 TIter next1(fOutputs);
637 AliAnalysisDataContainer *output;
638 while ((output=(AliAnalysisDataContainer*)next1())) {
639 // Special outputs have the files already closed and written.
640 if (output->IsSpecialOutput()) continue;
641 const char *filename = output->GetFileName();
642 if (!(strcmp(filename, "default"))) {
643 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
644 TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
646 if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
651 if (!strlen(filename)) continue;
652 if (!output->GetData()) continue;
653 TFile *file = output->GetFile();
654 TDirectory *opwd = gDirectory;
655 file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
656 if (!file) file = new TFile(filename, "RECREATE");
657 if (file->IsZombie()) continue;
658 output->SetFile(file);
660 if (fDebug > 1) printf(" writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
661 if (output->GetData()->InheritsFrom(TCollection::Class())) {
662 // If data is a collection, we set the name of the collection
663 // as the one of the container and we save as a single key.
664 TCollection *coll = (TCollection*)output->GetData();
665 coll->SetName(output->GetName());
666 coll->Write(output->GetName(), TObject::kSingleKey);
668 if (output->GetData()->InheritsFrom(TTree::Class())) {
669 TTree *tree = (TTree*)output->GetData();
670 tree->SetDirectory(file);
673 output->GetData()->Write();
676 if (opwd) opwd->cd();
679 while ((output=(AliAnalysisDataContainer*)next1())) {
680 // Close all files at output
681 TDirectory *opwd = gDirectory;
682 if (output->GetFile()) output->GetFile()->Close();
683 if (opwd) opwd->cd();
686 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
687 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
688 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
690 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
692 TDirectory *cdir = gDirectory;
693 TFile f("syswatch.root", "RECREATE");
695 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
696 tree->SetMarkerStyle(kCircle);
697 tree->SetMarkerColor(kBlue);
698 tree->SetMarkerSize(0.5);
699 if (!gROOT->IsBatch()) {
700 tree->SetAlias("event", "id0");
701 tree->SetAlias("memUSED", "pI.fMemVirtual");
702 tree->SetAlias("userCPU", "pI.fCpuUser");
703 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
704 c->Divide(2,1,0.01,0.01);
706 tree->Draw("memUSED:event","","", 1234567890, 0);
708 tree->Draw("userCPU:event","","", 1234567890, 0);
714 if (cdir) cdir->cd();
716 if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
719 //______________________________________________________________________________
720 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
722 // Adds a user task to the global list of tasks.
723 if (fTasks->FindObject(task)) {
724 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
727 task->SetActive(kFALSE);
731 //______________________________________________________________________________
732 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
734 // Retreive task by name.
735 if (!fTasks) return NULL;
736 return (AliAnalysisTask*)fTasks->FindObject(name);
739 //______________________________________________________________________________
740 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
741 TClass *datatype, EAliAnalysisContType type, const char *filename)
743 // Create a data container of a certain type. Types can be:
744 // kExchangeContainer = 0, used to exchange date between tasks
745 // kInputContainer = 1, used to store input data
746 // kOutputContainer = 2, used for posting results
747 if (fContainers->FindObject(name)) {
748 Error("CreateContainer","A container named %s already defined !\n",name);
751 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
752 fContainers->Add(cont);
754 case kInputContainer:
757 case kOutputContainer:
759 if (filename && strlen(filename)) {
760 cont->SetFileName(filename);
761 cont->SetDataOwned(kFALSE); // data owned by the file
764 case kExchangeContainer:
770 //______________________________________________________________________________
771 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
772 AliAnalysisDataContainer *cont)
774 // Connect input of an existing task to a data container.
775 if (!fTasks->FindObject(task)) {
777 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
779 Bool_t connected = task->ConnectInput(islot, cont);
783 //______________________________________________________________________________
784 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
785 AliAnalysisDataContainer *cont)
787 // Connect output of an existing task to a data container.
788 if (!fTasks->FindObject(task)) {
790 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
792 Bool_t connected = task->ConnectOutput(islot, cont);
796 //______________________________________________________________________________
797 void AliAnalysisManager::CleanContainers()
799 // Clean data from all containers that have already finished all client tasks.
800 TIter next(fContainers);
801 AliAnalysisDataContainer *cont;
802 while ((cont=(AliAnalysisDataContainer *)next())) {
803 if (cont->IsOwnedData() &&
804 cont->IsDataReady() &&
805 cont->ClientsExecuted()) cont->DeleteData();
809 //______________________________________________________________________________
810 Bool_t AliAnalysisManager::InitAnalysis()
812 // Initialization of analysis chain of tasks. Should be called after all tasks
813 // and data containers are properly connected
814 // Check for input/output containers
816 // Check for top tasks (depending only on input data containers)
817 if (!fTasks->First()) {
818 Error("InitAnalysis", "Analysis has no tasks !");
822 AliAnalysisTask *task;
823 AliAnalysisDataContainer *cont;
826 Bool_t iszombie = kFALSE;
827 Bool_t istop = kTRUE;
829 while ((task=(AliAnalysisTask*)next())) {
832 Int_t ninputs = task->GetNinputs();
833 for (i=0; i<ninputs; i++) {
834 cont = task->GetInputSlot(i)->GetContainer();
842 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
845 if (iszombie) continue;
846 // Check if cont is an input container
847 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
848 // Connect to parent task
852 fTopTasks->Add(task);
856 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
859 // Check now if there are orphan tasks
860 for (i=0; i<ntop; i++) {
861 task = (AliAnalysisTask*)fTopTasks->At(i);
866 while ((task=(AliAnalysisTask*)next())) {
867 if (!task->IsUsed()) {
869 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
872 // Check the task hierarchy (no parent task should depend on data provided
873 // by a daughter task)
874 for (i=0; i<ntop; i++) {
875 task = (AliAnalysisTask*)fTopTasks->At(i);
876 if (task->CheckCircularDeps()) {
877 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
882 // Check that all containers feeding post-event loop tasks are in the outputs list
883 TIter nextcont(fContainers); // loop over all containers
884 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
885 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
886 if (cont->HasConsumers()) {
887 // Check if one of the consumers is post event loop
888 TIter nextconsumer(cont->GetConsumers());
889 while ((task=(AliAnalysisTask*)nextconsumer())) {
890 if (task->IsPostEventLoop()) {
898 // Check if all special output containers have a file name provided
899 TIter nextout(fOutputs);
900 while ((cont=(AliAnalysisDataContainer*)nextout())) {
901 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
902 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
910 //______________________________________________________________________________
911 void AliAnalysisManager::PrintStatus(Option_t *option) const
913 // Print task hierarchy.
915 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
918 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
920 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
921 TIter next(fTopTasks);
922 AliAnalysisTask *task;
923 while ((task=(AliAnalysisTask*)next()))
924 task->PrintTask(option);
927 //______________________________________________________________________________
928 void AliAnalysisManager::ResetAnalysis()
930 // Reset all execution flags and clean containers.
934 //______________________________________________________________________________
935 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
937 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
938 // MIX. Process nentries starting from firstentry
940 Error("StartAnalysis","Analysis manager was not initialized !");
943 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
944 TString anaType = type;
946 fMode = kLocalAnalysis;
947 if (anaType.Contains("proof")) fMode = kProofAnalysis;
948 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
949 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
951 if (fMode == kGridAnalysis) {
952 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
953 fMode = kLocalAnalysis;
956 SetEventLoop(kFALSE);
957 // Enable event loop mode if a tree was provided
958 if (tree || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
961 TString ttype = "TTree";
962 if (tree && tree->IsA() == TChain::Class()) {
963 chain = (TChain*)tree;
964 if (!chain || !chain->GetListOfFiles()->First()) {
965 Error("StartAnalysis", "Cannot process null or empty chain...");
971 // Initialize locally all tasks (happens for all modes)
973 AliAnalysisTask *task;
974 while ((task=(AliAnalysisTask*)next())) {
982 // Call CreateOutputObjects for all tasks
983 while ((task=(AliAnalysisTask*)nextT())) {
984 TDirectory *curdir = gDirectory;
985 task->CreateOutputObjects();
986 if (curdir) curdir->cd();
992 // Run tree-based analysis via AliAnalysisSelector
993 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
994 fSelector = new AliAnalysisSelector(this);
995 tree->Process(fSelector, "", nentries, firstentry);
998 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
999 printf("StartAnalysis: no PROOF!!!\n");
1002 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1003 gROOT->ProcessLine(line);
1006 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1007 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1009 printf("StartAnalysis: no chain\n");
1014 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
1016 case kMixingAnalysis:
1017 // Run event mixing analysis
1019 Error("StartAnalysis", "Cannot run event mixing without event pool");
1022 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1023 fSelector = new AliAnalysisSelector(this);
1024 while ((chain=fEventPool->GetNextChain())) {
1026 // Call NotifyBinChange for all tasks
1027 while ((task=(AliAnalysisTask*)next()))
1028 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1029 chain->Process(fSelector);
1031 PackOutput(fSelector->GetOutputList());
1036 //______________________________________________________________________________
1037 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1039 // Start analysis for this manager on a given dataset. Analysis task can be:
1040 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1042 Error("StartAnalysis","Analysis manager was not initialized !");
1045 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
1046 TString anaType = type;
1048 if (!anaType.Contains("proof")) {
1049 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1052 fMode = kProofAnalysis;
1054 SetEventLoop(kTRUE);
1055 // Set the dataset flag
1056 TObject::SetBit(kUseDataSet);
1059 // Initialize locally all tasks
1061 AliAnalysisTask *task;
1062 while ((task=(AliAnalysisTask*)next())) {
1066 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1067 printf("StartAnalysis: no PROOF!!!\n");
1070 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1071 gROOT->ProcessLine(line);
1072 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
1073 if (!gROOT->ProcessLine(line)) {
1074 Error("StartAnalysis", "Dataset %s not found", dataset);
1077 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1078 dataset, nentries, firstentry);
1079 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1080 gROOT->ProcessLine(line);
1083 //______________________________________________________________________________
1084 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
1086 // Opens a special output file used in PROOF.
1088 if (fMode!=kProofAnalysis || !fSelector) {
1089 Error("OpenProofFile","Cannot open PROOF file %s",filename);
1092 sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
1093 if (fDebug > 1) printf("=== %s\n", line);
1094 gROOT->ProcessLine(line);
1095 sprintf(line, "pf->OpenFile(\"%s\");", option);
1096 gROOT->ProcessLine(line);
1098 gROOT->ProcessLine("pf->Print()");
1099 printf(" == proof file name: %s\n", gFile->GetName());
1101 sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1102 if (fDebug > 1) printf("=== %s\n", line);
1103 gROOT->ProcessLine(line);
1107 //______________________________________________________________________________
1108 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1110 // Execute analysis.
1111 static Long64_t ncalls = 0;
1112 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1113 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1116 Error("ExecAnalysis", "Analysis manager was not initialized !");
1119 AliAnalysisTask *task;
1120 // Check if the top tree is active.
1123 // De-activate all tasks
1124 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1125 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1127 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1130 cont->SetData(fTree); // This will notify all consumers
1131 Long64_t entry = fTree->GetTree()->GetReadEntry();
1134 // Call BeginEvent() for optional input/output and MC services
1135 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1136 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1137 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1139 // Execute the tasks
1140 // TIter next1(cont->GetConsumers());
1141 TIter next1(fTopTasks);
1142 while ((task=(AliAnalysisTask*)next1())) {
1144 cout << " Executing task " << task->GetName() << endl;
1147 task->ExecuteTask(option);
1150 // Call FinishEvent() for optional output and MC services
1151 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1152 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1153 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1154 // Gather system information if requested
1155 if (getsysInfo && ((ncalls%fNSysInfo)==0))
1156 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1159 // The event loop is not controlled by TSelector
1161 // Call BeginEvent() for optional input/output and MC services
1162 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
1163 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
1164 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1165 TIter next2(fTopTasks);
1166 while ((task=(AliAnalysisTask*)next2())) {
1167 task->SetActive(kTRUE);
1169 cout << " Executing task " << task->GetName() << endl;
1171 task->ExecuteTask(option);
1174 // Call FinishEvent() for optional output and MC services
1175 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1176 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1177 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1180 //______________________________________________________________________________
1181 void AliAnalysisManager::FinishAnalysis()