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"
35 #include "AliAnalysisManager.h"
36 #include "AliAnalysisTask.h"
37 #include "AliAnalysisDataContainer.h"
38 #include "AliAnalysisDataSlot.h"
40 ClassImp(AliAnalysisManager)
42 //______________________________________________________________________________
43 AliAnalysisManager::AliAnalysisManager() : TSelector(),
53 // Default constructor.
54 if (TClass::IsCallingNew() != TClass::kDummyNew) {
55 fContainers = new TObjArray();
56 fInputs = new TObjArray();
57 fOutputs = new TObjArray();
58 fTasks = new TObjArray();
59 fTopTasks = new TObjArray();
60 fZombies = new TObjArray();
61 // fStatus = new AliAnalysisInfo(this);
66 //______________________________________________________________________________
67 AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
79 fInitOK = other.fInitOK;
80 fContainers = new TObjArray(*other.fContainers);
81 fInputs = new TObjArray(*other.fInputs);
82 fOutputs = new TObjArray(*other.fOutputs);
83 fTasks = new TObjArray(*other.fTasks);
84 fTopTasks = new TObjArray(*other.fTopTasks);
85 fZombies = new TObjArray(*other.fZombies);
86 // fStatus = new AliAnalysisInfo(this);
89 //______________________________________________________________________________
90 AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other)
94 TSelector::operator=(other);
96 fInitOK = other.fInitOK;
97 fContainers = new TObjArray(*other.fContainers);
98 fInputs = new TObjArray(*other.fInputs);
99 fOutputs = new TObjArray(*other.fOutputs);
100 fTasks = new TObjArray(*other.fTasks);
101 fTopTasks = new TObjArray(*other.fTopTasks);
102 fZombies = new TObjArray(*other.fZombies);
103 // fStatus = new AliAnalysisInfo(this);
109 //______________________________________________________________________________
110 AliAnalysisManager::~AliAnalysisManager()
113 if (fContainers) {fContainers->Delete(); delete fContainers;}
114 if (fInputs) delete fInputs;
115 if (fOutputs) delete fOutputs;
116 if (fTasks) {fTasks->Delete(); delete fTasks;}
117 if (fTopTasks) delete fTopTasks;
118 if (fZombies) delete fZombies;
120 //______________________________________________________________________________
121 Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall)
123 // Read one entry of the tree or a whole branch.
124 return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0;
127 //______________________________________________________________________________
128 void AliAnalysisManager::Init(TTree *tree)
130 // The Init() function is called when the selector needs to initialize
131 // a new tree or chain. Typically here the branch addresses of the tree
132 // will be set. It is normaly not necessary to make changes to the
133 // generated code, but the routine can be extended by the user if needed.
134 // Init() will be called many times when running with PROOF.
135 printf("AliAnalysisManager::Init(%s)\n", tree->GetName());
137 cout<<"You have to call InitAnalysis first"<<endl;
138 //AliError("You have to call InitAnalysis first");
143 AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0);
147 //______________________________________________________________________________
148 void AliAnalysisManager::Begin(TTree *tree)
150 // The Begin() function is called at the start of the query.
151 // When running with PROOF Begin() is only called on the client.
152 // The tree argument is deprecated (on PROOF 0 is passed).
153 printf("AliAnalysisManager::Begin(%s)\n", tree->GetName());
157 //______________________________________________________________________________
158 void AliAnalysisManager::SlaveBegin(TTree *tree)
160 // The SlaveBegin() function is called after the Begin() function.
161 // When running with PROOF SlaveBegin() is called on each slave server.
162 // The tree argument is deprecated (on PROOF 0 is passed).
163 printf("AliAnalysisManager::SlaveBegin(%s)\n", tree->GetName());
167 //______________________________________________________________________________
168 Bool_t AliAnalysisManager::Notify()
170 // The Notify() function is called when a new file is opened. This
171 // can be either for a new TTree in a TChain or when when a new TTree
172 // is started when using PROOF. It is normaly not necessary to make changes
173 // to the generated code, but the routine can be extended by the
174 // user if needed. The return value is currently not used.
176 TFile *curfile = fTree->GetCurrentFile();
177 if (curfile) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName());
182 //______________________________________________________________________________
183 Bool_t AliAnalysisManager::Process(Long64_t entry)
185 // The Process() function is called for each entry in the tree (or possibly
186 // keyed object in the case of PROOF) to be processed. The entry argument
187 // specifies which entry in the currently loaded tree is to be processed.
188 // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
189 // to read either all or the required parts of the data. When processing
190 // keyed objects with PROOF, the object is already loaded and is available
191 // via the fObject pointer.
193 // This function should contain the "body" of the analysis. It can contain
194 // simple or elaborate selection criteria, run algorithms on the data
195 // of the event and typically fill histograms.
197 // WARNING when a selector is used with a TChain, you must use
198 // the pointer to the current TTree to call GetEntry(entry).
199 // The entry is always the local entry number in the current tree.
200 // Assuming that fChain is the pointer to the TChain being processed,
201 // use fChain->GetTree()->GetEntry(entry).
203 // printf("AliAnalysisManager::Process(%lld)\n", entry);
209 //______________________________________________________________________________
210 void AliAnalysisManager::SlaveTerminate()
212 // The SlaveTerminate() function is called after all entries or objects
213 // have been processed. When running with PROOF SlaveTerminate() is called
214 // on each slave server.
216 printf("AliAnalysisManager::SlaveTerminate()\n");
219 cout<<"ERROR: Output list not initialized."<<endl;
220 //AliError("ERROR: Output list not initialized.");
223 TIter next(fOutputs);
224 AliAnalysisDataContainer *output;
225 while ((output=(AliAnalysisDataContainer *)next())) {
226 output->SetDataOwned(kFALSE);
227 fOutput->Add(output->GetData());
231 //______________________________________________________________________________
232 void AliAnalysisManager::Terminate()
234 // The Terminate() function is the last function to be called during
235 // a query. It always runs on the client, it can be used to present
236 // the results graphically or save the results to file.
237 printf("AliAnalysisManager::Terminate()\n");
238 AliAnalysisDataContainer *output;
239 AliAnalysisTask *task;
240 TIter next(fOutputs);
241 while ((output=(AliAnalysisDataContainer *)next())) output->WriteData();
243 // Call Terminate() for tasks
244 while ((task=(AliAnalysisTask*)next1())) task->Terminate();
247 //______________________________________________________________________________
248 void AliAnalysisManager::AddTask(AliAnalysisTask *task)
250 // Adds a user task to the global list of tasks.
251 task->SetActive(kFALSE);
255 //______________________________________________________________________________
256 AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const
258 // Retreive task by name.
259 if (!fTasks) return NULL;
260 return (AliAnalysisTask*)fTasks->FindObject(name);
263 //______________________________________________________________________________
264 AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name,
265 TClass *datatype, EAliAnalysisContType type)
267 // Create a data container of a certain type. Types can be:
268 // kNormalContainer = 0, used to exchange date between tasks
269 // kInputContainer = 1, used to store input data
270 // kOutputContainer = 2, used for posting results
271 AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype);
272 fContainers->Add(cont);
274 case kInputContainer:
277 case kOutputContainer:
280 case kNormalContainer:
286 //______________________________________________________________________________
287 Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot,
288 AliAnalysisDataContainer *cont)
290 // Connect input of an existing task to a data container.
291 if (!fTasks->FindObject(task)) {
293 cout<<"Task "<<task->GetName()<<" not registered. Now owned by analysis manager"<<endl;
294 //AliInfo(Form("Task %s not registered. Now owned by analysis manager", task->GetName()));
296 Bool_t connected = task->ConnectInput(islot, cont);
300 //______________________________________________________________________________
301 Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot,
302 AliAnalysisDataContainer *cont)
304 // Connect output of an existing task to a data container.
305 if (!fTasks->FindObject(task)) {
307 cout<<"Task "<<task->GetName()<<"not registered. Now owned by analysis manager"<<endl;
308 //AliInfo(Form("Task %s not registered. Now owned by analysis manager", task->GetName()));
310 Bool_t connected = task->ConnectOutput(islot, cont);
314 //______________________________________________________________________________
315 void AliAnalysisManager::CleanContainers()
317 // Clean data from all containers that have already finished all client tasks.
318 TIter next(fContainers);
319 AliAnalysisDataContainer *cont;
320 while ((cont=(AliAnalysisDataContainer *)next())) {
321 if (cont->IsOwnedData() &&
322 cont->IsDataReady() &&
323 cont->ClientsExecuted()) cont->DeleteData();
327 //______________________________________________________________________________
328 Bool_t AliAnalysisManager::InitAnalysis()
330 // Initialization of analysis chain of tasks. Should be called after all tasks
331 // and data containers are properly connected
332 // Check for input/output containers
334 // if (!fInputs->GetEntriesFast()) {
335 // AliError("No input container defined. At least one container should store input data");
338 // if (!fOutputs->GetEntriesFast()) {
339 // AliError("No output container defined. At least one container should store output data");
342 // Check for top tasks (depending only on input data containers)
343 if (!fTasks->First()) {
344 cout<<"Analysis has no tasks !"<<endl;
345 //AliError("Analysis have no tasks !");
349 AliAnalysisTask *task;
350 AliAnalysisDataContainer *cont;
353 Bool_t iszombie = kFALSE;
354 Bool_t istop = kTRUE;
356 while ((task=(AliAnalysisTask*)next())) {
359 Int_t ninputs = task->GetNinputs();
361 // task->SetZombie();
362 // fZombies->Add(task);
364 // AliWarning(Form("Task %s has no input slots defined ! Declared zombie...",task->GetName()));
367 for (i=0; i<ninputs; i++) {
368 cont = task->GetInputSlot(i)->GetContainer();
376 cout<<"Input slot "<<i<<" of task "<<task->GetName()<<" has no container connected ! Declared zombie..."<<endl;
377 //AliWarning(Form("Input slot %i of task %s has no container connected ! Declared zombie...",i,task->GetName()));
379 if (iszombie) continue;
380 // Check if cont is an input container
381 if (istop && !fInputs->FindObject(cont)) istop=kFALSE;
382 // Connect to parent task
386 fTopTasks->Add(task);
390 cout<<"No top task defined. At least one task should be connected only to input containers"<<endl;
391 //AliError("No top task defined. At least one task should be connected only to input containers");
394 // Check now if there are orphan tasks
395 for (i=0; i<ntop; i++) {
396 task = (AliAnalysisTask*)fTopTasks->At(i);
401 while ((task=(AliAnalysisTask*)next())) {
402 if (!task->IsUsed()) {
404 cout<<"Task "<<task->GetName()<<" is orphan"<<endl;
405 //AliWarning(Form("Task %s is orphan",task->GetName()));
408 // Check the task hierarchy (no parent task should depend on data provided
409 // by a daughter task)
410 for (i=0; i<ntop; i++) {
411 task = (AliAnalysisTask*)fTopTasks->At(i);
412 if (task->CheckCircularDeps()) {
413 cout<<"Found illegal circular dependencies between following tasks:"<<endl;
414 //AliError("Found illegal circular dependencies between following tasks:");
423 //______________________________________________________________________________
424 void AliAnalysisManager::PrintStatus(Option_t *option) const
426 // Print task hierarchy.
427 TIter next(fTopTasks);
428 AliAnalysisTask *task;
429 while ((task=(AliAnalysisTask*)next()))
430 task->PrintTask(option);
433 //______________________________________________________________________________
434 void AliAnalysisManager::ResetAnalysis()
436 // Reset all execution flags and clean containers.
440 //______________________________________________________________________________
441 void AliAnalysisManager::ExecAnalysis(Option_t *option)
445 cout<<"Analysis manager was not initialized !"<<endl;
446 //AliError("Analysis manager was not initialized !");
449 AliAnalysisTask *task;
450 // Check if the top tree is active.
453 // De-activate all tasks
454 while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE);
455 AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0);
457 cout<<"Cannot execute analysis in TSelector mode without at least one top container"<<endl;
458 //AliError("Cannot execute analysis in TSelector mode without at least one top container");
461 cont->SetData(fTree); // This will notify all consumers
462 TIter next1(cont->GetConsumers());
463 while ((task=(AliAnalysisTask*)next1())) {
464 // task->SetActive(kTRUE);
465 task->ExecuteTask(option);
469 // The event loop is not controlled by TSelector
470 TIter next2(fTopTasks);
471 while ((task=(AliAnalysisTask*)next2())) {
472 task->SetActive(kTRUE);
473 printf("executing %s\n", task->GetName());
474 task->ExecuteTask(option);
478 //______________________________________________________________________________
479 void AliAnalysisManager::FinishAnalysis()