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());
502 TString outFilename = file->GetName();
504 // Restore current directory
505 if (opwd) opwd->cd();
506 // Check if a special output location was provided or the output files have to be merged
507 if (strlen(fSpecialOutputLocation.Data())) {
508 TString remote = fSpecialOutputLocation;
510 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
511 remote += Form("%s_%d_", gSystem->HostName(), gid);
512 remote += output->GetFileName();
513 TFile::Cp ( outFilename.Data(), remote.Data() );
515 // No special location specified-> use TProofOutputFile as merging utility
516 // The file at this output slot must be opened in CreateOutputObjects
517 if (fDebug > 1) printf(" File %s to be merged...\n", output->GetFileName());
522 if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
525 //______________________________________________________________________________
526 void AliAnalysisManager::ImportWrappers(TList *source)
528 // Import data in output containers from wrappers coming in source.
529 if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
530 TIter next(fOutputs);
531 AliAnalysisDataContainer *cont;
532 AliAnalysisDataWrapper *wrap;
534 while ((cont=(AliAnalysisDataContainer*)next())) {
536 if (cont->GetProducer()->IsPostEventLoop()) continue;
537 const char *filename = cont->GetFileName();
538 Bool_t isManagedByHandler = kFALSE;
539 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
540 isManagedByHandler = kTRUE;
541 filename = fOutputEventHandler->GetOutputFileName();
543 if (cont->IsSpecialOutput()) {
544 if (strlen(fSpecialOutputLocation.Data()) && !isManagedByHandler) continue;
545 // Copy merged file from PROOF scratch space
547 TObject *pof = source->FindObject(filename);
548 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
549 Error("ImportWrappers", "TProofOutputFile object not found in output list for container %s", cont->GetName());
552 gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", ((TProofOutputFile*)0x%lx)->GetOutputFileName();)", full_path, pof));
554 printf(" Copying file %s from PROOF scratch space\n", full_path);
555 Bool_t gotit = TFile::Cp(full_path, filename);
557 Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
559 // Normally we should connect data from the copied file to the
560 // corresponding output container, but it is not obvious how to do this
561 // automatically if several objects in file...
562 TFile *f = new TFile(filename, "READ");
564 if (!isManagedByHandler) obj = f->Get(cont->GetName());
565 if (!obj && !isManagedByHandler) {
566 Error("ImportWrappers", "Could not find object %s in file %s", cont->GetName(), filename);
569 wrap = new AliAnalysisDataWrapper(obj);
570 wrap->SetDeleteData(kFALSE);
572 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
574 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
579 printf(" Importing data for container %s", cont->GetName());
580 if (strlen(filename)) printf(" -> file %s\n", cont->GetFileName());
583 cont->ImportData(wrap);
585 if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
588 //______________________________________________________________________________
589 void AliAnalysisManager::UnpackOutput(TList *source)
591 // Called by AliAnalysisSelector::Terminate only on the client.
592 if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
594 Error("UnpackOutput", "No target. Aborting.");
597 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
599 if (fMode == kProofAnalysis) ImportWrappers(source);
601 TIter next(fOutputs);
602 AliAnalysisDataContainer *output;
603 while ((output=(AliAnalysisDataContainer*)next())) {
604 if (!output->GetData()) continue;
605 // Check if there are client tasks that run post event loop
606 if (output->HasConsumers()) {
607 // Disable event loop semaphore
608 output->SetPostEventLoop(kTRUE);
609 TObjArray *list = output->GetConsumers();
610 Int_t ncons = list->GetEntriesFast();
611 for (Int_t i=0; i<ncons; i++) {
612 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
613 task->CheckNotify(kTRUE);
614 // If task is active, execute it
615 if (task->IsPostEventLoop() && task->IsActive()) {
616 if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
622 if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
625 //______________________________________________________________________________
626 void AliAnalysisManager::Terminate()
628 // The Terminate() function is the last function to be called during
629 // a query. It always runs on the client, it can be used to present
630 // the results graphically.
631 if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
632 AliAnalysisTask *task;
634 // Call Terminate() for tasks
635 while ((task=(AliAnalysisTask*)next())) task->Terminate();
637 TIter next1(fOutputs);
638 AliAnalysisDataContainer *output;
639 while ((output=(AliAnalysisDataContainer*)next1())) {
640 // Special outputs have the files already closed and written.
641 if (output->IsSpecialOutput()&&(fMode == kProofAnalysis)) continue;
642 const char *filename = output->GetFileName();
643 if (!(strcmp(filename, "default"))) {
644 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
645 TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
647 if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
652 if (!strlen(filename)) continue;
653 if (!output->GetData()) continue;
654 TFile *file = output->GetFile();
655 TDirectory *opwd = gDirectory;
656 file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
657 if (!file) file = new TFile(filename, "RECREATE");
658 if (file->IsZombie()) continue;
659 output->SetFile(file);
661 if (fDebug > 1) printf(" writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
662 if (output->GetData()->InheritsFrom(TCollection::Class())) {
663 // If data is a collection, we set the name of the collection
664 // as the one of the container and we save as a single key.
665 TCollection *coll = (TCollection*)output->GetData();
666 coll->SetName(output->GetName());
667 coll->Write(output->GetName(), TObject::kSingleKey);
669 if (output->GetData()->InheritsFrom(TTree::Class())) {
670 TTree *tree = (TTree*)output->GetData();
671 tree->SetDirectory(file);
674 output->GetData()->Write();
677 if (opwd) opwd->cd();
680 while ((output=(AliAnalysisDataContainer*)next1())) {
681 // Close all files at output
682 TDirectory *opwd = gDirectory;
683 if (output->GetFile()) output->GetFile()->Close();
684 if (opwd) opwd->cd();
687 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
688 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
689 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
691 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
693 TDirectory *cdir = gDirectory;
694 TFile f("syswatch.root", "RECREATE");
696 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
697 tree->SetMarkerStyle(kCircle);
698 tree->SetMarkerColor(kBlue);
699 tree->SetMarkerSize(0.5);
700 if (!gROOT->IsBatch()) {
701 tree->SetAlias("event", "id0");
702 tree->SetAlias("memUSED", "pI.fMemVirtual");
703 tree->SetAlias("userCPU", "pI.fCpuUser");
704 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
705 c->Divide(2,1,0.01,0.01);
707 tree->Draw("memUSED:event","","", 1234567890, 0);
709 tree->Draw("userCPU:event","","", 1234567890, 0);
715 if (cdir) cdir->cd();
717 if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
720 //______________________________________________________________________________
721 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
723 // Adds a user task to the global list of tasks.
724 if (fTasks->FindObject(task)) {
725 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
728 task->SetActive(kFALSE);
732 //______________________________________________________________________________
733 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
735 // Retreive task by name.
736 if (!fTasks) return NULL;
737 return (AliAnalysisTask*)fTasks->FindObject(name);
740 //______________________________________________________________________________
741 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
742 TClass *datatype, EAliAnalysisContType type, const char *filename)
744 // Create a data container of a certain type. Types can be:
745 // kExchangeContainer = 0, used to exchange date between tasks
746 // kInputContainer = 1, used to store input data
747 // kOutputContainer = 2, used for posting results
748 if (fContainers->FindObject(name)) {
749 Error("CreateContainer","A container named %s already defined !\n",name);
752 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
753 fContainers->Add(cont);
755 case kInputContainer:
758 case kOutputContainer:
760 if (filename && strlen(filename)) {
761 cont->SetFileName(filename);
762 cont->SetDataOwned(kFALSE); // data owned by the file
765 case kExchangeContainer:
771 //______________________________________________________________________________
772 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
773 AliAnalysisDataContainer *cont)
775 // Connect input of an existing task to a data container.
776 if (!fTasks->FindObject(task)) {
778 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
780 Bool_t connected = task->ConnectInput(islot, cont);
784 //______________________________________________________________________________
785 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
786 AliAnalysisDataContainer *cont)
788 // Connect output of an existing task to a data container.
789 if (!fTasks->FindObject(task)) {
791 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
793 Bool_t connected = task->ConnectOutput(islot, cont);
797 //______________________________________________________________________________
798 void AliAnalysisManager::CleanContainers()
800 // Clean data from all containers that have already finished all client tasks.
801 TIter next(fContainers);
802 AliAnalysisDataContainer *cont;
803 while ((cont=(AliAnalysisDataContainer *)next())) {
804 if (cont->IsOwnedData() &&
805 cont->IsDataReady() &&
806 cont->ClientsExecuted()) cont->DeleteData();
810 //______________________________________________________________________________
811 Bool_t AliAnalysisManager::InitAnalysis()
813 // Initialization of analysis chain of tasks. Should be called after all tasks
814 // and data containers are properly connected
815 // Check for input/output containers
817 // Check for top tasks (depending only on input data containers)
818 if (!fTasks->First()) {
819 Error("InitAnalysis", "Analysis has no tasks !");
823 AliAnalysisTask *task;
824 AliAnalysisDataContainer *cont;
827 Bool_t iszombie = kFALSE;
828 Bool_t istop = kTRUE;
830 while ((task=(AliAnalysisTask*)next())) {
833 Int_t ninputs = task->GetNinputs();
834 for (i=0; i<ninputs; i++) {
835 cont = task->GetInputSlot(i)->GetContainer();
843 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
846 if (iszombie) continue;
847 // Check if cont is an input container
848 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
849 // Connect to parent task
853 fTopTasks->Add(task);
857 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
860 // Check now if there are orphan tasks
861 for (i=0; i<ntop; i++) {
862 task = (AliAnalysisTask*)fTopTasks->At(i);
867 while ((task=(AliAnalysisTask*)next())) {
868 if (!task->IsUsed()) {
870 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
873 // Check the task hierarchy (no parent task should depend on data provided
874 // by a daughter task)
875 for (i=0; i<ntop; i++) {
876 task = (AliAnalysisTask*)fTopTasks->At(i);
877 if (task->CheckCircularDeps()) {
878 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
883 // Check that all containers feeding post-event loop tasks are in the outputs list
884 TIter nextcont(fContainers); // loop over all containers
885 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
886 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
887 if (cont->HasConsumers()) {
888 // Check if one of the consumers is post event loop
889 TIter nextconsumer(cont->GetConsumers());
890 while ((task=(AliAnalysisTask*)nextconsumer())) {
891 if (task->IsPostEventLoop()) {
899 // Check if all special output containers have a file name provided
900 TIter nextout(fOutputs);
901 while ((cont=(AliAnalysisDataContainer*)nextout())) {
902 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
903 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
911 //______________________________________________________________________________
912 void AliAnalysisManager::PrintStatus(Option_t *option) const
914 // Print task hierarchy.
916 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
919 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
921 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
922 TIter next(fTopTasks);
923 AliAnalysisTask *task;
924 while ((task=(AliAnalysisTask*)next()))
925 task->PrintTask(option);
928 //______________________________________________________________________________
929 void AliAnalysisManager::ResetAnalysis()
931 // Reset all execution flags and clean containers.
935 //______________________________________________________________________________
936 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
938 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
939 // MIX. Process nentries starting from firstentry
941 Error("StartAnalysis","Analysis manager was not initialized !");
944 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
945 TString anaType = type;
947 fMode = kLocalAnalysis;
948 if (anaType.Contains("proof")) fMode = kProofAnalysis;
949 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
950 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
952 if (fMode == kGridAnalysis) {
953 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
954 fMode = kLocalAnalysis;
957 SetEventLoop(kFALSE);
958 // Enable event loop mode if a tree was provided
959 if (tree || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
962 TString ttype = "TTree";
963 if (tree && tree->IsA() == TChain::Class()) {
964 chain = (TChain*)tree;
965 if (!chain || !chain->GetListOfFiles()->First()) {
966 Error("StartAnalysis", "Cannot process null or empty chain...");
972 // Initialize locally all tasks (happens for all modes)
974 AliAnalysisTask *task;
975 while ((task=(AliAnalysisTask*)next())) {
983 // Call CreateOutputObjects for all tasks
984 while ((task=(AliAnalysisTask*)nextT())) {
985 TDirectory *curdir = gDirectory;
986 task->CreateOutputObjects();
987 if (curdir) curdir->cd();
993 // Run tree-based analysis via AliAnalysisSelector
994 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
995 fSelector = new AliAnalysisSelector(this);
996 tree->Process(fSelector, "", nentries, firstentry);
999 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1000 printf("StartAnalysis: no PROOF!!!\n");
1003 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1004 gROOT->ProcessLine(line);
1007 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1008 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1010 printf("StartAnalysis: no chain\n");
1015 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
1017 case kMixingAnalysis:
1018 // Run event mixing analysis
1020 Error("StartAnalysis", "Cannot run event mixing without event pool");
1023 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1024 fSelector = new AliAnalysisSelector(this);
1025 while ((chain=fEventPool->GetNextChain())) {
1027 // Call NotifyBinChange for all tasks
1028 while ((task=(AliAnalysisTask*)next()))
1029 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1030 chain->Process(fSelector);
1032 PackOutput(fSelector->GetOutputList());
1037 //______________________________________________________________________________
1038 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1040 // Start analysis for this manager on a given dataset. Analysis task can be:
1041 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1043 Error("StartAnalysis","Analysis manager was not initialized !");
1046 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
1047 TString anaType = type;
1049 if (!anaType.Contains("proof")) {
1050 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1053 fMode = kProofAnalysis;
1055 SetEventLoop(kTRUE);
1056 // Set the dataset flag
1057 TObject::SetBit(kUseDataSet);
1060 // Initialize locally all tasks
1062 AliAnalysisTask *task;
1063 while ((task=(AliAnalysisTask*)next())) {
1067 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1068 printf("StartAnalysis: no PROOF!!!\n");
1071 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1072 gROOT->ProcessLine(line);
1073 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
1074 if (!gROOT->ProcessLine(line)) {
1075 Error("StartAnalysis", "Dataset %s not found", dataset);
1078 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1079 dataset, nentries, firstentry);
1080 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1081 gROOT->ProcessLine(line);
1084 //______________________________________________________________________________
1085 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
1087 // Opens a special output file used in PROOF.
1089 if (fMode!=kProofAnalysis || !fSelector) {
1090 Error("OpenProofFile","Cannot open PROOF file %s",filename);
1093 sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
1094 if (fDebug > 1) printf("=== %s\n", line);
1095 gROOT->ProcessLine(line);
1096 sprintf(line, "pf->OpenFile(\"%s\");", option);
1097 gROOT->ProcessLine(line);
1099 gROOT->ProcessLine("pf->Print()");
1100 printf(" == proof file name: %s\n", gFile->GetName());
1102 sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1103 if (fDebug > 1) printf("=== %s\n", line);
1104 gROOT->ProcessLine(line);
1108 //______________________________________________________________________________
1109 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1111 // Execute analysis.
1112 static Long64_t ncalls = 0;
1113 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1114 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1117 Error("ExecAnalysis", "Analysis manager was not initialized !");
1120 AliAnalysisTask *task;
1121 // Check if the top tree is active.
1124 // De-activate all tasks
1125 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1126 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1128 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1131 cont->SetData(fTree); // This will notify all consumers
1132 Long64_t entry = fTree->GetTree()->GetReadEntry();
1135 // Call BeginEvent() for optional input/output and MC services
1136 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1137 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1138 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1140 // Execute the tasks
1141 // TIter next1(cont->GetConsumers());
1142 TIter next1(fTopTasks);
1143 while ((task=(AliAnalysisTask*)next1())) {
1145 cout << " Executing task " << task->GetName() << endl;
1148 task->ExecuteTask(option);
1151 // Call FinishEvent() for optional output and MC services
1152 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1153 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1154 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1155 // Gather system information if requested
1156 if (getsysInfo && ((ncalls%fNSysInfo)==0))
1157 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1160 // The event loop is not controlled by TSelector
1162 // Call BeginEvent() for optional input/output and MC services
1163 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
1164 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
1165 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1166 TIter next2(fTopTasks);
1167 while ((task=(AliAnalysisTask*)next2())) {
1168 task->SetActive(kTRUE);
1170 cout << " Executing task " << task->GetName() << endl;
1172 task->ExecuteTask(option);
1175 // Call FinishEvent() for optional output and MC services
1176 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1177 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1178 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1181 //______________________________________________________________________________
1182 void AliAnalysisManager::FinishAnalysis()