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