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 |
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 | //============================================================================== |
25 | // |
26 | //============================================================================== |
27 | |
28 | #include "TClass.h" |
29 | #include "AliLog.h" |
30 | |
31 | #include "AliAnalysisManager.h" |
32 | #include "AliAnalysisTask.h" |
33 | #include "AliAnalysisDataContainer.h" |
34 | #include "AliAnalysisDataSlot.h" |
35 | |
36 | ClassImp(AliAnalysisManager) |
37 | |
38 | //______________________________________________________________________________ |
39 | AliAnalysisManager::AliAnalysisManager() : TSelector() |
40 | { |
41 | // Default constructor. |
42 | if (TClass::IsCallingNew() == TClass::kDummyNew) { |
43 | fInitOK = kFALSE; |
44 | fContainers = 0; |
45 | fInputs = 0; |
46 | fOutputs = 0; |
47 | fTasks = 0; |
48 | fTopTasks = 0; |
49 | fZombies = 0; |
50 | fStatus = 0; |
51 | } else { |
52 | fInitOK = kFALSE; |
53 | fContainers = new TObjArray(); |
54 | fInputs = new TObjArray(); |
55 | fOutputs = new TObjArray(); |
56 | fTasks = new TObjArray(); |
57 | fTopTasks = new TObjArray(); |
58 | fZombies = new TObjArray(); |
59 | // fStatus = new AliAnalysisInfo(this); |
60 | } |
61 | } |
62 | |
63 | //______________________________________________________________________________ |
64 | AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other) |
65 | :TSelector(other) |
66 | { |
67 | // Copy constructor. |
68 | fInitOK = other.IsInitialized(); |
69 | fContainers = new TObjArray(*other.GetContainers()); |
70 | fInputs = new TObjArray(*other.GetInputs()); |
71 | fOutputs = new TObjArray(*other.GetOutputs()); |
72 | fTasks = new TObjArray(*other.GetTasks()); |
73 | fTopTasks = new TObjArray(*other.GetTopTasks()); |
74 | fZombies = new TObjArray(*other.GetZombieTasks()); |
75 | // fStatus = new AliAnalysisInfo(this); |
76 | } |
77 | |
78 | //______________________________________________________________________________ |
79 | AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& other) |
80 | { |
81 | // Assignment |
82 | if (&other != this) { |
83 | TSelector::operator=(other); |
84 | fInitOK = other.IsInitialized(); |
85 | fContainers = new TObjArray(*other.GetContainers()); |
86 | fInputs = new TObjArray(*other.GetInputs()); |
87 | fOutputs = new TObjArray(*other.GetOutputs()); |
88 | fTasks = new TObjArray(*other.GetTasks()); |
89 | fTopTasks = new TObjArray(*other.GetTopTasks()); |
90 | fZombies = new TObjArray(*other.GetZombieTasks()); |
91 | // fStatus = new AliAnalysisInfo(this); |
92 | } |
93 | return *this; |
94 | } |
95 | |
96 | //______________________________________________________________________________ |
97 | AliAnalysisManager::~AliAnalysisManager() |
98 | { |
99 | // Destructor. |
100 | if (fContainers) {fContainers->Delete(); delete fContainers;} |
101 | if (fInputs) delete fInputs; |
102 | if (fOutputs) delete fOutputs; |
103 | if (fTasks) {fTasks->Delete(); delete fTasks;} |
104 | if (fTopTasks) delete fTopTasks; |
105 | if (fZombies) delete fZombies; |
106 | } |
107 | |
108 | //______________________________________________________________________________ |
109 | void AliAnalysisManager::Init(TTree */*tree*/) |
110 | { |
111 | // The Init() function is called when the selector needs to initialize |
112 | // a new tree or chain. Typically here the branch addresses of the tree |
113 | // will be set. It is normaly not necessary to make changes to the |
114 | // generated code, but the routine can be extended by the user if needed. |
115 | // Init() will be called many times when running with PROOF. |
116 | } |
117 | |
118 | //______________________________________________________________________________ |
119 | void AliAnalysisManager::Begin(TTree */*tree*/) |
120 | { |
121 | // The Begin() function is called at the start of the query. |
122 | // When running with PROOF Begin() is only called on the client. |
123 | // The tree argument is deprecated (on PROOF 0 is passed). |
124 | } |
125 | |
126 | //______________________________________________________________________________ |
127 | void AliAnalysisManager::SlaveBegin(TTree */*tree*/) |
128 | { |
129 | // The SlaveBegin() function is called after the Begin() function. |
130 | // When running with PROOF SlaveBegin() is called on each slave server. |
131 | // The tree argument is deprecated (on PROOF 0 is passed). |
132 | } |
133 | |
134 | //______________________________________________________________________________ |
135 | Bool_t AliAnalysisManager::Process(Long64_t /*entry*/) |
136 | { |
137 | // The Process() function is called for each entry in the tree (or possibly |
138 | // keyed object in the case of PROOF) to be processed. The entry argument |
139 | // specifies which entry in the currently loaded tree is to be processed. |
140 | // It can be passed to either TTree::GetEntry() or TBranch::GetEntry() |
141 | // to read either all or the required parts of the data. When processing |
142 | // keyed objects with PROOF, the object is already loaded and is available |
143 | // via the fObject pointer. |
144 | // |
145 | // This function should contain the "body" of the analysis. It can contain |
146 | // simple or elaborate selection criteria, run algorithms on the data |
147 | // of the event and typically fill histograms. |
148 | |
149 | // WARNING when a selector is used with a TChain, you must use |
150 | // the pointer to the current TTree to call GetEntry(entry). |
151 | // The entry is always the local entry number in the current tree. |
152 | // Assuming that fChain is the pointer to the TChain being processed, |
153 | // use fChain->GetTree()->GetEntry(entry). |
154 | return kFALSE; |
155 | } |
156 | |
157 | //______________________________________________________________________________ |
158 | void AliAnalysisManager::SlaveTerminate() |
159 | { |
160 | // The SlaveTerminate() function is called after all entries or objects |
161 | // have been processed. When running with PROOF SlaveTerminate() is called |
162 | // on each slave server. |
163 | } |
164 | |
165 | //______________________________________________________________________________ |
166 | void AliAnalysisManager::Terminate() |
167 | { |
168 | // The Terminate() function is the last function to be called during |
169 | // a query. It always runs on the client, it can be used to present |
170 | // the results graphically or save the results to file. |
171 | } |
172 | |
173 | //______________________________________________________________________________ |
174 | void AliAnalysisManager::AddTask(AliAnalysisTask *task) |
175 | { |
176 | // Adds a user task to the global list of tasks. |
177 | task->SetActive(kFALSE); |
178 | fTasks->Add(task); |
179 | } |
180 | |
181 | //______________________________________________________________________________ |
182 | AliAnalysisTask *AliAnalysisManager::GetTask(const char *name) const |
183 | { |
184 | // Retreive task by name. |
185 | if (!fTasks) return NULL; |
186 | return (AliAnalysisTask*)fTasks->FindObject(name); |
187 | } |
188 | |
189 | //______________________________________________________________________________ |
190 | AliAnalysisDataContainer *AliAnalysisManager::CreateContainer(const char *name, |
191 | TClass *datatype, EAliAnalysisContType type) |
192 | { |
193 | // Create a data container of a certain type. Types can be: |
194 | // kNormalContainer = 0, used to exchange date between tasks |
195 | // kInputContainer = 1, used to store input data |
196 | // kOutputContainer = 2, used for posting results |
197 | AliAnalysisDataContainer *cont = new AliAnalysisDataContainer(name, datatype); |
198 | fContainers->Add(cont); |
199 | switch (type) { |
200 | case kInputContainer: |
201 | fInputs->Add(cont); |
202 | break; |
203 | case kOutputContainer: |
204 | fOutputs->Add(cont); |
205 | break; |
206 | case kNormalContainer: |
207 | break; |
208 | } |
209 | return cont; |
210 | } |
211 | |
212 | //______________________________________________________________________________ |
213 | Bool_t AliAnalysisManager::ConnectInput(AliAnalysisTask *task, Int_t islot, |
214 | AliAnalysisDataContainer *cont) |
215 | { |
216 | // Connect input of an existing task to a data container. |
217 | if (!fTasks->FindObject(task)) { |
218 | AddTask(task); |
219 | AliInfo(Form("Task %s not registered. Now owned by analysis manager", task->GetName())); |
220 | } |
221 | Bool_t connected = task->ConnectInput(islot, cont); |
222 | return connected; |
223 | } |
224 | |
225 | //______________________________________________________________________________ |
226 | Bool_t AliAnalysisManager::ConnectOutput(AliAnalysisTask *task, Int_t islot, |
227 | AliAnalysisDataContainer *cont) |
228 | { |
229 | // Connect output of an existing task to a data container. |
230 | if (!fTasks->FindObject(task)) { |
231 | AddTask(task); |
232 | AliInfo(Form("Task %s not registered. Now owned by analysis manager", task->GetName())); |
233 | } |
234 | Bool_t connected = task->ConnectOutput(islot, cont); |
235 | return connected; |
236 | } |
237 | |
238 | //______________________________________________________________________________ |
239 | void AliAnalysisManager::CleanContainers() |
240 | { |
241 | // Clean data from all containers that have already finished all client tasks. |
242 | TIter next(fContainers); |
243 | AliAnalysisDataContainer *cont; |
244 | while ((cont=(AliAnalysisDataContainer *)next())) { |
245 | if (cont->IsOwnedData() && |
246 | cont->IsDataReady() && |
247 | cont->ClientsExecuted()) cont->DeleteData(); |
248 | } |
249 | } |
250 | |
251 | //______________________________________________________________________________ |
252 | Bool_t AliAnalysisManager::InitAnalysis() |
253 | { |
254 | // Initialization of analysis chain of tasks. Should be called after all tasks |
255 | // and data containers are properly connected |
256 | // Check for input/output containers |
257 | fInitOK = kFALSE; |
258 | if (!fInputs->GetEntriesFast()) { |
259 | AliError("No input container defined. At least one container should store input data"); |
260 | return kFALSE; |
261 | } |
262 | if (!fOutputs->GetEntriesFast()) { |
263 | AliError("No output container defined. At least one container should store output data"); |
264 | return kFALSE; |
265 | } |
266 | // Check for top tasks (depending only on input data containers) |
267 | if (!fTasks->First()) { |
268 | AliError("Analysis have no tasks !"); |
269 | return kFALSE; |
270 | } |
271 | TIter next(fTasks); |
272 | AliAnalysisTask *task; |
273 | AliAnalysisDataContainer *cont; |
274 | Int_t ntop = 0; |
275 | Int_t nzombies = 0; |
276 | Bool_t is_zombie = kFALSE; |
277 | Bool_t is_top = kTRUE; |
278 | Int_t i; |
279 | while ((task=(AliAnalysisTask*)next())) { |
280 | is_top = kTRUE; |
281 | is_zombie = kFALSE; |
282 | Int_t ninputs = task->GetNinputs(); |
283 | if (!ninputs) { |
284 | task->SetZombie(); |
285 | fZombies->Add(task); |
286 | nzombies++; |
287 | AliWarning(Form("Task %s has no input slots defined ! Declared zombie...",task->GetName())); |
288 | continue; |
289 | } |
290 | for (i=0; i<ninputs; i++) { |
291 | cont = task->GetInputSlot(i)->GetContainer(); |
292 | if (!cont) { |
293 | if (!is_zombie) { |
294 | task->SetZombie(); |
295 | fZombies->Add(task); |
296 | nzombies++; |
297 | is_zombie = kTRUE; |
298 | } |
299 | AliWarning(Form("Input slot %i of task %s has no container connected ! Declared zombie...", |
300 | i,task->GetName())); |
301 | } |
302 | if (is_zombie) continue; |
303 | // Check if cont is an input container |
304 | if (is_top && !fInputs->FindObject(cont)) is_top=kFALSE; |
305 | // Connect to parent task |
306 | } |
307 | if (is_top) { |
308 | ntop++; |
309 | fTopTasks->Add(task); |
310 | } |
311 | } |
312 | if (!ntop) { |
313 | AliError("No top task defined. At least one task should be connected only to input containers"); |
314 | return kFALSE; |
315 | } |
316 | // Check now if there are orphan tasks |
317 | for (i=0; i<ntop; i++) { |
318 | task = (AliAnalysisTask*)fTopTasks->At(i); |
319 | task->SetUsed(); |
320 | } |
321 | Int_t norphans = 0; |
322 | next.Reset(); |
323 | while ((task=(AliAnalysisTask*)next())) { |
324 | if (!task->IsUsed()) { |
325 | norphans++; |
326 | AliWarning(Form("Task %s is orphan",task->GetName())); |
327 | } |
328 | } |
329 | // Check the task hierarchy (no parent task should depend on data provided |
330 | // by a daughter task) |
331 | for (i=0; i<ntop; i++) { |
332 | task = (AliAnalysisTask*)fTopTasks->At(i); |
333 | if (task->CheckCircularDeps()) { |
334 | AliError("Found illegal circular dependencies between following tasks:"); |
335 | PrintStatus("dep"); |
336 | return kFALSE; |
337 | } |
338 | } |
339 | return kTRUE; |
340 | } |
341 | |
342 | //______________________________________________________________________________ |
343 | void AliAnalysisManager::PrintStatus(Option_t *option) const |
344 | { |
345 | // Print task hierarchy. |
346 | TIter next(fTopTasks); |
347 | AliAnalysisTask *task; |
348 | while ((task=(AliAnalysisTask*)next())) |
349 | task->PrintTask(option); |
350 | } |
351 | |
352 | //______________________________________________________________________________ |
353 | void AliAnalysisManager::ResetAnalysis() |
354 | { |
355 | // Reset all execution flags and clean containers. |
356 | CleanContainers(); |
357 | } |
358 | |
359 | //______________________________________________________________________________ |
360 | void AliAnalysisManager::ExecAnalysis(Option_t *option) |
361 | { |
362 | // Execute analysis. |
363 | TIter next(fTopTasks); |
364 | AliAnalysisTask *task; |
365 | while ((task=(AliAnalysisTask*)next())) |
366 | task->ExecuteTask(option); |
367 | } |
368 | |
369 | //______________________________________________________________________________ |
370 | void AliAnalysisManager::FinishAnalysis() |
371 | { |
372 | // Finish analysis. |
373 | } |