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 if (output->GetProducer()->IsPostEventLoop()) continue;
392 const char *filename = output->GetFileName();
393 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
394 isManagedByHandler = kTRUE;
395 filename = fOutputEventHandler->GetOutputFileName();
397 // Check if data was posted to this container. If not, issue an error.
398 if (!output->GetData() && !isManagedByHandler) {
399 Error("PackOutput", "No data for output container %s. Forgot to PostData ?\n", output->GetName());
402 if (!output->IsSpecialOutput()) {
404 if (strlen(filename) && !isManagedByHandler) {
405 // File resident outputs
406 TFile *file = output->GetFile();
407 // Backup current folder
408 TDirectory *opwd = gDirectory;
409 // Create file if not existing and register to container.
410 if (file) file->cd();
411 else file = new TFile(filename, "RECREATE");
412 if (file->IsZombie()) {
413 Fatal("PackOutput", "Could not recreate file %s\n", filename);
416 output->SetFile(file);
417 // Clear file list to release object ownership to user.
419 // Save data to file, then close.
420 if (output->GetData()->InheritsFrom(TCollection::Class())) {
421 // If data is a collection, we set the name of the collection
422 // as the one of the container and we save as a single key.
423 TCollection *coll = (TCollection*)output->GetData();
424 coll->SetName(output->GetName());
425 coll->Write(output->GetName(), TObject::kSingleKey);
427 output->GetData()->Write();
429 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
431 printf(" file %s listing content:\n", filename);
435 // Restore current directory
436 if (opwd) opwd->cd();
438 // Memory-resident outputs
439 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
441 AliAnalysisDataWrapper *wrap = 0;
442 if (isManagedByHandler) {
443 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
444 wrap->SetName(output->GetName());
446 else wrap =output->ExportData();
447 // Output wrappers must delete data after merging (AG 13/11/07)
448 wrap->SetDeleteData(kTRUE);
452 TDirectory *opwd = gDirectory;
453 TFile *file = output->GetFile();
454 if (isManagedByHandler) {
455 // Terminate IO for files managed by the output handler
456 if (file) file->Write();
457 fOutputEventHandler->TerminateIO();
462 AliAnalysisTask *producer = output->GetProducer();
464 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
465 output->GetFileName(), output->GetName(), producer->ClassName());
469 // Release object ownership to users after writing data to file
470 if (output->GetData()->InheritsFrom(TCollection::Class())) {
471 // If data is a collection, we set the name of the collection
472 // as the one of the container and we save as a single key.
473 TCollection *coll = (TCollection*)output->GetData();
474 coll->SetName(output->GetName());
475 coll->Write(output->GetName(), TObject::kSingleKey);
477 output->GetData()->Write();
481 printf(" file %s listing content:\n", output->GetFileName());
485 // Restore current directory
486 if (opwd) opwd->cd();
487 // Check if a special output location was provided or the output files have to be merged
488 if (strlen(fSpecialOutputLocation.Data())) {
489 TString remote = fSpecialOutputLocation;
491 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
492 remote += Form("%s_%d_", gSystem->HostName(), gid);
493 remote += output->GetFileName();
494 TFile::Cp(output->GetFileName(), remote.Data());
496 // No special location specified-> use TProofOutputFile as merging utility
497 // The file at this output slot must be opened in CreateOutputObjects
498 if (fDebug > 1) printf(" File %s to be merged...\n", output->GetFileName());
503 if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
506 //______________________________________________________________________________
507 void AliAnalysisManager::ImportWrappers(TList *source)
509 // Import data in output containers from wrappers coming in source.
510 if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
511 TIter next(fOutputs);
512 AliAnalysisDataContainer *cont;
513 AliAnalysisDataWrapper *wrap;
515 while ((cont=(AliAnalysisDataContainer*)next())) {
517 if (cont->GetProducer()->IsPostEventLoop()) continue;
518 const char *filename = cont->GetFileName();
519 Bool_t isManagedByHandler = kFALSE;
520 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
521 isManagedByHandler = kTRUE;
522 filename = fOutputEventHandler->GetOutputFileName();
524 if (cont->IsSpecialOutput()) {
525 if (strlen(fSpecialOutputLocation.Data()) && !isManagedByHandler) continue;
526 // Copy merged file from PROOF scratch space
528 TObject *pof = source->FindObject(filename);
529 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
530 Error("ImportWrappers", "TProofOutputFile object not found in output list for container %s", cont->GetName());
533 gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", ((TProofOutputFile*)0x%lx)->GetOutputFileName();)", full_path, pof));
535 printf(" Copying file %s from PROOF scratch space\n", full_path);
536 Bool_t gotit = TFile::Cp(full_path, filename);
538 Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
540 // Normally we should connect data from the copied file to the
541 // corresponding output container, but it is not obvious how to do this
542 // automatically if several objects in file...
543 TFile *f = new TFile(filename, "READ");
544 TObject *obj = f->Get(cont->GetName());
546 Error("ImportWrappers", "Could not find object %s in file %s", cont->GetName(), filename);
549 wrap = new AliAnalysisDataWrapper(obj);
550 wrap->SetDeleteData(kFALSE);
552 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
554 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
559 printf(" Importing data for container %s", cont->GetName());
560 if (strlen(filename)) printf(" -> file %s\n", cont->GetFileName());
563 cont->ImportData(wrap);
565 if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
568 //______________________________________________________________________________
569 void AliAnalysisManager::UnpackOutput(TList *source)
571 // Called by AliAnalysisSelector::Terminate only on the client.
572 if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
574 Error("UnpackOutput", "No target. Aborting.");
577 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
579 if (fMode == kProofAnalysis) ImportWrappers(source);
581 TIter next(fOutputs);
582 AliAnalysisDataContainer *output;
583 while ((output=(AliAnalysisDataContainer*)next())) {
584 if (!output->GetData()) continue;
585 // Check if there are client tasks that run post event loop
586 if (output->HasConsumers()) {
587 // Disable event loop semaphore
588 output->SetPostEventLoop(kTRUE);
589 TObjArray *list = output->GetConsumers();
590 Int_t ncons = list->GetEntriesFast();
591 for (Int_t i=0; i<ncons; i++) {
592 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
593 task->CheckNotify(kTRUE);
594 // If task is active, execute it
595 if (task->IsPostEventLoop() && task->IsActive()) {
596 if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
602 if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
605 //______________________________________________________________________________
606 void AliAnalysisManager::Terminate()
608 // The Terminate() function is the last function to be called during
609 // a query. It always runs on the client, it can be used to present
610 // the results graphically.
611 if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
612 AliAnalysisTask *task;
614 // Call Terminate() for tasks
615 while ((task=(AliAnalysisTask*)next())) task->Terminate();
617 TIter next1(fOutputs);
618 AliAnalysisDataContainer *output;
619 while ((output=(AliAnalysisDataContainer*)next1())) {
620 // Close all files at output
621 // Special outputs have the files already closed and written.
622 if (output->IsSpecialOutput()) continue;
623 const char *filename = output->GetFileName();
624 if (!(strcmp(filename, "default"))) {
625 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
626 TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
628 if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
633 if (!strlen(filename)) continue;
634 if (!output->GetData()) continue;
635 TFile *file = output->GetFile();
636 TDirectory *opwd = gDirectory;
640 file = new TFile(filename, "RECREATE");
641 if (file->IsZombie()) continue;
642 output->SetFile(file);
644 if (fDebug > 1) printf(" writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
645 if (output->GetData()->InheritsFrom(TCollection::Class())) {
646 // If data is a collection, we set the name of the collection
647 // as the one of the container and we save as a single key.
648 TCollection *coll = (TCollection*)output->GetData();
649 coll->SetName(output->GetName());
650 coll->Write(output->GetName(), TObject::kSingleKey);
652 output->GetData()->Write();
655 if (opwd) opwd->cd();
658 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
659 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
660 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
662 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
664 TDirectory *cdir = gDirectory;
665 TFile f("syswatch.root", "RECREATE");
667 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
668 tree->SetMarkerStyle(kCircle);
669 tree->SetMarkerColor(kBlue);
670 tree->SetMarkerSize(0.5);
671 if (!gROOT->IsBatch()) {
672 tree->SetAlias("event", "id0");
673 tree->SetAlias("memUSED", "pI.fMemVirtual");
674 tree->SetAlias("userCPU", "pI.fCpuUser");
675 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
676 c->Divide(2,1,0.01,0.01);
678 tree->Draw("memUSED:event","","", 1234567890, 0);
680 tree->Draw("userCPU:event","","", 1234567890, 0);
686 if (cdir) cdir->cd();
688 if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
691 //______________________________________________________________________________
692 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
694 // Adds a user task to the global list of tasks.
695 if (fTasks->FindObject(task)) {
696 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
699 task->SetActive(kFALSE);
703 //______________________________________________________________________________
704 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
706 // Retreive task by name.
707 if (!fTasks) return NULL;
708 return (AliAnalysisTask*)fTasks->FindObject(name);
711 //______________________________________________________________________________
712 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
713 TClass *datatype, EAliAnalysisContType type, const char *filename)
715 // Create a data container of a certain type. Types can be:
716 // kExchangeContainer = 0, used to exchange date between tasks
717 // kInputContainer = 1, used to store input data
718 // kOutputContainer = 2, used for posting results
719 if (fContainers->FindObject(name)) {
720 Error("CreateContainer","A container named %s already defined !\n",name);
723 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
724 fContainers->Add(cont);
726 case kInputContainer:
729 case kOutputContainer:
731 if (filename && strlen(filename)) {
732 cont->SetFileName(filename);
733 cont->SetDataOwned(kFALSE); // data owned by the file
736 case kExchangeContainer:
742 //______________________________________________________________________________
743 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
744 AliAnalysisDataContainer *cont)
746 // Connect input of an existing task to a data container.
747 if (!fTasks->FindObject(task)) {
749 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
751 Bool_t connected = task->ConnectInput(islot, cont);
755 //______________________________________________________________________________
756 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
757 AliAnalysisDataContainer *cont)
759 // Connect output of an existing task to a data container.
760 if (!fTasks->FindObject(task)) {
762 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
764 Bool_t connected = task->ConnectOutput(islot, cont);
768 //______________________________________________________________________________
769 void AliAnalysisManager::CleanContainers()
771 // Clean data from all containers that have already finished all client tasks.
772 TIter next(fContainers);
773 AliAnalysisDataContainer *cont;
774 while ((cont=(AliAnalysisDataContainer *)next())) {
775 if (cont->IsOwnedData() &&
776 cont->IsDataReady() &&
777 cont->ClientsExecuted()) cont->DeleteData();
781 //______________________________________________________________________________
782 Bool_t AliAnalysisManager::InitAnalysis()
784 // Initialization of analysis chain of tasks. Should be called after all tasks
785 // and data containers are properly connected
786 // Check for input/output containers
788 // Check for top tasks (depending only on input data containers)
789 if (!fTasks->First()) {
790 Error("InitAnalysis", "Analysis has no tasks !");
794 AliAnalysisTask *task;
795 AliAnalysisDataContainer *cont;
798 Bool_t iszombie = kFALSE;
799 Bool_t istop = kTRUE;
801 while ((task=(AliAnalysisTask*)next())) {
804 Int_t ninputs = task->GetNinputs();
805 for (i=0; i<ninputs; i++) {
806 cont = task->GetInputSlot(i)->GetContainer();
814 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
817 if (iszombie) continue;
818 // Check if cont is an input container
819 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
820 // Connect to parent task
824 fTopTasks->Add(task);
828 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
831 // Check now if there are orphan tasks
832 for (i=0; i<ntop; i++) {
833 task = (AliAnalysisTask*)fTopTasks->At(i);
838 while ((task=(AliAnalysisTask*)next())) {
839 if (!task->IsUsed()) {
841 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
844 // Check the task hierarchy (no parent task should depend on data provided
845 // by a daughter task)
846 for (i=0; i<ntop; i++) {
847 task = (AliAnalysisTask*)fTopTasks->At(i);
848 if (task->CheckCircularDeps()) {
849 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
854 // Check that all containers feeding post-event loop tasks are in the outputs list
855 TIter nextcont(fContainers); // loop over all containers
856 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
857 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
858 if (cont->HasConsumers()) {
859 // Check if one of the consumers is post event loop
860 TIter nextconsumer(cont->GetConsumers());
861 while ((task=(AliAnalysisTask*)nextconsumer())) {
862 if (task->IsPostEventLoop()) {
870 // Check if all special output containers have a file name provided
871 TIter nextout(fOutputs);
872 while ((cont=(AliAnalysisDataContainer*)nextout())) {
873 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
874 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
882 //______________________________________________________________________________
883 void AliAnalysisManager::PrintStatus(Option_t *option) const
885 // Print task hierarchy.
887 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
890 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
892 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
893 TIter next(fTopTasks);
894 AliAnalysisTask *task;
895 while ((task=(AliAnalysisTask*)next()))
896 task->PrintTask(option);
899 //______________________________________________________________________________
900 void AliAnalysisManager::ResetAnalysis()
902 // Reset all execution flags and clean containers.
906 //______________________________________________________________________________
907 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
909 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
910 // MIX. Process nentries starting from firstentry
912 Error("StartAnalysis","Analysis manager was not initialized !");
915 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
916 TString anaType = type;
918 fMode = kLocalAnalysis;
919 if (anaType.Contains("proof")) fMode = kProofAnalysis;
920 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
921 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
923 if (fMode == kGridAnalysis) {
924 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
925 fMode = kLocalAnalysis;
928 SetEventLoop(kFALSE);
929 // Enable event loop mode if a tree was provided
930 if (tree || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
933 TString ttype = "TTree";
934 if (tree && tree->IsA() == TChain::Class()) {
935 chain = (TChain*)tree;
936 if (!chain || !chain->GetListOfFiles()->First()) {
937 Error("StartAnalysis", "Cannot process null or empty chain...");
943 // Initialize locally all tasks (happens for all modes)
945 AliAnalysisTask *task;
946 while ((task=(AliAnalysisTask*)next())) {
954 // Call CreateOutputObjects for all tasks
955 while ((task=(AliAnalysisTask*)nextT())) {
956 TDirectory *curdir = gDirectory;
957 task->CreateOutputObjects();
958 if (curdir) curdir->cd();
964 // Run tree-based analysis via AliAnalysisSelector
965 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
966 fSelector = new AliAnalysisSelector(this);
967 tree->Process(fSelector, "", nentries, firstentry);
970 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
971 printf("StartAnalysis: no PROOF!!!\n");
974 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
975 gROOT->ProcessLine(line);
978 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
979 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
981 printf("StartAnalysis: no chain\n");
986 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
988 case kMixingAnalysis:
989 // Run event mixing analysis
991 Error("StartAnalysis", "Cannot run event mixing without event pool");
994 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
995 fSelector = new AliAnalysisSelector(this);
996 while ((chain=fEventPool->GetNextChain())) {
998 // Call NotifyBinChange for all tasks
999 while ((task=(AliAnalysisTask*)next()))
1000 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1001 chain->Process(fSelector);
1003 PackOutput(fSelector->GetOutputList());
1008 //______________________________________________________________________________
1009 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1011 // Start analysis for this manager on a given dataset. Analysis task can be:
1012 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1014 Error("StartAnalysis","Analysis manager was not initialized !");
1017 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
1018 TString anaType = type;
1020 if (!anaType.Contains("proof")) {
1021 Error("Cannot process datasets in %s mode. Try PROOF.", type);
1024 fMode = kProofAnalysis;
1026 SetEventLoop(kTRUE);
1027 // Set the dataset flag
1028 TObject::SetBit(kUseDataSet);
1031 // Initialize locally all tasks
1033 AliAnalysisTask *task;
1034 while ((task=(AliAnalysisTask*)next())) {
1038 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1039 printf("StartAnalysis: no PROOF!!!\n");
1042 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1043 gROOT->ProcessLine(line);
1044 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
1045 if (!gROOT->ProcessLine(line)) {
1046 Error("StartAnalysis", "Dataset %s not found", dataset);
1049 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1050 dataset, nentries, firstentry);
1051 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1052 gROOT->ProcessLine(line);
1055 //______________________________________________________________________________
1056 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
1058 // Opens a special output file used in PROOF.
1060 if (fMode!=kProofAnalysis || !fSelector) {
1061 Error("OpenProofFile","Cannot open PROOF file %s",filename);
1064 sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
1065 if (fDebug > 1) printf("=== %s\n", line);
1066 gROOT->ProcessLine(line);
1067 sprintf(line, "pf->OpenFile(\"%s\");", option);
1068 gROOT->ProcessLine(line);
1070 gROOT->ProcessLine("pf->Print()");
1071 printf(" == proof file name: %s\n", gFile->GetName());
1073 sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1074 if (fDebug > 1) printf("=== %s\n", line);
1075 gROOT->ProcessLine(line);
1079 //______________________________________________________________________________
1080 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1082 // Execute analysis.
1083 static Long64_t ncalls = 0;
1084 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1085 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1088 Error("ExecAnalysis", "Analysis manager was not initialized !");
1091 AliAnalysisTask *task;
1092 // Check if the top tree is active.
1095 // De-activate all tasks
1096 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1097 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1099 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1102 cont->SetData(fTree); // This will notify all consumers
1103 Long64_t entry = fTree->GetTree()->GetReadEntry();
1106 // Call BeginEvent() for optional input/output and MC services
1107 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1108 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1109 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1111 // Execute the tasks
1112 // TIter next1(cont->GetConsumers());
1113 TIter next1(fTopTasks);
1114 while ((task=(AliAnalysisTask*)next1())) {
1116 cout << " Executing task " << task->GetName() << endl;
1119 task->ExecuteTask(option);
1122 // Call FinishEvent() for optional output and MC services
1123 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1124 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1125 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1126 // Gather system information if requested
1127 if (getsysInfo && ((ncalls%fNSysInfo)==0))
1128 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1131 // The event loop is not controlled by TSelector
1133 // Call BeginEvent() for optional input/output and MC services
1134 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
1135 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
1136 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1137 TIter next2(fTopTasks);
1138 while ((task=(AliAnalysisTask*)next2())) {
1139 task->SetActive(kTRUE);
1141 cout << " Executing task " << task->GetName() << endl;
1143 task->ExecuteTask(option);
1146 // Call FinishEvent() for optional output and MC services
1147 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1148 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1149 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1152 //______________________________________________________________________________
1153 void AliAnalysisManager::FinishAnalysis()