d3106602 |
1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * |
3 | * * |
4 | * Author: The ALICE Off-line Project. * |
5 | * Contributors are mentioned in the code where appropriate. * |
6 | * * |
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 | **************************************************************************/ |
15 | |
16 | /* $Id$ */ |
17 | // Author: Andrei Gheata, 31/05/2006 |
18 | |
19 | //============================================================================== |
20 | // AliAnalysysManager - Manager analysis class. Allows creation of several |
37153431 |
21 | // analysis tasks and data containers storing their input/output. Allows |
d3106602 |
22 | // connecting/chaining tasks via shared data containers. Serializes the current |
23 | // event for all tasks depending only on initial input data. |
24 | //============================================================================== |
25 | // |
26 | //============================================================================== |
27 | |
c52c2132 |
28 | #include <Riostream.h> |
11026a80 |
29 | |
c52c2132 |
30 | #include <TClass.h> |
31 | #include <TFile.h> |
32 | #include <TMethodCall.h> |
33 | #include <TChain.h> |
34 | #include <TSystem.h> |
35 | #include <TROOT.h> |
d3106602 |
36 | |
d3106602 |
37 | #include "AliAnalysisTask.h" |
38 | #include "AliAnalysisDataContainer.h" |
39 | #include "AliAnalysisDataSlot.h" |
d2f1d9ef |
40 | #include "AliVEventHandler.h" |
c52c2132 |
41 | #include "AliAnalysisManager.h" |
d3106602 |
42 | |
43 | ClassImp(AliAnalysisManager) |
44 | |
c52c2132 |
45 | AliAnalysisManager *AliAnalysisManager::fgAnalysisManager = NULL; |
46 | |
d3106602 |
47 | //______________________________________________________________________________ |
c52c2132 |
48 | AliAnalysisManager::AliAnalysisManager() |
49 | :TNamed(), |
327eaf46 |
50 | fTree(NULL), |
54cff064 |
51 | fInputEventHandler(NULL), |
6bb2b24f |
52 | fOutputEventHandler(NULL), |
53 | fMCtruthEventHandler(NULL), |
c52c2132 |
54 | fCurrentEntry(-1), |
55 | fMode(kLocalAnalysis), |
37a26056 |
56 | fInitOK(kFALSE), |
6bb2b24f |
57 | fDebug(0), |
c52c2132 |
58 | fTasks(NULL), |
59 | fTopTasks(NULL), |
60 | fZombies(NULL), |
37a26056 |
61 | fContainers(NULL), |
62 | fInputs(NULL), |
c52c2132 |
63 | fOutputs(NULL) |
64 | { |
65 | // Dummy constructor. |
66 | fgAnalysisManager = this; |
b1310ef5 |
67 | SetEventLoop(kTRUE); |
c52c2132 |
68 | } |
69 | |
70 | //______________________________________________________________________________ |
71 | AliAnalysisManager::AliAnalysisManager(const char *name, const char *title) |
72 | :TNamed(name,title), |
73 | fTree(NULL), |
54cff064 |
74 | fInputEventHandler(NULL), |
6bb2b24f |
75 | fOutputEventHandler(NULL), |
76 | fMCtruthEventHandler(NULL), |
c52c2132 |
77 | fCurrentEntry(-1), |
78 | fMode(kLocalAnalysis), |
79 | fInitOK(kFALSE), |
80 | fDebug(0), |
37a26056 |
81 | fTasks(NULL), |
82 | fTopTasks(NULL), |
c52c2132 |
83 | fZombies(NULL), |
84 | fContainers(NULL), |
85 | fInputs(NULL), |
86 | fOutputs(NULL) |
d3106602 |
87 | { |
88 | // Default constructor. |
c52c2132 |
89 | fgAnalysisManager = this; |
90 | fTasks = new TObjArray(); |
91 | fTopTasks = new TObjArray(); |
92 | fZombies = new TObjArray(); |
93 | fContainers = new TObjArray(); |
94 | fInputs = new TObjArray(); |
37153431 |
95 | fOutputs = new TObjArray(); |
b1310ef5 |
96 | SetEventLoop(kTRUE); |
d3106602 |
97 | } |
98 | |
99 | //______________________________________________________________________________ |
100 | AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other) |
c52c2132 |
101 | :TNamed(other), |
327eaf46 |
102 | fTree(NULL), |
54cff064 |
103 | fInputEventHandler(NULL), |
6bb2b24f |
104 | fOutputEventHandler(NULL), |
105 | fMCtruthEventHandler(NULL), |
c52c2132 |
106 | fCurrentEntry(-1), |
107 | fMode(other.fMode), |
108 | fInitOK(other.fInitOK), |
109 | fDebug(other.fDebug), |
37a26056 |
110 | fTasks(NULL), |
111 | fTopTasks(NULL), |
c52c2132 |
112 | fZombies(NULL), |
113 | fContainers(NULL), |
114 | fInputs(NULL), |
115 | fOutputs(NULL) |
d3106602 |
116 | { |
117 | // Copy constructor. |
37a26056 |
118 | fTasks = new TObjArray(*other.fTasks); |
119 | fTopTasks = new TObjArray(*other.fTopTasks); |
120 | fZombies = new TObjArray(*other.fZombies); |
c52c2132 |
121 | fContainers = new TObjArray(*other.fContainers); |
122 | fInputs = new TObjArray(*other.fInputs); |
123 | fOutputs = new TObjArray(*other.fOutputs); |
124 | fgAnalysisManager = this; |
d3106602 |
125 | } |
126 | |
127 | //______________________________________________________________________________ |
128 | AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other) |
129 | { |
130 | // Assignment |
131 | if (&other != this) { |
c52c2132 |
132 | TNamed::operator=(other); |
54cff064 |
133 | fInputEventHandler = other.fInputEventHandler; |
6bb2b24f |
134 | fOutputEventHandler = other.fOutputEventHandler; |
135 | fMCtruthEventHandler = other.fMCtruthEventHandler; |
c52c2132 |
136 | fTree = NULL; |
137 | fCurrentEntry = -1; |
138 | fMode = other.fMode; |
37a26056 |
139 | fInitOK = other.fInitOK; |
c52c2132 |
140 | fDebug = other.fDebug; |
37a26056 |
141 | fTasks = new TObjArray(*other.fTasks); |
142 | fTopTasks = new TObjArray(*other.fTopTasks); |
143 | fZombies = new TObjArray(*other.fZombies); |
c52c2132 |
144 | fContainers = new TObjArray(*other.fContainers); |
145 | fInputs = new TObjArray(*other.fInputs); |
146 | fOutputs = new TObjArray(*other.fOutputs); |
147 | fgAnalysisManager = this; |
d3106602 |
148 | } |
149 | return *this; |
150 | } |
151 | |
152 | //______________________________________________________________________________ |
153 | AliAnalysisManager::~AliAnalysisManager() |
154 | { |
155 | // Destructor. |
d3106602 |
156 | if (fTasks) {fTasks->Delete(); delete fTasks;} |
157 | if (fTopTasks) delete fTopTasks; |
158 | if (fZombies) delete fZombies; |
c52c2132 |
159 | if (fContainers) {fContainers->Delete(); delete fContainers;} |
160 | if (fInputs) delete fInputs; |
161 | if (fOutputs) delete fOutputs; |
162 | if (fgAnalysisManager==this) fgAnalysisManager = NULL; |
d3106602 |
163 | } |
c52c2132 |
164 | |
d3106602 |
165 | //______________________________________________________________________________ |
327eaf46 |
166 | Int_t AliAnalysisManager::GetEntry(Long64_t entry, Int_t getall) |
167 | { |
168 | // Read one entry of the tree or a whole branch. |
c52c2132 |
169 | if (fDebug > 1) { |
170 | cout << "== AliAnalysisManager::GetEntry()" << endl; |
171 | } |
172 | fCurrentEntry = entry; |
327eaf46 |
173 | return fTree ? fTree->GetTree()->GetEntry(entry, getall) : 0; |
174 | } |
175 | |
176 | //______________________________________________________________________________ |
177 | void AliAnalysisManager::Init(TTree *tree) |
d3106602 |
178 | { |
179 | // The Init() function is called when the selector needs to initialize |
180 | // a new tree or chain. Typically here the branch addresses of the tree |
181 | // will be set. It is normaly not necessary to make changes to the |
182 | // generated code, but the routine can be extended by the user if needed. |
183 | // Init() will be called many times when running with PROOF. |
327eaf46 |
184 | if (!tree) return; |
c52c2132 |
185 | if (fDebug > 1) { |
981f2614 |
186 | printf("->AliAnalysisManager::Init(%s)\n", tree->GetName()); |
c52c2132 |
187 | } |
188 | if (!fInitOK) InitAnalysis(); |
189 | if (!fInitOK) return; |
327eaf46 |
190 | fTree = tree; |
191 | AliAnalysisDataContainer *top = (AliAnalysisDataContainer*)fInputs->At(0); |
c52c2132 |
192 | if (!top) { |
193 | cout<<"Error: No top input container !" <<endl; |
194 | return; |
37153431 |
195 | } |
327eaf46 |
196 | top->SetData(tree); |
981f2614 |
197 | if (fDebug > 1) { |
198 | printf("<-AliAnalysisManager::Init(%s)\n", tree->GetName()); |
199 | } |
d3106602 |
200 | } |
201 | |
202 | //______________________________________________________________________________ |
327eaf46 |
203 | void AliAnalysisManager::Begin(TTree *tree) |
d3106602 |
204 | { |
205 | // The Begin() function is called at the start of the query. |
206 | // When running with PROOF Begin() is only called on the client. |
207 | // The tree argument is deprecated (on PROOF 0 is passed). |
c52c2132 |
208 | if (fDebug > 1) { |
209 | cout << "AliAnalysisManager::Begin()" << endl; |
210 | } |
327eaf46 |
211 | Init(tree); |
d3106602 |
212 | } |
213 | |
214 | //______________________________________________________________________________ |
327eaf46 |
215 | void AliAnalysisManager::SlaveBegin(TTree *tree) |
d3106602 |
216 | { |
217 | // The SlaveBegin() function is called after the Begin() function. |
218 | // When running with PROOF SlaveBegin() is called on each slave server. |
219 | // The tree argument is deprecated (on PROOF 0 is passed). |
c52c2132 |
220 | if (fDebug > 1) { |
981f2614 |
221 | cout << "->AliAnalysisManager::SlaveBegin()" << endl; |
37153431 |
222 | } |
8c9485b2 |
223 | // Call InitIO of EventHandler |
6bb2b24f |
224 | if (fOutputEventHandler) { |
8ca08916 |
225 | if (fMode == kProofAnalysis) { |
6bb2b24f |
226 | fOutputEventHandler->InitIO("proof"); |
8ca08916 |
227 | } else { |
6bb2b24f |
228 | fOutputEventHandler->InitIO("local"); |
8ca08916 |
229 | } |
bcce695f |
230 | } |
54cff064 |
231 | if (fInputEventHandler) { |
232 | fInputEventHandler->SetInputTree(tree); |
233 | fInputEventHandler->InitIO(""); |
234 | } |
235 | |
8c9485b2 |
236 | // |
c52c2132 |
237 | TIter next(fTasks); |
238 | AliAnalysisTask *task; |
239 | // Call CreateOutputObjects for all tasks |
c5a87c56 |
240 | while ((task=(AliAnalysisTask*)next())) { |
241 | TDirectory *curdir = gDirectory; |
c52c2132 |
242 | task->CreateOutputObjects(); |
c5a87c56 |
243 | if (curdir) curdir->cd(); |
244 | } |
c52c2132 |
245 | if (fMode == kLocalAnalysis) Init(tree); |
981f2614 |
246 | if (fDebug > 1) { |
247 | cout << "<-AliAnalysisManager::SlaveBegin()" << endl; |
248 | } |
d3106602 |
249 | } |
250 | |
251 | //______________________________________________________________________________ |
327eaf46 |
252 | Bool_t AliAnalysisManager::Notify() |
253 | { |
254 | // The Notify() function is called when a new file is opened. This |
255 | // can be either for a new TTree in a TChain or when when a new TTree |
256 | // is started when using PROOF. It is normaly not necessary to make changes |
257 | // to the generated code, but the routine can be extended by the |
258 | // user if needed. The return value is currently not used. |
53faeca4 |
259 | if (fTree) { |
260 | TFile *curfile = fTree->GetCurrentFile(); |
c52c2132 |
261 | if (curfile && fDebug>1) printf("AliAnalysisManager::Notify() file: %s\n", curfile->GetName()); |
981f2614 |
262 | TIter next(fTasks); |
263 | AliAnalysisTask *task; |
264 | // Call Notify for all tasks |
265 | while ((task=(AliAnalysisTask*)next())) |
53faeca4 |
266 | task->Notify(); |
890126ab |
267 | |
6bb2b24f |
268 | // Call Notify of the MC truth handler |
890126ab |
269 | if (fMCtruthEventHandler) { |
53faeca4 |
270 | fMCtruthEventHandler->Notify(curfile->GetName()); |
890126ab |
271 | } |
327eaf46 |
272 | } |
273 | return kTRUE; |
274 | } |
275 | |
276 | //______________________________________________________________________________ |
277 | Bool_t AliAnalysisManager::Process(Long64_t entry) |
d3106602 |
278 | { |
279 | // The Process() function is called for each entry in the tree (or possibly |
280 | // keyed object in the case of PROOF) to be processed. The entry argument |
281 | // specifies which entry in the currently loaded tree is to be processed. |
282 | // It can be passed to either TTree::GetEntry() or TBranch::GetEntry() |
283 | // to read either all or the required parts of the data. When processing |
284 | // keyed objects with PROOF, the object is already loaded and is available |
285 | // via the fObject pointer. |
286 | // |
287 | // This function should contain the "body" of the analysis. It can contain |
288 | // simple or elaborate selection criteria, run algorithms on the data |
289 | // of the event and typically fill histograms. |
290 | |
291 | // WARNING when a selector is used with a TChain, you must use |
292 | // the pointer to the current TTree to call GetEntry(entry). |
293 | // The entry is always the local entry number in the current tree. |
294 | // Assuming that fChain is the pointer to the TChain being processed, |
295 | // use fChain->GetTree()->GetEntry(entry). |
c52c2132 |
296 | if (fDebug > 1) { |
981f2614 |
297 | cout << "->AliAnalysisManager::Process()" << endl; |
37153431 |
298 | } |
6bb2b24f |
299 | |
300 | if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(); |
301 | if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(); |
302 | |
327eaf46 |
303 | GetEntry(entry); |
304 | ExecAnalysis(); |
981f2614 |
305 | if (fDebug > 1) { |
306 | cout << "<-AliAnalysisManager::Process()" << endl; |
307 | } |
327eaf46 |
308 | return kTRUE; |
d3106602 |
309 | } |
310 | |
311 | //______________________________________________________________________________ |
c52c2132 |
312 | void AliAnalysisManager::PackOutput(TList *target) |
d3106602 |
313 | { |
981f2614 |
314 | // Pack all output data containers in the output list. Called at SlaveTerminate |
315 | // stage in PROOF case for each slave. |
c52c2132 |
316 | if (fDebug > 1) { |
981f2614 |
317 | cout << "->AliAnalysisManager::PackOutput()" << endl; |
c52c2132 |
318 | } |
319 | if (!target) { |
320 | Error("PackOutput", "No target. Aborting."); |
321 | return; |
37153431 |
322 | } |
c52c2132 |
323 | |
6bb2b24f |
324 | if (fOutputEventHandler) fOutputEventHandler ->Terminate(); |
325 | if (fMCtruthEventHandler) fMCtruthEventHandler->Terminate(); |
8c9485b2 |
326 | |
c52c2132 |
327 | if (fMode == kProofAnalysis) { |
328 | TIter next(fOutputs); |
329 | AliAnalysisDataContainer *output; |
330 | while ((output=(AliAnalysisDataContainer*)next())) { |
981f2614 |
331 | if (fDebug > 1) printf(" Packing container %s...\n", output->GetName()); |
332 | if (output->GetData()) target->Add(output->ExportData()); |
c52c2132 |
333 | } |
334 | } |
c52c2132 |
335 | if (fDebug > 1) { |
981f2614 |
336 | printf("<-AliAnalysisManager::PackOutput: output list contains %d containers\n", target->GetSize()); |
37153431 |
337 | } |
c52c2132 |
338 | } |
339 | |
340 | //______________________________________________________________________________ |
981f2614 |
341 | void AliAnalysisManager::ImportWrappers(TList *source) |
c52c2132 |
342 | { |
981f2614 |
343 | // Import data in output containers from wrappers coming in source. |
344 | if (fDebug > 1) { |
345 | cout << "->AliAnalysisManager::ImportWrappers()" << endl; |
346 | } |
327eaf46 |
347 | TIter next(fOutputs); |
981f2614 |
348 | AliAnalysisDataContainer *cont; |
349 | AliAnalysisDataWrapper *wrap; |
350 | Int_t icont = 0; |
c52c2132 |
351 | while ((cont=(AliAnalysisDataContainer*)next())) { |
981f2614 |
352 | wrap = (AliAnalysisDataWrapper*)source->FindObject(cont->GetName()); |
353 | if (!wrap && fDebug>1) { |
354 | printf("(WW) ImportWrappers: container %s not found in analysis output !\n", cont->GetName()); |
c52c2132 |
355 | continue; |
356 | } |
981f2614 |
357 | icont++; |
358 | if (fDebug > 1) printf(" Importing data for container %s\n", wrap->GetName()); |
359 | if (cont->GetFileName()) printf(" -> %s\n", cont->GetFileName()); |
360 | cont->ImportData(wrap); |
c52c2132 |
361 | } |
981f2614 |
362 | if (fDebug > 1) { |
363 | cout << "<-AliAnalysisManager::ImportWrappers(): "<< icont << " containers imported" << endl; |
364 | } |
c52c2132 |
365 | } |
366 | |
367 | //______________________________________________________________________________ |
368 | void AliAnalysisManager::UnpackOutput(TList *source) |
369 | { |
370 | // Called by AliAnalysisSelector::Terminate. Output containers should |
371 | // be in source in the same order as in fOutputs. |
981f2614 |
372 | if (fDebug > 1) { |
373 | cout << "->AliAnalysisManager::UnpackOutput()" << endl; |
374 | } |
c52c2132 |
375 | if (!source) { |
981f2614 |
376 | Error("UnpackOutput", "No target. Aborting."); |
c52c2132 |
377 | return; |
378 | } |
379 | if (fDebug > 1) { |
c52c2132 |
380 | printf(" Source list contains %d containers\n", source->GetSize()); |
327eaf46 |
381 | } |
c52c2132 |
382 | |
981f2614 |
383 | if (fMode == kProofAnalysis) ImportWrappers(source); |
37153431 |
384 | |
981f2614 |
385 | TIter next(fOutputs); |
c52c2132 |
386 | AliAnalysisDataContainer *output; |
387 | while ((output=(AliAnalysisDataContainer*)next())) { |
c52c2132 |
388 | if (!output->GetData()) continue; |
b1310ef5 |
389 | // Check if there are client tasks that run post event loop |
390 | if (output->HasConsumers()) { |
391 | // Disable event loop semaphore |
392 | output->SetPostEventLoop(kTRUE); |
393 | TObjArray *list = output->GetConsumers(); |
394 | Int_t ncons = list->GetEntriesFast(); |
395 | for (Int_t i=0; i<ncons; i++) { |
396 | AliAnalysisTask *task = (AliAnalysisTask*)list->At(i); |
397 | task->CheckNotify(kTRUE); |
398 | // If task is active, execute it |
399 | if (task->IsPostEventLoop() && task->IsActive()) { |
400 | if (fDebug > 1) { |
401 | cout << "== Executing post event loop task " << task->GetName() << endl; |
402 | } |
403 | task->ExecuteTask(); |
404 | } |
405 | } |
406 | } |
c52c2132 |
407 | // Check if the output need to be written to a file. |
408 | const char *filename = output->GetFileName(); |
8ca08916 |
409 | if (!(strcmp(filename, "default"))) { |
6bb2b24f |
410 | if (fOutputEventHandler) filename = fOutputEventHandler->GetOutputFileName(); |
8ca08916 |
411 | } |
412 | |
c52c2132 |
413 | if (!filename || !strlen(filename)) continue; |
414 | TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename); |
415 | if (file) file->cd(); |
416 | else file = new TFile(filename, "RECREATE"); |
417 | if (file->IsZombie()) continue; |
418 | // Reparent data to this file |
419 | TMethodCall callEnv; |
420 | if (output->GetData()->IsA()) |
421 | callEnv.InitWithPrototype(output->GetData()->IsA(), "SetDirectory", "TDirectory*"); |
422 | if (callEnv.IsValid()) { |
423 | callEnv.SetParam((Long_t) file); |
424 | callEnv.Execute(output->GetData()); |
425 | } |
b1109411 |
426 | output->GetData()->Write(); |
c52c2132 |
427 | } |
981f2614 |
428 | if (fDebug > 1) { |
429 | cout << "<-AliAnalysisManager::UnpackOutput()" << endl; |
430 | } |
d3106602 |
431 | } |
432 | |
433 | //______________________________________________________________________________ |
434 | void AliAnalysisManager::Terminate() |
435 | { |
436 | // The Terminate() function is the last function to be called during |
437 | // a query. It always runs on the client, it can be used to present |
c52c2132 |
438 | // the results graphically. |
439 | if (fDebug > 1) { |
981f2614 |
440 | cout << "->AliAnalysisManager::Terminate()" << endl; |
c52c2132 |
441 | } |
327eaf46 |
442 | AliAnalysisTask *task; |
c52c2132 |
443 | TIter next(fTasks); |
327eaf46 |
444 | // Call Terminate() for tasks |
c52c2132 |
445 | while ((task=(AliAnalysisTask*)next())) task->Terminate(); |
981f2614 |
446 | if (fDebug > 1) { |
447 | cout << "<-AliAnalysisManager::Terminate()" << endl; |
448 | } |
8c9485b2 |
449 | // |
6bb2b24f |
450 | if (fOutputEventHandler) fOutputEventHandler->TerminateIO(); |
d3106602 |
451 | } |
452 | |
453 | //______________________________________________________________________________ |
454 | void AliAnalysisManager::AddTask(AliAnalysisTask *task) |
455 | { |
456 | // Adds a user task to the global list of tasks. |
457 | task->SetActive(kFALSE); |
458 | fTasks->Add(task); |
459 | } |
460 | |
461 | //______________________________________________________________________________ |
462 | AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const |
463 | { |
464 | // Retreive task by name. |
465 | if (!fTasks) return NULL; |
466 | return (AliAnalysisTask*)fTasks->FindObject(name); |
467 | } |
468 | |
469 | //______________________________________________________________________________ |
470 | AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, |
c52c2132 |
471 | TClass *datatype, EAliAnalysisContType type, const char *filename) |
d3106602 |
472 | { |
473 | // Create a data container of a certain type. Types can be: |
c52c2132 |
474 | // kExchangeContainer = 0, used to exchange date between tasks |
d3106602 |
475 | // kInputContainer = 1, used to store input data |
476 | // kOutputContainer = 2, used for posting results |
b1310ef5 |
477 | if (fContainers->FindObject(name)) { |
478 | Error("CreateContainer","A container named %s already defined !\n",name); |
479 | return NULL; |
480 | } |
d3106602 |
481 | AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype); |
482 | fContainers->Add(cont); |
483 | switch (type) { |
484 | case kInputContainer: |
485 | fInputs->Add(cont); |
486 | break; |
487 | case kOutputContainer: |
488 | fOutputs->Add(cont); |
c52c2132 |
489 | if (filename && strlen(filename)) cont->SetFileName(filename); |
d3106602 |
490 | break; |
c52c2132 |
491 | case kExchangeContainer: |
d3106602 |
492 | break; |
493 | } |
494 | return cont; |
495 | } |
496 | |
497 | //______________________________________________________________________________ |
498 | Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot, |
499 | AliAnalysisDataContainer *cont) |
500 | { |
501 | // Connect input of an existing task to a data container. |
502 | if (!fTasks->FindObject(task)) { |
503 | AddTask(task); |
c52c2132 |
504 | Warning("ConnectInput", "Task %s not registered. Now owned by analysis manager", task->GetName()); |
d3106602 |
505 | } |
506 | Bool_t connected = task->ConnectInput(islot, cont); |
507 | return connected; |
508 | } |
509 | |
510 | //______________________________________________________________________________ |
511 | Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot, |
512 | AliAnalysisDataContainer *cont) |
513 | { |
514 | // Connect output of an existing task to a data container. |
515 | if (!fTasks->FindObject(task)) { |
516 | AddTask(task); |
c52c2132 |
517 | Warning("ConnectOutput", "Task %s not registered. Now owned by analysis manager", task->GetName()); |
d3106602 |
518 | } |
519 | Bool_t connected = task->ConnectOutput(islot, cont); |
520 | return connected; |
521 | } |
522 | |
523 | //______________________________________________________________________________ |
524 | void AliAnalysisManager::CleanContainers() |
525 | { |
526 | // Clean data from all containers that have already finished all client tasks. |
527 | TIter next(fContainers); |
528 | AliAnalysisDataContainer *cont; |
529 | while ((cont=(AliAnalysisDataContainer *)next())) { |
530 | if (cont->IsOwnedData() && |
531 | cont->IsDataReady() && |
532 | cont->ClientsExecuted()) cont->DeleteData(); |
533 | } |
534 | } |
535 | |
536 | //______________________________________________________________________________ |
537 | Bool_t AliAnalysisManager::InitAnalysis() |
538 | { |
539 | // Initialization of analysis chain of tasks. Should be called after all tasks |
540 | // and data containers are properly connected |
541 | // Check for input/output containers |
542 | fInitOK = kFALSE; |
d3106602 |
543 | // Check for top tasks (depending only on input data containers) |
544 | if (!fTasks->First()) { |
c52c2132 |
545 | Error("InitAnalysis", "Analysis has no tasks !"); |
d3106602 |
546 | return kFALSE; |
547 | } |
548 | TIter next(fTasks); |
549 | AliAnalysisTask *task; |
550 | AliAnalysisDataContainer *cont; |
551 | Int_t ntop = 0; |
552 | Int_t nzombies = 0; |
327eaf46 |
553 | Bool_t iszombie = kFALSE; |
554 | Bool_t istop = kTRUE; |
d3106602 |
555 | Int_t i; |
556 | while ((task=(AliAnalysisTask*)next())) { |
327eaf46 |
557 | istop = kTRUE; |
558 | iszombie = kFALSE; |
d3106602 |
559 | Int_t ninputs = task->GetNinputs(); |
d3106602 |
560 | for (i=0; i<ninputs; i++) { |
561 | cont = task->GetInputSlot(i)->GetContainer(); |
562 | if (!cont) { |
327eaf46 |
563 | if (!iszombie) { |
d3106602 |
564 | task->SetZombie(); |
565 | fZombies->Add(task); |
566 | nzombies++; |
327eaf46 |
567 | iszombie = kTRUE; |
d3106602 |
568 | } |
c52c2132 |
569 | Error("InitAnalysis", "Input slot %d of task %s has no container connected ! Declared zombie...", |
570 | i, task->GetName()); |
d3106602 |
571 | } |
327eaf46 |
572 | if (iszombie) continue; |
d3106602 |
573 | // Check if cont is an input container |
327eaf46 |
574 | if (istop && !fInputs->FindObject(cont)) istop=kFALSE; |
d3106602 |
575 | // Connect to parent task |
576 | } |
327eaf46 |
577 | if (istop) { |
d3106602 |
578 | ntop++; |
579 | fTopTasks->Add(task); |
580 | } |
581 | } |
582 | if (!ntop) { |
c52c2132 |
583 | Error("InitAnalysis", "No top task defined. At least one task should be connected only to input containers"); |
d3106602 |
584 | return kFALSE; |
585 | } |
586 | // Check now if there are orphan tasks |
587 | for (i=0; i<ntop; i++) { |
588 | task = (AliAnalysisTask*)fTopTasks->At(i); |
589 | task->SetUsed(); |
590 | } |
591 | Int_t norphans = 0; |
592 | next.Reset(); |
593 | while ((task=(AliAnalysisTask*)next())) { |
594 | if (!task->IsUsed()) { |
595 | norphans++; |
c52c2132 |
596 | Warning("InitAnalysis", "Task %s is orphan", task->GetName()); |
d3106602 |
597 | } |
598 | } |
599 | // Check the task hierarchy (no parent task should depend on data provided |
600 | // by a daughter task) |
601 | for (i=0; i<ntop; i++) { |
602 | task = (AliAnalysisTask*)fTopTasks->At(i); |
603 | if (task->CheckCircularDeps()) { |
c52c2132 |
604 | Error("InitAnalysis", "Found illegal circular dependencies between following tasks:"); |
d3106602 |
605 | PrintStatus("dep"); |
606 | return kFALSE; |
607 | } |
608 | } |
b1310ef5 |
609 | // Check that all containers feeding post-event loop tasks are in the outputs list |
610 | TIter nextcont(fContainers); // loop over all containers |
611 | while ((cont=(AliAnalysisDataContainer*)nextcont())) { |
612 | if (!cont->IsPostEventLoop() && !fOutputs->FindObject(cont)) { |
613 | if (cont->HasConsumers()) { |
614 | // Check if one of the consumers is post event loop |
615 | TIter nextconsumer(cont->GetConsumers()); |
616 | while ((task=(AliAnalysisTask*)nextconsumer())) { |
617 | if (task->IsPostEventLoop()) { |
618 | fOutputs->Add(cont); |
619 | break; |
620 | } |
621 | } |
622 | } |
623 | } |
624 | } |
327eaf46 |
625 | fInitOK = kTRUE; |
d3106602 |
626 | return kTRUE; |
627 | } |
628 | |
629 | //______________________________________________________________________________ |
630 | void AliAnalysisManager::PrintStatus(Option_t *option) const |
631 | { |
632 | // Print task hierarchy. |
633 | TIter next(fTopTasks); |
634 | AliAnalysisTask *task; |
635 | while ((task=(AliAnalysisTask*)next())) |
636 | task->PrintTask(option); |
637 | } |
638 | |
639 | //______________________________________________________________________________ |
640 | void AliAnalysisManager::ResetAnalysis() |
641 | { |
642 | // Reset all execution flags and clean containers. |
643 | CleanContainers(); |
644 | } |
645 | |
c52c2132 |
646 | //______________________________________________________________________________ |
647 | void AliAnalysisManager::StartAnalysis(const char *type, TTree *tree) |
648 | { |
649 | // Start analysis for this manager. Analysis task can be: LOCAL, PROOF or GRID. |
650 | if (!fInitOK) { |
651 | Error("StartAnalysis","Analysis manager was not initialized !"); |
652 | return; |
653 | } |
654 | if (fDebug>1) { |
655 | cout << "StartAnalysis: " << GetName() << endl; |
656 | } |
657 | TString anaType = type; |
658 | anaType.ToLower(); |
659 | fMode = kLocalAnalysis; |
660 | if (tree) { |
661 | if (anaType.Contains("proof")) fMode = kProofAnalysis; |
662 | else if (anaType.Contains("grid")) fMode = kGridAnalysis; |
663 | } |
664 | if (fMode == kGridAnalysis) { |
665 | Warning("StartAnalysis", "GRID analysis mode not implemented. Running local."); |
981f2614 |
666 | fMode = kLocalAnalysis; |
667 | } |
efd53803 |
668 | char line[128]; |
669 | SetEventLoop(kFALSE); |
b1310ef5 |
670 | // Disable all branches if requested and set event loop mode |
efd53803 |
671 | if (tree) { |
b1310ef5 |
672 | if (TestBit(kDisableBranches)) { |
673 | printf("Disabling all branches...\n"); |
674 | // tree->SetBranchStatus("*",0); // not yet working |
675 | } |
efd53803 |
676 | SetEventLoop(kTRUE); |
677 | } |
efd53803 |
678 | |
c52c2132 |
679 | TChain *chain = dynamic_cast<TChain*>(tree); |
9b33830a |
680 | |
681 | // Initialize locally all tasks |
682 | TIter next(fTasks); |
683 | AliAnalysisTask *task; |
efd53803 |
684 | while ((task=(AliAnalysisTask*)next())) { |
efd53803 |
685 | task->LocalInit(); |
686 | } |
687 | |
c52c2132 |
688 | switch (fMode) { |
689 | case kLocalAnalysis: |
690 | if (!tree) { |
981f2614 |
691 | TIter next(fTasks); |
692 | AliAnalysisTask *task; |
693 | // Call CreateOutputObjects for all tasks |
c5a87c56 |
694 | while ((task=(AliAnalysisTask*)next())) { |
695 | TDirectory *curdir = gDirectory; |
696 | task->CreateOutputObjects(); |
697 | if (curdir) curdir->cd(); |
698 | } |
c52c2132 |
699 | ExecAnalysis(); |
981f2614 |
700 | Terminate(); |
c52c2132 |
701 | return; |
702 | } |
703 | // Run tree-based analysis via AliAnalysisSelector |
8eedd442 |
704 | // gROOT->ProcessLine(".L $ALICE_ROOT/ANALYSIS/AliAnalysisSelector.cxx+"); |
c52c2132 |
705 | cout << "===== RUNNING LOCAL ANALYSIS " << GetName() << " ON TREE " << tree->GetName() << endl; |
706 | sprintf(line, "AliAnalysisSelector *selector = new AliAnalysisSelector((AliAnalysisManager*)0x%lx);",(ULong_t)this); |
707 | gROOT->ProcessLine(line); |
708 | sprintf(line, "((TTree*)0x%lx)->Process(selector);",(ULong_t)tree); |
709 | gROOT->ProcessLine(line); |
710 | break; |
711 | case kProofAnalysis: |
712 | if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) { |
713 | printf("StartAnalysis: no PROOF!!!\n"); |
714 | return; |
715 | } |
716 | sprintf(line, "gProof->AddInput((TObject*)0x%lx);", (ULong_t)this); |
717 | gROOT->ProcessLine(line); |
718 | if (chain) { |
719 | chain->SetProof(); |
720 | cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON CHAIN " << chain->GetName() << endl; |
8eedd442 |
721 | chain->Process("AliAnalysisSelector"); |
c52c2132 |
722 | } else { |
723 | printf("StartAnalysis: no chain\n"); |
724 | return; |
725 | } |
726 | break; |
727 | case kGridAnalysis: |
728 | Warning("StartAnalysis", "GRID analysis mode not implemented. Running local."); |
729 | } |
730 | } |
731 | |
d3106602 |
732 | //______________________________________________________________________________ |
733 | void AliAnalysisManager::ExecAnalysis(Option_t *option) |
734 | { |
735 | // Execute analysis. |
327eaf46 |
736 | if (!fInitOK) { |
c52c2132 |
737 | Error("ExecAnalysis", "Analysis manager was not initialized !"); |
327eaf46 |
738 | return; |
739 | } |
d3106602 |
740 | AliAnalysisTask *task; |
327eaf46 |
741 | // Check if the top tree is active. |
742 | if (fTree) { |
c52c2132 |
743 | if (fDebug>1) { |
744 | printf("AliAnalysisManager::ExecAnalysis\n"); |
745 | } |
327eaf46 |
746 | TIter next(fTasks); |
747 | // De-activate all tasks |
748 | while ((task=(AliAnalysisTask*)next())) task->SetActive(kFALSE); |
749 | AliAnalysisDataContainer *cont = (AliAnalysisDataContainer*)fInputs->At(0); |
750 | if (!cont) { |
c52c2132 |
751 | Error("ExecAnalysis","Cannot execute analysis in TSelector mode without at least one top container"); |
327eaf46 |
752 | return; |
753 | } |
754 | cont->SetData(fTree); // This will notify all consumers |
6bb2b24f |
755 | // |
756 | // Call BeginEvent() for optional output and MC services |
757 | if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(); |
758 | if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(); |
759 | // |
760 | // Execute the tasks |
327eaf46 |
761 | TIter next1(cont->GetConsumers()); |
762 | while ((task=(AliAnalysisTask*)next1())) { |
c52c2132 |
763 | if (fDebug >1) { |
764 | cout << " Executing task " << task->GetName() << endl; |
765 | } |
6bb2b24f |
766 | |
327eaf46 |
767 | task->ExecuteTask(option); |
768 | } |
6bb2b24f |
769 | // |
770 | // Call FinishEvent() for optional output and MC services |
771 | if (fOutputEventHandler) fOutputEventHandler ->FinishEvent(); |
772 | if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent(); |
773 | // |
327eaf46 |
774 | return; |
775 | } |
776 | // The event loop is not controlled by TSelector |
6bb2b24f |
777 | // |
778 | // Call BeginEvent() for optional output and MC services |
779 | if (fOutputEventHandler) fOutputEventHandler ->BeginEvent(); |
780 | if (fMCtruthEventHandler) fMCtruthEventHandler->BeginEvent(); |
327eaf46 |
781 | TIter next2(fTopTasks); |
782 | while ((task=(AliAnalysisTask*)next2())) { |
783 | task->SetActive(kTRUE); |
c52c2132 |
784 | if (fDebug > 1) { |
785 | cout << " Executing task " << task->GetName() << endl; |
786 | } |
d3106602 |
787 | task->ExecuteTask(option); |
327eaf46 |
788 | } |
6bb2b24f |
789 | // |
790 | // Call FinishEvent() for optional output and MC services |
791 | if (fOutputEventHandler) fOutputEventHandler->FinishEvent(); |
792 | if (fMCtruthEventHandler) fMCtruthEventHandler->FinishEvent(); |
d3106602 |
793 | } |
794 | |
795 | //______________________________________________________________________________ |
796 | void AliAnalysisManager::FinishAnalysis() |
797 | { |
798 | // Finish analysis. |
799 | } |