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 output->GetData()->Write();
430 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
432 printf(" file %s listing content:\n", filename);
436 // Restore current directory
437 if (opwd) opwd->cd();
439 // Memory-resident outputs
440 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
442 AliAnalysisDataWrapper *wrap = 0;
443 if (isManagedByHandler) {
444 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
445 wrap->SetName(output->GetName());
447 else wrap =output->ExportData();
448 // Output wrappers must delete data after merging (AG 13/11/07)
449 wrap->SetDeleteData(kTRUE);
453 TDirectory *opwd = gDirectory;
454 TFile *file = output->GetFile();
455 if (fDebug > 1 && file) printf("PackOutput %s: file merge, special output\n", output->GetName());
456 if (isManagedByHandler) {
457 // Terminate IO for files managed by the output handler
458 if (file) file->Write();
459 if (file && fDebug > 2) {
460 printf(" handled file %s listing content:\n", file->GetName());
463 fOutputEventHandler->TerminateIO();
468 AliAnalysisTask *producer = output->GetProducer();
470 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
471 output->GetFileName(), output->GetName(), producer->ClassName());
475 // Release object ownership to users after writing data to file
476 if (output->GetData()->InheritsFrom(TCollection::Class())) {
477 // If data is a collection, we set the name of the collection
478 // as the one of the container and we save as a single key.
479 TCollection *coll = (TCollection*)output->GetData();
480 coll->SetName(output->GetName());
481 coll->Write(output->GetName(), TObject::kSingleKey);
483 output->GetData()->Write();
487 printf(" file %s listing content:\n", output->GetFileName());
491 // Restore current directory
492 if (opwd) opwd->cd();
493 // Check if a special output location was provided or the output files have to be merged
494 if (strlen(fSpecialOutputLocation.Data())) {
495 TString remote = fSpecialOutputLocation;
497 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
498 remote += Form("%s_%d_", gSystem->HostName(), gid);
499 remote += output->GetFileName();
500 TFile::Cp(output->GetFileName(), remote.Data());
502 // No special location specified-> use TProofOutputFile as merging utility
503 // The file at this output slot must be opened in CreateOutputObjects
504 if (fDebug > 1) printf(" File %s to be merged...\n", output->GetFileName());
509 if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
512 //______________________________________________________________________________
513 void AliAnalysisManager::ImportWrappers(TList *source)
515 // Import data in output containers from wrappers coming in source.
516 if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
517 TIter next(fOutputs);
518 AliAnalysisDataContainer *cont;
519 AliAnalysisDataWrapper *wrap;
521 while ((cont=(AliAnalysisDataContainer*)next())) {
523 if (cont->GetProducer()->IsPostEventLoop()) continue;
524 const char *filename = cont->GetFileName();
525 Bool_t isManagedByHandler = kFALSE;
526 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
527 isManagedByHandler = kTRUE;
528 filename = fOutputEventHandler->GetOutputFileName();
530 if (cont->IsSpecialOutput()) {
531 if (strlen(fSpecialOutputLocation.Data()) && !isManagedByHandler) continue;
532 // Copy merged file from PROOF scratch space
534 TObject *pof = source->FindObject(filename);
535 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
536 Error("ImportWrappers", "TProofOutputFile object not found in output list for container %s", cont->GetName());
539 gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", ((TProofOutputFile*)0x%lx)->GetOutputFileName();)", full_path, pof));
541 printf(" Copying file %s from PROOF scratch space\n", full_path);
542 Bool_t gotit = TFile::Cp(full_path, filename);
544 Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
546 // Normally we should connect data from the copied file to the
547 // corresponding output container, but it is not obvious how to do this
548 // automatically if several objects in file...
549 TFile *f = new TFile(filename, "READ");
551 if (!isManagedByHandler) obj = f->Get(cont->GetName());
552 if (!obj && !isManagedByHandler) {
553 Error("ImportWrappers", "Could not find object %s in file %s", cont->GetName(), filename);
556 wrap = new AliAnalysisDataWrapper(obj);
557 wrap->SetDeleteData(kFALSE);
559 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
561 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
566 printf(" Importing data for container %s", cont->GetName());
567 if (strlen(filename)) printf(" -> file %s\n", cont->GetFileName());
570 cont->ImportData(wrap);
572 if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
575 //______________________________________________________________________________
576 void AliAnalysisManager::UnpackOutput(TList *source)
578 // Called by AliAnalysisSelector::Terminate only on the client.
579 if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
581 Error("UnpackOutput", "No target. Aborting.");
584 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
586 if (fMode == kProofAnalysis) ImportWrappers(source);
588 TIter next(fOutputs);
589 AliAnalysisDataContainer *output;
590 while ((output=(AliAnalysisDataContainer*)next())) {
591 if (!output->GetData()) continue;
592 // Check if there are client tasks that run post event loop
593 if (output->HasConsumers()) {
594 // Disable event loop semaphore
595 output->SetPostEventLoop(kTRUE);
596 TObjArray *list = output->GetConsumers();
597 Int_t ncons = list->GetEntriesFast();
598 for (Int_t i=0; i<ncons; i++) {
599 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
600 task->CheckNotify(kTRUE);
601 // If task is active, execute it
602 if (task->IsPostEventLoop() && task->IsActive()) {
603 if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
609 if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
612 //______________________________________________________________________________
613 void AliAnalysisManager::Terminate()
615 // The Terminate() function is the last function to be called during
616 // a query. It always runs on the client, it can be used to present
617 // the results graphically.
618 if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
619 AliAnalysisTask *task;
621 // Call Terminate() for tasks
622 while ((task=(AliAnalysisTask*)next())) task->Terminate();
624 TIter next1(fOutputs);
625 AliAnalysisDataContainer *output;
626 while ((output=(AliAnalysisDataContainer*)next1())) {
627 // Close all files at output
628 // Special outputs have the files already closed and written.
629 if (output->IsSpecialOutput()) continue;
630 const char *filename = output->GetFileName();
631 if (!(strcmp(filename, "default"))) {
632 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
633 TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
635 if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
640 if (!strlen(filename)) continue;
641 if (!output->GetData()) continue;
642 TFile *file = output->GetFile();
643 TDirectory *opwd = gDirectory;
647 file = new TFile(filename, "RECREATE");
648 if (file->IsZombie()) continue;
649 output->SetFile(file);
651 if (fDebug > 1) printf(" writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
652 if (output->GetData()->InheritsFrom(TCollection::Class())) {
653 // If data is a collection, we set the name of the collection
654 // as the one of the container and we save as a single key.
655 TCollection *coll = (TCollection*)output->GetData();
656 coll->SetName(output->GetName());
657 coll->Write(output->GetName(), TObject::kSingleKey);
659 output->GetData()->Write();
662 if (opwd) opwd->cd();
665 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
666 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
667 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
669 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
671 TDirectory *cdir = gDirectory;
672 TFile f("syswatch.root", "RECREATE");
674 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
675 tree->SetMarkerStyle(kCircle);
676 tree->SetMarkerColor(kBlue);
677 tree->SetMarkerSize(0.5);
678 if (!gROOT->IsBatch()) {
679 tree->SetAlias("event", "id0");
680 tree->SetAlias("memUSED", "pI.fMemVirtual");
681 tree->SetAlias("userCPU", "pI.fCpuUser");
682 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
683 c->Divide(2,1,0.01,0.01);
685 tree->Draw("memUSED:event","","", 1234567890, 0);
687 tree->Draw("userCPU:event","","", 1234567890, 0);
693 if (cdir) cdir->cd();
695 if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
698 //______________________________________________________________________________
699 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
701 // Adds a user task to the global list of tasks.
702 if (fTasks->FindObject(task)) {
703 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
706 task->SetActive(kFALSE);
710 //______________________________________________________________________________
711 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
713 // Retreive task by name.
714 if (!fTasks) return NULL;
715 return (AliAnalysisTask*)fTasks->FindObject(name);
718 //______________________________________________________________________________
719 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
720 TClass *datatype, EAliAnalysisContType type, const char *filename)
722 // Create a data container of a certain type. Types can be:
723 // kExchangeContainer = 0, used to exchange date between tasks
724 // kInputContainer = 1, used to store input data
725 // kOutputContainer = 2, used for posting results
726 if (fContainers->FindObject(name)) {
727 Error("CreateContainer","A container named %s already defined !\n",name);
730 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
731 fContainers->Add(cont);
733 case kInputContainer:
736 case kOutputContainer:
738 if (filename && strlen(filename)) {
739 cont->SetFileName(filename);
740 cont->SetDataOwned(kFALSE); // data owned by the file
743 case kExchangeContainer:
749 //______________________________________________________________________________
750 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
751 AliAnalysisDataContainer *cont)
753 // Connect input of an existing task to a data container.
754 if (!fTasks->FindObject(task)) {
756 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
758 Bool_t connected = task->ConnectInput(islot, cont);
762 //______________________________________________________________________________
763 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
764 AliAnalysisDataContainer *cont)
766 // Connect output of an existing task to a data container.
767 if (!fTasks->FindObject(task)) {
769 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
771 Bool_t connected = task->ConnectOutput(islot, cont);
775 //______________________________________________________________________________
776 void AliAnalysisManager::CleanContainers()
778 // Clean data from all containers that have already finished all client tasks.
779 TIter next(fContainers);
780 AliAnalysisDataContainer *cont;
781 while ((cont=(AliAnalysisDataContainer *)next())) {
782 if (cont->IsOwnedData() &&
783 cont->IsDataReady() &&
784 cont->ClientsExecuted()) cont->DeleteData();
788 //______________________________________________________________________________
789 Bool_t AliAnalysisManager::InitAnalysis()
791 // Initialization of analysis chain of tasks. Should be called after all tasks
792 // and data containers are properly connected
793 // Check for input/output containers
795 // Check for top tasks (depending only on input data containers)
796 if (!fTasks->First()) {
797 Error("InitAnalysis", "Analysis has no tasks !");
801 AliAnalysisTask *task;
802 AliAnalysisDataContainer *cont;
805 Bool_t iszombie = kFALSE;
806 Bool_t istop = kTRUE;
808 while ((task=(AliAnalysisTask*)next())) {
811 Int_t ninputs = task->GetNinputs();
812 for (i=0; i<ninputs; i++) {
813 cont = task->GetInputSlot(i)->GetContainer();
821 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
824 if (iszombie) continue;
825 // Check if cont is an input container
826 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
827 // Connect to parent task
831 fTopTasks->Add(task);
835 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
838 // Check now if there are orphan tasks
839 for (i=0; i<ntop; i++) {
840 task = (AliAnalysisTask*)fTopTasks->At(i);
845 while ((task=(AliAnalysisTask*)next())) {
846 if (!task->IsUsed()) {
848 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
851 // Check the task hierarchy (no parent task should depend on data provided
852 // by a daughter task)
853 for (i=0; i<ntop; i++) {
854 task = (AliAnalysisTask*)fTopTasks->At(i);
855 if (task->CheckCircularDeps()) {
856 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
861 // Check that all containers feeding post-event loop tasks are in the outputs list
862 TIter nextcont(fContainers); // loop over all containers
863 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
864 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
865 if (cont->HasConsumers()) {
866 // Check if one of the consumers is post event loop
867 TIter nextconsumer(cont->GetConsumers());
868 while ((task=(AliAnalysisTask*)nextconsumer())) {
869 if (task->IsPostEventLoop()) {
877 // Check if all special output containers have a file name provided
878 TIter nextout(fOutputs);
879 while ((cont=(AliAnalysisDataContainer*)nextout())) {
880 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
881 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
889 //______________________________________________________________________________
890 void AliAnalysisManager::PrintStatus(Option_t *option) const
892 // Print task hierarchy.
894 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
897 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
899 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
900 TIter next(fTopTasks);
901 AliAnalysisTask *task;
902 while ((task=(AliAnalysisTask*)next()))
903 task->PrintTask(option);
906 //______________________________________________________________________________
907 void AliAnalysisManager::ResetAnalysis()
909 // Reset all execution flags and clean containers.
913 //______________________________________________________________________________
914 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
916 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
917 // MIX. Process nentries starting from firstentry
919 Error("StartAnalysis","Analysis manager was not initialized !");
922 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
923 TString anaType = type;
925 fMode = kLocalAnalysis;
926 if (anaType.Contains("proof")) fMode = kProofAnalysis;
927 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
928 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
930 if (fMode == kGridAnalysis) {
931 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
932 fMode = kLocalAnalysis;
935 SetEventLoop(kFALSE);
936 // Enable event loop mode if a tree was provided
937 if (tree || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
940 TString ttype = "TTree";
941 if (tree && tree->IsA() == TChain::Class()) {
942 chain = (TChain*)tree;
943 if (!chain || !chain->GetListOfFiles()->First()) {
944 Error("StartAnalysis", "Cannot process null or empty chain...");
950 // Initialize locally all tasks (happens for all modes)
952 AliAnalysisTask *task;
953 while ((task=(AliAnalysisTask*)next())) {
961 // Call CreateOutputObjects for all tasks
962 while ((task=(AliAnalysisTask*)nextT())) {
963 TDirectory *curdir = gDirectory;
964 task->CreateOutputObjects();
965 if (curdir) curdir->cd();
971 // Run tree-based analysis via AliAnalysisSelector
972 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
973 fSelector = new AliAnalysisSelector(this);
974 tree->Process(fSelector, "", nentries, firstentry);
977 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
978 printf("StartAnalysis: no PROOF!!!\n");
981 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
982 gROOT->ProcessLine(line);
985 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
986 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
988 printf("StartAnalysis: no chain\n");
993 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
995 case kMixingAnalysis:
996 // Run event mixing analysis
998 Error("StartAnalysis", "Cannot run event mixing without event pool");
1001 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1002 fSelector = new AliAnalysisSelector(this);
1003 while ((chain=fEventPool->GetNextChain())) {
1005 // Call NotifyBinChange for all tasks
1006 while ((task=(AliAnalysisTask*)next()))
1007 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1008 chain->Process(fSelector);
1010 PackOutput(fSelector->GetOutputList());
1015 //______________________________________________________________________________
1016 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1018 // Start analysis for this manager on a given dataset. Analysis task can be:
1019 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1021 Error("StartAnalysis","Analysis manager was not initialized !");
1024 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
1025 TString anaType = type;
1027 if (!anaType.Contains("proof")) {
1028 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1031 fMode = kProofAnalysis;
1033 SetEventLoop(kTRUE);
1034 // Set the dataset flag
1035 TObject::SetBit(kUseDataSet);
1038 // Initialize locally all tasks
1040 AliAnalysisTask *task;
1041 while ((task=(AliAnalysisTask*)next())) {
1045 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1046 printf("StartAnalysis: no PROOF!!!\n");
1049 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1050 gROOT->ProcessLine(line);
1051 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
1052 if (!gROOT->ProcessLine(line)) {
1053 Error("StartAnalysis", "Dataset %s not found", dataset);
1056 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1057 dataset, nentries, firstentry);
1058 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1059 gROOT->ProcessLine(line);
1062 //______________________________________________________________________________
1063 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
1065 // Opens a special output file used in PROOF.
1067 if (fMode!=kProofAnalysis || !fSelector) {
1068 Error("OpenProofFile","Cannot open PROOF file %s",filename);
1071 sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
1072 if (fDebug > 1) printf("=== %s\n", line);
1073 gROOT->ProcessLine(line);
1074 sprintf(line, "pf->OpenFile(\"%s\");", option);
1075 gROOT->ProcessLine(line);
1077 gROOT->ProcessLine("pf->Print()");
1078 printf(" == proof file name: %s\n", gFile->GetName());
1080 sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1081 if (fDebug > 1) printf("=== %s\n", line);
1082 gROOT->ProcessLine(line);
1086 //______________________________________________________________________________
1087 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1089 // Execute analysis.
1090 static Long64_t ncalls = 0;
1091 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1092 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1095 Error("ExecAnalysis", "Analysis manager was not initialized !");
1098 AliAnalysisTask *task;
1099 // Check if the top tree is active.
1102 // De-activate all tasks
1103 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1104 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1106 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1109 cont->SetData(fTree); // This will notify all consumers
1110 Long64_t entry = fTree->GetTree()->GetReadEntry();
1113 // Call BeginEvent() for optional input/output and MC services
1114 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1115 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1116 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1118 // Execute the tasks
1119 // TIter next1(cont->GetConsumers());
1120 TIter next1(fTopTasks);
1121 while ((task=(AliAnalysisTask*)next1())) {
1123 cout << " Executing task " << task->GetName() << endl;
1126 task->ExecuteTask(option);
1129 // Call FinishEvent() for optional output and MC services
1130 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1131 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1132 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1133 // Gather system information if requested
1134 if (getsysInfo && ((ncalls%fNSysInfo)==0))
1135 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1138 // The event loop is not controlled by TSelector
1140 // Call BeginEvent() for optional input/output and MC services
1141 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
1142 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
1143 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1144 TIter next2(fTopTasks);
1145 while ((task=(AliAnalysisTask*)next2())) {
1146 task->SetActive(kTRUE);
1148 cout << " Executing task " << task->GetName() << endl;
1150 task->ExecuteTask(option);
1153 // Call FinishEvent() for optional output and MC services
1154 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1155 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1156 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1159 //______________________________________________________________________________
1160 void AliAnalysisManager::FinishAnalysis()