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>
33 #include <TMethodCall.h>
39 #include "AliAnalysisSelector.h"
40 #include "AliAnalysisGrid.h"
41 #include "AliAnalysisTask.h"
42 #include "AliAnalysisDataContainer.h"
43 #include "AliAnalysisDataSlot.h"
44 #include "AliVEventHandler.h"
45 #include "AliVEventPool.h"
46 #include "AliSysInfo.h"
47 #include "AliAnalysisManager.h"
49 ClassImp(AliAnalysisManager)
51 AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL;
53 //______________________________________________________________________________
54 AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
57 fInputEventHandler(NULL),
58 fOutputEventHandler(NULL),
59 fMCtruthEventHandler(NULL),
63 fMode(kLocalAnalysis),
66 fSpecialOutputLocation(""),
76 // Default constructor.
77 fgAnalysisManager = this;
78 fTasks = new TObjArray();
79 fTopTasks = new TObjArray();
80 fZombies = new TObjArray();
81 fContainers = new TObjArray();
82 fInputs = new TObjArray();
83 fOutputs = new TObjArray();
87 //______________________________________________________________________________
88 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
91 fInputEventHandler(NULL),
92 fOutputEventHandler(NULL),
93 fMCtruthEventHandler(NULL),
98 fInitOK(other.fInitOK),
100 fSpecialOutputLocation(""),
111 fTasks = new TObjArray(*other.fTasks);
112 fTopTasks = new TObjArray(*other.fTopTasks);
113 fZombies = new TObjArray(*other.fZombies);
114 fContainers = new TObjArray(*other.fContainers);
115 fInputs = new TObjArray(*other.fInputs);
116 fOutputs = new TObjArray(*other.fOutputs);
117 fgAnalysisManager = this;
120 //______________________________________________________________________________
121 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
124 if (&other != this) {
125 TNamed::operator=(other);
126 fInputEventHandler = other.fInputEventHandler;
127 fOutputEventHandler = other.fOutputEventHandler;
128 fMCtruthEventHandler = other.fMCtruthEventHandler;
129 fEventPool = other.fEventPool;
132 fNSysInfo = other.fNSysInfo;
134 fInitOK = other.fInitOK;
135 fDebug = other.fDebug;
136 fTasks = new TObjArray(*other.fTasks);
137 fTopTasks = new TObjArray(*other.fTopTasks);
138 fZombies = new TObjArray(*other.fZombies);
139 fContainers = new TObjArray(*other.fContainers);
140 fInputs = new TObjArray(*other.fInputs);
141 fOutputs = new TObjArray(*other.fOutputs);
144 fgAnalysisManager = this;
149 //______________________________________________________________________________
150 AliAnalysisManager::~AliAnalysisManager()
153 if (fTasks) {fTasks->Delete(); delete fTasks;}
154 if (fTopTasks) delete fTopTasks;
155 if (fZombies) delete fZombies;
156 if (fContainers) {fContainers->Delete(); delete fContainers;}
157 if (fInputs) delete fInputs;
158 if (fOutputs) delete fOutputs;
159 if (fGridHandler) delete fGridHandler;
160 if (fgAnalysisManager==this) fgAnalysisManager = NULL;
163 //______________________________________________________________________________
164 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
166 // Read one entry of the tree or a whole branch.
167 if (fDebug > 0) printf("== AliAnalysisManager::GetEntry(%lld)\n", entry);
168 fCurrentEntry = entry;
169 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
172 //______________________________________________________________________________
173 void AliAnalysisManager::Init(TTree *tree)
175 // The Init() function is called when the selector needs to initialize
176 // a new tree or chain. Typically here the branch addresses of the tree
177 // will be set. It is normaly not necessary to make changes to the
178 // generated code, but the routine can be extended by the user if needed.
179 // Init() will be called many times when running with PROOF.
182 printf("->AliAnalysisManager::Init(%s)\n", tree->GetName());
185 // Call InitTree of EventHandler
186 if (fOutputEventHandler) {
187 if (fMode == kProofAnalysis) {
188 fOutputEventHandler->Init(0x0, "proof");
190 fOutputEventHandler->Init(0x0, "local");
194 if (fInputEventHandler) {
195 if (fMode == kProofAnalysis) {
196 fInputEventHandler->Init(tree, "proof");
198 fInputEventHandler->Init(tree, "local");
201 // If no input event handler we need to get the tree once
203 if(!tree->GetTree()) tree->LoadTree(0);
207 if (fMCtruthEventHandler) {
208 if (fMode == kProofAnalysis) {
209 fMCtruthEventHandler->Init(0x0, "proof");
211 fMCtruthEventHandler->Init(0x0, "local");
215 if (!fInitOK) InitAnalysis();
216 if (!fInitOK) return;
218 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
220 Error("Init","No top input container !");
225 printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName());
229 //______________________________________________________________________________
230 void AliAnalysisManager::SlaveBegin(TTree *tree)
232 // The SlaveBegin() function is called after the Begin() function.
233 // When running with PROOF SlaveBegin() is called on each slave server.
234 // The tree argument is deprecated (on PROOF 0 is passed).
235 if (fDebug > 0) printf("->AliAnalysisManager::SlaveBegin()\n");
236 static Bool_t isCalled = kFALSE;
237 TDirectory *curdir = gDirectory;
238 // Call SlaveBegin only once in case of mixing
239 if (isCalled && fMode==kMixingAnalysis) return;
240 // Call Init of EventHandler
241 if (fOutputEventHandler) {
242 if (fMode == kProofAnalysis) {
243 TIter nextout(fOutputs);
244 AliAnalysisDataContainer *c_aod;
245 while ((c_aod=(AliAnalysisDataContainer*)nextout())) if (!strcmp(c_aod->GetFileName(),"default")) break;
246 if (c_aod && c_aod->IsSpecialOutput()) {
248 if (fDebug > 1) printf(" Initializing special output file %s...\n", fOutputEventHandler->GetOutputFileName());
249 OpenProofFile(fOutputEventHandler->GetOutputFileName(), "RECREATE");
250 c_aod->SetFile(gFile);
251 fOutputEventHandler->Init("proofspecial");
254 fOutputEventHandler->Init("proof");
257 fOutputEventHandler->Init("local");
261 if (fInputEventHandler) {
262 fInputEventHandler->SetInputTree(tree);
263 if (fMode == kProofAnalysis) {
264 fInputEventHandler->Init("proof");
266 fInputEventHandler->Init("local");
270 if (fMCtruthEventHandler) {
271 if (fMode == kProofAnalysis) {
272 fMCtruthEventHandler->Init("proof");
274 fMCtruthEventHandler->Init("local");
277 if (curdir) curdir->cd();
280 AliAnalysisTask *task;
281 // Call CreateOutputObjects for all tasks
282 while ((task=(AliAnalysisTask*)next())) {
284 task->CreateOutputObjects();
285 if (curdir) curdir->cd();
289 if (fDebug > 0) printf("<-AliAnalysisManager::SlaveBegin()\n");
292 //______________________________________________________________________________
293 Bool_t AliAnalysisManager::Notify()
295 // The Notify() function is called when a new file is opened. This
296 // can be either for a new TTree in a TChain or when when a new TTree
297 // is started when using PROOF. It is normaly not necessary to make changes
298 // to the generated code, but the routine can be extended by the
299 // user if needed. The return value is currently not used.
301 Error("Notify","No current tree.");
304 TFile *curfile = fTree->GetCurrentFile();
306 Error("Notify","No current file");
310 if (fDebug > 0) printf("->AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
312 AliAnalysisTask *task;
313 // Call Notify for all tasks
314 while ((task=(AliAnalysisTask*)next()))
317 // Call Notify of the event handlers
318 if (fInputEventHandler) {
319 fInputEventHandler->Notify(curfile->GetName());
322 if (fOutputEventHandler) {
323 fOutputEventHandler->Notify(curfile->GetName());
326 if (fMCtruthEventHandler) {
327 fMCtruthEventHandler->Notify(curfile->GetName());
329 if (fDebug > 0) printf("<-AliAnalysisManager::Notify()\n");
333 //______________________________________________________________________________
334 Bool_t AliAnalysisManager::Process(Long64_t entry)
336 // The Process() function is called for each entry in the tree (or possibly
337 // keyed object in the case of PROOF) to be processed. The entry argument
338 // specifies which entry in the currently loaded tree is to be processed.
339 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
340 // to read either all or the required parts of the data. When processing
341 // keyed objects with PROOF, the object is already loaded and is available
342 // via the fObject pointer.
344 // This function should contain the "body" of the analysis. It can contain
345 // simple or elaborate selection criteria, run algorithms on the data
346 // of the event and typically fill histograms.
348 // WARNING when a selector is used with a TChain, you must use
349 // the pointer to the current TTree to call GetEntry(entry).
350 // The entry is always the local entry number in the current tree.
351 // Assuming that fChain is the pointer to the TChain being processed,
352 // use fChain->GetTree()->GetEntry(entry).
353 if (fDebug > 0) printf("->AliAnalysisManager::Process(%lld)\n", entry);
355 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
356 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
357 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
361 if (fDebug > 0) printf("<-AliAnalysisManager::Process()\n");
365 //______________________________________________________________________________
366 void AliAnalysisManager::PackOutput(TList *target)
368 // Pack all output data containers in the output list. Called at SlaveTerminate
369 // stage in PROOF case for each slave.
370 if (fDebug > 0) printf("->AliAnalysisManager::PackOutput()\n");
372 Error("PackOutput", "No target. Aborting.");
375 if (fInputEventHandler) fInputEventHandler ->Terminate();
376 if (fOutputEventHandler) fOutputEventHandler ->Terminate();
377 if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate();
379 // Call FinishTaskOutput() for each event loop task (not called for
380 // post-event loop tasks - use Terminate() fo those)
381 TIter nexttask(fTasks);
382 AliAnalysisTask *task;
383 while ((task=(AliAnalysisTask*)nexttask())) {
384 if (!task->IsPostEventLoop()) {
385 if (fDebug > 0) printf("->FinishTaskOutput: task %s\n", task->GetName());
386 task->FinishTaskOutput();
387 if (fDebug > 0) printf("<-FinishTaskOutput: task %s\n", task->GetName());
391 if (fMode == kProofAnalysis) {
392 TIter next(fOutputs);
393 AliAnalysisDataContainer *output;
394 Bool_t isManagedByHandler = kFALSE;
395 while ((output=(AliAnalysisDataContainer*)next())) {
396 // Do not consider outputs of post event loop tasks
397 isManagedByHandler = kFALSE;
398 if (output->GetProducer()->IsPostEventLoop()) continue;
399 const char *filename = output->GetFileName();
400 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
401 isManagedByHandler = kTRUE;
402 filename = fOutputEventHandler->GetOutputFileName();
404 // Check if data was posted to this container. If not, issue an error.
405 if (!output->GetData() && !isManagedByHandler) {
406 Error("PackOutput", "No data for output container %s. Forgot to PostData ?\n", output->GetName());
409 if (!output->IsSpecialOutput()) {
411 if (strlen(filename) && !isManagedByHandler) {
412 // File resident outputs
413 TFile *file = output->GetFile();
414 // Backup current folder
415 TDirectory *opwd = gDirectory;
416 // Create file if not existing and register to container.
417 if (file) file->cd();
418 else file = new TFile(filename, "RECREATE");
419 if (file->IsZombie()) {
420 Fatal("PackOutput", "Could not recreate file %s\n", filename);
423 output->SetFile(file);
424 // Clear file list to release object ownership to user.
426 // Save data to file, then close.
427 if (output->GetData()->InheritsFrom(TCollection::Class())) {
428 // If data is a collection, we set the name of the collection
429 // as the one of the container and we save as a single key.
430 TCollection *coll = (TCollection*)output->GetData();
431 coll->SetName(output->GetName());
432 coll->Write(output->GetName(), TObject::kSingleKey);
434 if (output->GetData()->InheritsFrom(TTree::Class())) {
435 TTree *tree = (TTree*)output->GetData();
436 tree->SetDirectory(file);
439 output->GetData()->Write();
442 if (fDebug > 1) printf("PackOutput %s: memory merge, file resident output\n", output->GetName());
444 printf(" file %s listing content:\n", filename);
448 // Restore current directory
449 if (opwd) opwd->cd();
451 // Memory-resident outputs
452 if (fDebug > 1) printf("PackOutput %s: memory merge memory resident output\n", filename);
454 AliAnalysisDataWrapper *wrap = 0;
455 if (isManagedByHandler) {
456 wrap = new AliAnalysisDataWrapper(fOutputEventHandler->GetTree());
457 wrap->SetName(output->GetName());
459 else wrap =output->ExportData();
460 // Output wrappers must NOT delete data after merging - the user owns them
461 wrap->SetDeleteData(kFALSE);
465 TDirectory *opwd = gDirectory;
466 TFile *file = output->GetFile();
467 if (fDebug > 1 && file) printf("PackOutput %s: file merge, special output\n", output->GetName());
468 if (isManagedByHandler) {
469 // Terminate IO for files managed by the output handler
470 if (file) file->Write();
471 if (file && fDebug > 2) {
472 printf(" handled file %s listing content:\n", file->GetName());
475 fOutputEventHandler->TerminateIO();
480 AliAnalysisTask *producer = output->GetProducer();
482 "File %s for special container %s was NOT opened in %s::CreateOutputObjects !!!",
483 output->GetFileName(), output->GetName(), producer->ClassName());
487 // Release object ownership to users after writing data to file
488 if (output->GetData()->InheritsFrom(TCollection::Class())) {
489 // If data is a collection, we set the name of the collection
490 // as the one of the container and we save as a single key.
491 TCollection *coll = (TCollection*)output->GetData();
492 coll->SetName(output->GetName());
493 coll->Write(output->GetName(), TObject::kSingleKey);
495 if (output->GetData()->InheritsFrom(TTree::Class())) {
496 TTree *tree = (TTree*)output->GetData();
497 tree->SetDirectory(file);
500 output->GetData()->Write();
505 printf(" file %s listing content:\n", output->GetFileName());
508 TString outFilename = file->GetName();
510 // Restore current directory
511 if (opwd) opwd->cd();
512 // Check if a special output location was provided or the output files have to be merged
513 if (strlen(fSpecialOutputLocation.Data())) {
514 TString remote = fSpecialOutputLocation;
516 Int_t gid = gROOT->ProcessLine("gProofServ->GetGroupId();");
517 remote += Form("%s_%d_", gSystem->HostName(), gid);
518 remote += output->GetFileName();
519 TFile::Cp ( outFilename.Data(), remote.Data() );
521 // No special location specified-> use TProofOutputFile as merging utility
522 // The file at this output slot must be opened in CreateOutputObjects
523 if (fDebug > 1) printf(" File %s to be merged...\n", output->GetFileName());
528 if (fDebug > 0) printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize());
531 //______________________________________________________________________________
532 void AliAnalysisManager::ImportWrappers(TList *source)
534 // Import data in output containers from wrappers coming in source.
535 if (fDebug > 0) printf("->AliAnalysisManager::ImportWrappers()\n");
536 TIter next(fOutputs);
537 AliAnalysisDataContainer *cont;
538 AliAnalysisDataWrapper *wrap;
540 Bool_t inGrid = (fMode == kGridAnalysis)?kTRUE:kFALSE;
541 while ((cont=(AliAnalysisDataContainer*)next())) {
543 if (cont->GetProducer()->IsPostEventLoop() && !inGrid) continue;
544 const char *filename = cont->GetFileName();
545 Bool_t isManagedByHandler = kFALSE;
546 if (!(strcmp(filename, "default")) && fOutputEventHandler) {
547 isManagedByHandler = kTRUE;
548 filename = fOutputEventHandler->GetOutputFileName();
550 if (cont->IsSpecialOutput() || inGrid) {
551 if (strlen(fSpecialOutputLocation.Data()) && !isManagedByHandler) continue;
552 // Copy merged file from PROOF scratch space.
553 // In case of grid the files are already in the current directory.
557 TObject *pof = source->FindObject(filename);
558 if (!pof || !pof->InheritsFrom("TProofOutputFile")) {
559 Error("ImportWrappers", "TProofOutputFile object not found in output list for container %s", cont->GetName());
562 gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", ((TProofOutputFile*)0x%lx)->GetOutputFileName();)", full_path, pof));
563 gROOT->ProcessLine(Form("sprintf((char*)0x%lx, \"%%s\", gProof->GetUrl();)", ch_url));
564 TString clientUrl(ch_url);
565 TString full_path_str(full_path);
566 if (clientUrl.Contains("localhost")){
567 TObjArray* array = full_path_str.Tokenize ( "//" );
568 TObjString *strobj = ( TObjString *)array->At(1);
569 TObjArray* arrayPort = strobj->GetString().Tokenize ( ":" );
570 TObjString *strobjPort = ( TObjString *) arrayPort->At(1);
571 full_path_str.ReplaceAll(strobj->GetString().Data(),"localhost:PORT");
572 full_path_str.ReplaceAll(":PORT",Form(":%s",strobjPort->GetString().Data()));
573 if (fDebug > 1) Info("ImportWrappers","Using tunnel from %s to %s",full_path_str.Data(),filename);
576 printf(" Copying file %s from PROOF scratch space\n", full_path_str.Data());
577 Bool_t gotit = TFile::Cp(full_path_str.Data(), filename);
579 Error("ImportWrappers", "Could not get file %s from proof scratch space", cont->GetFileName());
583 // Normally we should connect data from the copied file to the
584 // corresponding output container, but it is not obvious how to do this
585 // automatically if several objects in file...
586 TFile *f = TFile::Open(filename, "READ");
588 Error("ImportWrappers", "Cannot open file %s in read-only mode", filename);
592 // Try to fetch first a list object having the container name.
593 obj = f->Get(cont->GetName());
595 // Fetch first object from file having the container type.
596 TIter nextkey(f->GetListOfKeys());
598 while ((key=(TKey*)nextkey())) {
599 obj = f->Get(key->GetName());
600 if (obj && obj->IsA()->InheritsFrom(cont->GetType())) break;
604 Error("ImportWrappers", "Could not find object for container %s in file %s", cont->GetName(), filename);
607 wrap = new AliAnalysisDataWrapper(obj);
608 wrap->SetDeleteData(kFALSE);
610 if (!wrap) wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName());
612 Error("ImportWrappers","Container %s not found in analysis output !", cont->GetName());
617 printf(" Importing data for container %s", cont->GetName());
618 if (strlen(filename)) printf(" -> file %s\n", filename);
621 cont->ImportData(wrap);
623 if (fDebug > 0) printf("<-AliAnalysisManager::ImportWrappers(): %d containers imported\n", icont);
626 //______________________________________________________________________________
627 void AliAnalysisManager::UnpackOutput(TList *source)
629 // Called by AliAnalysisSelector::Terminate only on the client.
630 if (fDebug > 0) printf("->AliAnalysisManager::UnpackOutput()\n");
632 Error("UnpackOutput", "No target. Aborting.");
635 if (fDebug > 1) printf(" Source list contains %d containers\n", source->GetSize());
637 if (fMode == kProofAnalysis) ImportWrappers(source);
639 TIter next(fOutputs);
640 AliAnalysisDataContainer *output;
641 while ((output=(AliAnalysisDataContainer*)next())) {
642 if (!output->GetData()) continue;
643 // Check if there are client tasks that run post event loop
644 if (output->HasConsumers()) {
645 // Disable event loop semaphore
646 output->SetPostEventLoop(kTRUE);
647 TObjArray *list = output->GetConsumers();
648 Int_t ncons = list->GetEntriesFast();
649 for (Int_t i=0; i<ncons; i++) {
650 AliAnalysisTask *task = (AliAnalysisTask*)list->At(i);
651 task->CheckNotify(kTRUE);
652 // If task is active, execute it
653 if (task->IsPostEventLoop() && task->IsActive()) {
654 if (fDebug > 0) printf("== Executing post event loop task %s\n", task->GetName());
660 if (fDebug > 0) printf("<-AliAnalysisManager::UnpackOutput()\n");
663 //______________________________________________________________________________
664 void AliAnalysisManager::Terminate()
666 // The Terminate() function is the last function to be called during
667 // a query. It always runs on the client, it can be used to present
668 // the results graphically.
669 if (fDebug > 0) printf("->AliAnalysisManager::Terminate()\n");
670 AliAnalysisTask *task;
672 // Call Terminate() for tasks
673 while ((task=(AliAnalysisTask*)next())) task->Terminate();
675 TIter next1(fOutputs);
676 AliAnalysisDataContainer *output;
677 while ((output=(AliAnalysisDataContainer*)next1())) {
678 // Special outputs or grid files have the files already closed and written.
679 if (fMode == kGridAnalysis) continue;
680 if (output->IsSpecialOutput()&&(fMode == kProofAnalysis)) continue;
681 const char *filename = output->GetFileName();
682 if (!(strcmp(filename, "default"))) {
683 if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName();
684 TFile *aodfile = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
686 if (fDebug > 1) printf("Writing output handler file: %s\n", filename);
691 if (!strlen(filename)) continue;
692 if (!output->GetData()) continue;
693 TFile *file = output->GetFile();
694 TDirectory *opwd = gDirectory;
695 file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
696 if (!file) file = new TFile(filename, "RECREATE");
697 if (file->IsZombie()) continue;
698 output->SetFile(file);
700 if (fDebug > 1) printf(" writing output data %s to file %s\n", output->GetData()->GetName(), file->GetName());
701 if (output->GetData()->InheritsFrom(TCollection::Class())) {
702 // If data is a collection, we set the name of the collection
703 // as the one of the container and we save as a single key.
704 TCollection *coll = (TCollection*)output->GetData();
705 coll->SetName(output->GetName());
706 coll->Write(output->GetName(), TObject::kSingleKey);
708 if (output->GetData()->InheritsFrom(TTree::Class())) {
709 TTree *tree = (TTree*)output->GetData();
710 tree->SetDirectory(file);
713 output->GetData()->Write();
716 if (opwd) opwd->cd();
719 while ((output=(AliAnalysisDataContainer*)next1())) {
720 // Close all files at output
721 TDirectory *opwd = gDirectory;
722 if (output->GetFile()) output->GetFile()->Close();
723 if (opwd) opwd->cd();
726 if (fInputEventHandler) fInputEventHandler ->TerminateIO();
727 if (fOutputEventHandler) fOutputEventHandler ->TerminateIO();
728 if (fMCtruthEventHandler) fMCtruthEventHandler->TerminateIO();
730 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
732 TDirectory *cdir = gDirectory;
733 TFile f("syswatch.root", "RECREATE");
735 TTree *tree = AliSysInfo::MakeTree("syswatch.log");
736 tree->SetMarkerStyle(kCircle);
737 tree->SetMarkerColor(kBlue);
738 tree->SetMarkerSize(0.5);
739 if (!gROOT->IsBatch()) {
740 tree->SetAlias("event", "id0");
741 tree->SetAlias("memUSED", "pI.fMemVirtual");
742 tree->SetAlias("userCPU", "pI.fCpuUser");
743 TCanvas *c = new TCanvas("SysInfo","SysInfo",10,10,800,600);
744 c->Divide(2,1,0.01,0.01);
746 tree->Draw("memUSED:event","","", 1234567890, 0);
748 tree->Draw("userCPU:event","","", 1234567890, 0);
754 if (cdir) cdir->cd();
756 if (fDebug > 0) printf("<-AliAnalysisManager::Terminate()\n");
759 //______________________________________________________________________________
760 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
762 // Adds a user task to the global list of tasks.
763 if (fTasks->FindObject(task)) {
764 Warning("AddTask", "Task %s: the same object already added to the analysis manager. Not adding.", task->GetName());
767 task->SetActive(kFALSE);
771 //______________________________________________________________________________
772 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
774 // Retreive task by name.
775 if (!fTasks) return NULL;
776 return (AliAnalysisTask*)fTasks->FindObject(name);
779 //______________________________________________________________________________
780 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
781 TClass *datatype, EAliAnalysisContType type, const char *filename)
783 // Create a data container of a certain type. Types can be:
784 // kExchangeContainer = 0, used to exchange date between tasks
785 // kInputContainer = 1, used to store input data
786 // kOutputContainer = 2, used for posting results
787 if (fContainers->FindObject(name)) {
788 Error("CreateContainer","A container named %s already defined !\n",name);
791 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
792 fContainers->Add(cont);
794 case kInputContainer:
797 case kOutputContainer:
799 if (filename && strlen(filename)) {
800 cont->SetFileName(filename);
801 cont->SetDataOwned(kFALSE); // data owned by the file
804 case kExchangeContainer:
810 //______________________________________________________________________________
811 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
812 AliAnalysisDataContainer *cont)
814 // Connect input of an existing task to a data container.
815 if (!fTasks->FindObject(task)) {
817 Info("ConnectInput", "Task %s was not registered. Now owned by analysis manager", task->GetName());
819 Bool_t connected = task->ConnectInput(islot, cont);
823 //______________________________________________________________________________
824 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
825 AliAnalysisDataContainer *cont)
827 // Connect output of an existing task to a data container.
828 if (!fTasks->FindObject(task)) {
830 Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName());
832 Bool_t connected = task->ConnectOutput(islot, cont);
836 //______________________________________________________________________________
837 void AliAnalysisManager::CleanContainers()
839 // Clean data from all containers that have already finished all client tasks.
840 TIter next(fContainers);
841 AliAnalysisDataContainer *cont;
842 while ((cont=(AliAnalysisDataContainer *)next())) {
843 if (cont->IsOwnedData() &&
844 cont->IsDataReady() &&
845 cont->ClientsExecuted()) cont->DeleteData();
849 //______________________________________________________________________________
850 Bool_t AliAnalysisManager::InitAnalysis()
852 // Initialization of analysis chain of tasks. Should be called after all tasks
853 // and data containers are properly connected
854 // Check for input/output containers
856 // Check for top tasks (depending only on input data containers)
857 if (!fTasks->First()) {
858 Error("InitAnalysis", "Analysis has no tasks !");
862 AliAnalysisTask *task;
863 AliAnalysisDataContainer *cont;
866 Bool_t iszombie = kFALSE;
867 Bool_t istop = kTRUE;
869 while ((task=(AliAnalysisTask*)next())) {
872 Int_t ninputs = task->GetNinputs();
873 for (i=0; i<ninputs; i++) {
874 cont = task->GetInputSlot(i)->GetContainer();
882 Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...",
885 if (iszombie) continue;
886 // Check if cont is an input container
887 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
888 // Connect to parent task
892 fTopTasks->Add(task);
896 Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers");
899 // Check now if there are orphan tasks
900 for (i=0; i<ntop; i++) {
901 task = (AliAnalysisTask*)fTopTasks->At(i);
906 while ((task=(AliAnalysisTask*)next())) {
907 if (!task->IsUsed()) {
909 Warning("InitAnalysis", "Task %s is orphan", task->GetName());
912 // Check the task hierarchy (no parent task should depend on data provided
913 // by a daughter task)
914 for (i=0; i<ntop; i++) {
915 task = (AliAnalysisTask*)fTopTasks->At(i);
916 if (task->CheckCircularDeps()) {
917 Error("InitAnalysis", "Found illegal circular dependencies between following tasks:");
922 // Check that all containers feeding post-event loop tasks are in the outputs list
923 TIter nextcont(fContainers); // loop over all containers
924 while ((cont=(AliAnalysisDataContainer*)nextcont())) {
925 if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) {
926 if (cont->HasConsumers()) {
927 // Check if one of the consumers is post event loop
928 TIter nextconsumer(cont->GetConsumers());
929 while ((task=(AliAnalysisTask*)nextconsumer())) {
930 if (task->IsPostEventLoop()) {
938 // Check if all special output containers have a file name provided
939 TIter nextout(fOutputs);
940 while ((cont=(AliAnalysisDataContainer*)nextout())) {
941 if (cont->IsSpecialOutput() && !strlen(cont->GetFileName())) {
942 Error("InitAnalysis", "Wrong container %s : a file name MUST be provided for special outputs", cont->GetName());
950 //______________________________________________________________________________
951 void AliAnalysisManager::PrintStatus(Option_t *option) const
953 // Print task hierarchy.
955 Info("PrintStatus", "Analysis manager %s not initialized : call InitAnalysis() first", GetName());
958 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
960 Info("PrintStatus", "System information will be collected each %lld events", fNSysInfo);
961 TIter next(fTopTasks);
962 AliAnalysisTask *task;
963 while ((task=(AliAnalysisTask*)next()))
964 task->PrintTask(option);
967 //______________________________________________________________________________
968 void AliAnalysisManager::ResetAnalysis()
970 // Reset all execution flags and clean containers.
974 //______________________________________________________________________________
975 void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree, Long64_t nentries, Long64_t firstentry)
977 // Start analysis for this manager. Analysis task can be: LOCAL, PROOF, GRID or
978 // MIX. Process nentries starting from firstentry
980 Error("StartAnalysis","Analysis manager was not initialized !");
983 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
984 TString anaType = type;
986 fMode = kLocalAnalysis;
987 Bool_t runlocalinit = kTRUE;
988 if (anaType.Contains("file")) runlocalinit = kFALSE;
989 if (anaType.Contains("proof")) fMode = kProofAnalysis;
990 else if (anaType.Contains("grid")) fMode = kGridAnalysis;
991 else if (anaType.Contains("mix")) fMode = kMixingAnalysis;
993 if (fMode == kGridAnalysis) {
995 Error("StartAnalysis", "Cannot start grid analysis without a grid handler.");
996 Info("===", "Add an AliAnalysisAlien object as plugin for this manager and configure it.");
999 // Write analysis manager in the analysis file
1000 cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
1001 // run local task configuration
1002 TIter nextTask(fTasks);
1003 AliAnalysisTask *task;
1004 while ((task=(AliAnalysisTask*)nextTask())) {
1007 fGridHandler->StartAnalysis(nentries, firstentry);
1009 // Terminate grid analysis
1010 if (fGridHandler->GetRunMode() == AliAnalysisGrid::kOffline) return;
1011 cout << "===== MERGING OUTPUTS REGISTERED BY YOUR ANALYSIS JOB: " << GetName() << endl;
1012 if (!fGridHandler->MergeOutputs()) {
1013 // Return if outputs could not be merged or if it alien handler
1014 // was configured for offline mode or local testing.
1017 ImportWrappers(NULL);
1022 SetEventLoop(kFALSE);
1023 // Enable event loop mode if a tree was provided
1024 if (tree || fMode==kMixingAnalysis) SetEventLoop(kTRUE);
1027 TString ttype = "TTree";
1028 if (tree && tree->IsA() == TChain::Class()) {
1029 chain = (TChain*)tree;
1030 if (!chain || !chain->GetListOfFiles()->First()) {
1031 Error("StartAnalysis", "Cannot process null or empty chain...");
1037 // Initialize locally all tasks (happens for all modes)
1039 AliAnalysisTask *task;
1041 while ((task=(AliAnalysisTask*)next())) {
1047 case kLocalAnalysis:
1049 TIter nextT(fTasks);
1050 // Call CreateOutputObjects for all tasks
1051 while ((task=(AliAnalysisTask*)nextT())) {
1052 TDirectory *curdir = gDirectory;
1053 task->CreateOutputObjects();
1054 if (curdir) curdir->cd();
1060 // Run tree-based analysis via AliAnalysisSelector
1061 cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl;
1062 fSelector = new AliAnalysisSelector(this);
1063 tree->Process(fSelector, "", nentries, firstentry);
1065 case kProofAnalysis:
1066 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1067 printf("StartAnalysis: no PROOF!!!\n");
1070 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1071 gROOT->ProcessLine(line);
1074 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl;
1075 chain->Process("AliAnalysisSelector", "", nentries, firstentry);
1077 printf("StartAnalysis: no chain\n");
1082 Warning("StartAnalysis", "GRID analysis mode not implemented. Running local.");
1084 case kMixingAnalysis:
1085 // Run event mixing analysis
1087 Error("StartAnalysis", "Cannot run event mixing without event pool");
1090 cout << "===== RUNNING EVENT MIXING ANALYSIS " << GetName() << endl;
1091 fSelector = new AliAnalysisSelector(this);
1092 while ((chain=fEventPool->GetNextChain())) {
1094 // Call NotifyBinChange for all tasks
1095 while ((task=(AliAnalysisTask*)next()))
1096 if (!task->IsPostEventLoop()) task->NotifyBinChange();
1097 chain->Process(fSelector);
1099 PackOutput(fSelector->GetOutputList());
1104 //______________________________________________________________________________
1105 void AliAnalysisManager::StartAnalysis(const char *type, const char *dataset, Long64_t nentries, Long64_t firstentry)
1107 // Start analysis for this manager on a given dataset. Analysis task can be:
1108 // LOCAL, PROOF or GRID. Process nentries starting from firstentry.
1110 Error("StartAnalysis","Analysis manager was not initialized !");
1113 if (fDebug > 0) printf("StartAnalysis %s\n",GetName());
1114 TString anaType = type;
1116 if (!anaType.Contains("proof")) {
1117 Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
1120 fMode = kProofAnalysis;
1122 SetEventLoop(kTRUE);
1123 // Set the dataset flag
1124 TObject::SetBit(kUseDataSet);
1127 // Initialize locally all tasks
1129 AliAnalysisTask *task;
1130 while ((task=(AliAnalysisTask*)next())) {
1134 if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
1135 printf("StartAnalysis: no PROOF!!!\n");
1138 sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this);
1139 gROOT->ProcessLine(line);
1140 sprintf(line, "gProof->GetDataSet(\"%s\");", dataset);
1141 if (!gROOT->ProcessLine(line)) {
1142 Error("StartAnalysis", "Dataset %s not found", dataset);
1145 sprintf(line, "gProof->Process(\"%s\", \"AliAnalysisSelector\", \"\", %lld, %lld);",
1146 dataset, nentries, firstentry);
1147 cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dataset << endl;
1148 gROOT->ProcessLine(line);
1151 //______________________________________________________________________________
1152 TFile *AliAnalysisManager::OpenProofFile(const char *filename, const char *option)
1154 // Opens a special output file used in PROOF.
1156 if (fMode!=kProofAnalysis || !fSelector) {
1157 Error("OpenProofFile","Cannot open PROOF file %s",filename);
1160 sprintf(line, "TProofOutputFile *pf = new TProofOutputFile(\"%s\");", filename);
1161 if (fDebug > 1) printf("=== %s\n", line);
1162 gROOT->ProcessLine(line);
1163 sprintf(line, "pf->OpenFile(\"%s\");", option);
1164 gROOT->ProcessLine(line);
1166 gROOT->ProcessLine("pf->Print()");
1167 printf(" == proof file name: %s\n", gFile->GetName());
1169 sprintf(line, "((TList*)0x%lx)->Add(pf);",(ULong_t)fSelector->GetOutputList());
1170 if (fDebug > 1) printf("=== %s\n", line);
1171 gROOT->ProcessLine(line);
1175 //______________________________________________________________________________
1176 void AliAnalysisManager::ExecAnalysis(Option_t *option)
1178 // Execute analysis.
1179 static Long64_t ncalls = 0;
1180 Bool_t getsysInfo = ((fNSysInfo>0) && (fMode==kLocalAnalysis))?kTRUE:kFALSE;
1181 if (getsysInfo && ncalls==0) AliSysInfo::AddStamp("Start", (Int_t)ncalls);
1184 Error("ExecAnalysis", "Analysis manager was not initialized !");
1187 AliAnalysisTask *task;
1188 // Check if the top tree is active.
1191 // De-activate all tasks
1192 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
1193 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
1195 Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container");
1198 cont->SetData(fTree); // This will notify all consumers
1199 Long64_t entry = fTree->GetTree()->GetReadEntry();
1202 // Call BeginEvent() for optional input/output and MC services
1203 if (fInputEventHandler) fInputEventHandler ->BeginEvent(entry);
1204 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(entry);
1205 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(entry);
1207 // Execute the tasks
1208 // TIter next1(cont->GetConsumers());
1209 TIter next1(fTopTasks);
1210 while ((task=(AliAnalysisTask*)next1())) {
1212 cout << " Executing task " << task->GetName() << endl;
1215 task->ExecuteTask(option);
1218 // Call FinishEvent() for optional output and MC services
1219 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1220 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1221 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1222 // Gather system information if requested
1223 if (getsysInfo && ((ncalls%fNSysInfo)==0))
1224 AliSysInfo::AddStamp(Form("Event#%lld",ncalls),(Int_t)ncalls);
1227 // The event loop is not controlled by TSelector
1229 // Call BeginEvent() for optional input/output and MC services
1230 if (fInputEventHandler) fInputEventHandler ->BeginEvent(-1);
1231 if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(-1);
1232 if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(-1);
1233 TIter next2(fTopTasks);
1234 while ((task=(AliAnalysisTask*)next2())) {
1235 task->SetActive(kTRUE);
1237 cout << " Executing task " << task->GetName() << endl;
1239 task->ExecuteTask(option);
1242 // Call FinishEvent() for optional output and MC services
1243 if (fInputEventHandler) fInputEventHandler ->FinishEvent();
1244 if (fOutputEventHandler) fOutputEventHandler ->FinishEvent();
1245 if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent();
1248 //______________________________________________________________________________
1249 void AliAnalysisManager::FinishAnalysis()