1 | /* ======================================================================== *\
|
---|
2 | !
|
---|
3 | ! *
|
---|
4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
---|
5 | ! * Software. It is distributed to you in the hope that it can be a useful
|
---|
6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
---|
7 | ! * It is distributed WITHOUT ANY WARRANTY.
|
---|
8 | ! *
|
---|
9 | ! * Permission to use, copy, modify and distribute this software and its
|
---|
10 | ! * documentation for any purpose is hereby granted without fee,
|
---|
11 | ! * provided that the above copyright notice appear in all copies and
|
---|
12 | ! * that both that copyright notice and this permission notice appear
|
---|
13 | ! * in supporting documentation. It is provided "as is" without express
|
---|
14 | ! * or implied warranty.
|
---|
15 | ! *
|
---|
16 | !
|
---|
17 | !
|
---|
18 | ! Author(s): Thomas Bretz 12/2000 (tbretz@uni-sw.gwdg.de)
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2001
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | /////////////////////////////////////////////////////////////////////////////
|
---|
26 | // //
|
---|
27 | // MReadTree //
|
---|
28 | // //
|
---|
29 | // This tasks opens all branches in a specified tree and creates the //
|
---|
30 | // corresponding parameter containers if not already existing in the //
|
---|
31 | // parameter list. //
|
---|
32 | // //
|
---|
33 | // The Process function reads one events from the tree. To go through the //
|
---|
34 | // events of one tree make sure that the event number is increased from //
|
---|
35 | // outside. It makes also possible to go back by decreasing the number. //
|
---|
36 | // //
|
---|
37 | // If you don't want to start reading the first event you have to call //
|
---|
38 | // MReadTree::SetEventNum after instantiating your MReadTree-object. //
|
---|
39 | // //
|
---|
40 | // To make reading much faster (up to a factor of 10 to 20) you can //
|
---|
41 | // ensure that only the data you are really processing is enabled by //
|
---|
42 | // calling MReadTree::UseLeaf. //
|
---|
43 | // //
|
---|
44 | // Later we'll use TChain::SetNotify to notify MReadTree if the TChain //
|
---|
45 | // starts to read a new file. //
|
---|
46 | // //
|
---|
47 | /////////////////////////////////////////////////////////////////////////////
|
---|
48 |
|
---|
49 | #include "MReadTree.h"
|
---|
50 |
|
---|
51 | #include <fstream.h>
|
---|
52 |
|
---|
53 | #include <TFile.h>
|
---|
54 | #include <TChain.h>
|
---|
55 | #include <TArrayC.h>
|
---|
56 | #include <TObjArray.h>
|
---|
57 | #include <TChainElement.h>
|
---|
58 |
|
---|
59 | #include "MLog.h"
|
---|
60 | #include "MLogManip.h"
|
---|
61 |
|
---|
62 | #include "MTime.h"
|
---|
63 | #include "MParList.h"
|
---|
64 |
|
---|
65 | ClassImp(MReadTree);
|
---|
66 |
|
---|
67 | // --------------------------------------------------------------------------
|
---|
68 | //
|
---|
69 | // Default constructor. It creates an TChain instance which represents the
|
---|
70 | // the Tree you want to read and adds the given file (if you gave one).
|
---|
71 | // More files can be added using MReadTree::AddFile.
|
---|
72 | // Also an empty veto list is created. This list is used if you want to
|
---|
73 | // veto (disable or "don't enable") a branch in the tree.
|
---|
74 | //
|
---|
75 | MReadTree::MReadTree(const char *tname, const char *fname,
|
---|
76 | const char *name, const char *title)
|
---|
77 | : fNumEntry(0), fLeafEnabled(kFALSE)
|
---|
78 | {
|
---|
79 | *fName = name ? name : "MReadTree";
|
---|
80 | *fTitle = title ? title : "Task to loop over all events in one single tree";
|
---|
81 |
|
---|
82 | fVetoList = new TArrayC;
|
---|
83 | //
|
---|
84 | // open the input stream
|
---|
85 | //
|
---|
86 | fChain = new TChain(tname);
|
---|
87 |
|
---|
88 | if (fname)
|
---|
89 | fChain->Add(fname);
|
---|
90 | }
|
---|
91 |
|
---|
92 | // --------------------------------------------------------------------------
|
---|
93 | //
|
---|
94 | // Destructor. It deletes the TChain and veto list object
|
---|
95 | //
|
---|
96 | MReadTree::~MReadTree()
|
---|
97 | {
|
---|
98 | //
|
---|
99 | // Delete all the pointers to pointers to the objects where the
|
---|
100 | // branche data gets stored.
|
---|
101 | //
|
---|
102 | TIter Next(fChain->GetStatus());
|
---|
103 |
|
---|
104 | TChainElement *element = NULL;
|
---|
105 | while ((element=(TChainElement*)Next()))
|
---|
106 | delete (MParContainer**)element->GetBaddress();
|
---|
107 |
|
---|
108 | //
|
---|
109 | // Delete the chain and the veto list
|
---|
110 | //
|
---|
111 | delete fChain;
|
---|
112 | delete fVetoList;
|
---|
113 | }
|
---|
114 |
|
---|
115 | // --------------------------------------------------------------------------
|
---|
116 | //
|
---|
117 | // If you want to read the given tree over several files you must add
|
---|
118 | // the files here before PreProcess is called. Be careful: If the tree
|
---|
119 | // doesn't have the same contents (branches) it may confuse your
|
---|
120 | // program (trees which are are not existing in later files are not read
|
---|
121 | // anymore, tree wich are not existing in the first file are never read)
|
---|
122 | //
|
---|
123 | // Name may use the wildcarding notation, eg "xxx*.root" means all files
|
---|
124 | // starting with xxx in the current file system directory.
|
---|
125 | //
|
---|
126 | Int_t MReadTree::AddFile(const char *fname)
|
---|
127 | {
|
---|
128 | //
|
---|
129 | // FIXME! A check is missing whether the file already exists or not.
|
---|
130 | //
|
---|
131 | //
|
---|
132 | // returns the number of file which were added
|
---|
133 | //
|
---|
134 | return fChain->Add(fname);
|
---|
135 | }
|
---|
136 |
|
---|
137 | // --------------------------------------------------------------------------
|
---|
138 | //
|
---|
139 | // The first time this function is called all leafes/branches are disabled
|
---|
140 | // and the given branch/leaf is enabled. By enabling only the branches you
|
---|
141 | // are processing you can speed up your calculation many times (up to
|
---|
142 | // a factor of 10 or 20)
|
---|
143 | //
|
---|
144 | void MReadTree::UseLeaf(const char *name)
|
---|
145 | {
|
---|
146 | if (!fLeafEnabled)
|
---|
147 | {
|
---|
148 | *fLog << "Leaf choosing method enabled (only enabled leaves are read)." << endl;
|
---|
149 | fChain->SetBranchStatus("*", kFALSE);
|
---|
150 | fLeafEnabled = kTRUE;
|
---|
151 | }
|
---|
152 |
|
---|
153 | fChain->SetBranchStatus(name, kTRUE);
|
---|
154 | }
|
---|
155 |
|
---|
156 | // --------------------------------------------------------------------------
|
---|
157 | //
|
---|
158 | // The PreProcess loops (till now) over the branches in the given tree.
|
---|
159 | // It checks if the corresponding containers (containers with the same
|
---|
160 | // name than the branch name) are existing in the Parameter Container List.
|
---|
161 | // If not, a container of objec type 'branch-name' is created (everything
|
---|
162 | // after the last semicolon in the branch name is stripped). Only
|
---|
163 | // branches which don't have a veto (see VetoBranch) are enabled If the
|
---|
164 | // object isn't found in the root dictionary (a list of classes known by the
|
---|
165 | // root environment) the branch is skipped and an error message is printed
|
---|
166 | // out.
|
---|
167 | //
|
---|
168 | Bool_t MReadTree::PreProcess(MParList *pList)
|
---|
169 | {
|
---|
170 | //
|
---|
171 | // get number of events in this tree
|
---|
172 | //
|
---|
173 | fNumEntries = (UInt_t)fChain->GetEntries();
|
---|
174 |
|
---|
175 | if (!fNumEntries)
|
---|
176 | {
|
---|
177 | *fLog << dbginf << "No entries found in chain (file/s)." << endl;
|
---|
178 | return kFALSE;
|
---|
179 | }
|
---|
180 |
|
---|
181 | //
|
---|
182 | // output logging information
|
---|
183 | //
|
---|
184 | *fLog << fNumEntries << " entries found in file(s)." << endl;
|
---|
185 |
|
---|
186 | //
|
---|
187 | // Get all branches of this tree and
|
---|
188 | // create the Iterator to loop over all branches
|
---|
189 | //
|
---|
190 | TIter Next(fChain->GetListOfBranches());
|
---|
191 | TBranch *branch=NULL;
|
---|
192 |
|
---|
193 | Int_t num=0;
|
---|
194 | //
|
---|
195 | // loop over all tasks for processing
|
---|
196 | //
|
---|
197 | while ( (branch=(TBranch*)Next()) )
|
---|
198 | {
|
---|
199 | //
|
---|
200 | // Get Name of Branch
|
---|
201 | //
|
---|
202 | const char *name = branch->GetName();
|
---|
203 |
|
---|
204 | //
|
---|
205 | // Check if enabeling the branch is allowed
|
---|
206 | //
|
---|
207 | if (HasVeto(name))
|
---|
208 | {
|
---|
209 | *fLog << "Branch " << name << " has veto... skipped." << endl;
|
---|
210 | continue;
|
---|
211 | }
|
---|
212 |
|
---|
213 | //
|
---|
214 | // Create a pointer to the pointer to the object in which the
|
---|
215 | // branch data is stored. The pointers are stored in the TChain
|
---|
216 | // object and we get the pointers from there to delete it.
|
---|
217 | //
|
---|
218 | MParContainer **pcont = new MParContainer*;
|
---|
219 |
|
---|
220 | //
|
---|
221 | // check if object is existing in the list
|
---|
222 | //
|
---|
223 | *pcont=pList->FindCreateObj(name);
|
---|
224 |
|
---|
225 | if (!*pcont)
|
---|
226 | {
|
---|
227 | //
|
---|
228 | // if class is not existing in the (root) environment
|
---|
229 | // we cannot proceed reading this branch
|
---|
230 | //
|
---|
231 | *fLog << dbginf << "Warning: Class '" << name << "' not existing in dictionary. Branch skipped." << endl;
|
---|
232 | continue;
|
---|
233 | }
|
---|
234 |
|
---|
235 | //
|
---|
236 | // here pcont is a pointer the to container in which the data from
|
---|
237 | // the actual branch should be stored - enable branch.
|
---|
238 | //
|
---|
239 | fChain->SetBranchAddress(name, pcont);
|
---|
240 |
|
---|
241 | *fLog << "Branch " << name << " enabled for reading." << endl;
|
---|
242 |
|
---|
243 | num++;
|
---|
244 | }
|
---|
245 |
|
---|
246 | *fLog << "MReadTree will read " << num << " branches." << endl;
|
---|
247 |
|
---|
248 | return kTRUE;
|
---|
249 | }
|
---|
250 |
|
---|
251 | // --------------------------------------------------------------------------
|
---|
252 | //
|
---|
253 | // The Process-function reads one event from the tree (this contains all
|
---|
254 | // enabled branches) and increases the position in the file by one event.
|
---|
255 | // (Remark: The position can also be set by some member functions
|
---|
256 | // If the end of the file is reached the Eventloop is stopped.
|
---|
257 | //
|
---|
258 | Bool_t MReadTree::Process()
|
---|
259 | {
|
---|
260 | return fChain->GetEntry(fNumEntry++, 0) == 0 ? kFALSE : kTRUE;
|
---|
261 | }
|
---|
262 |
|
---|
263 | // --------------------------------------------------------------------------
|
---|
264 | //
|
---|
265 | // Get the Event with the current EventNumber fNumEntry
|
---|
266 | //
|
---|
267 | Bool_t MReadTree::GetEvent()
|
---|
268 | {
|
---|
269 | return fChain->GetEntry(fNumEntry) == 0 ? kFALSE : kTRUE;
|
---|
270 | }
|
---|
271 |
|
---|
272 | // --------------------------------------------------------------------------
|
---|
273 | //
|
---|
274 | // Decrease the number of the event which is read by Process() next
|
---|
275 | // by one or more
|
---|
276 | //
|
---|
277 | Bool_t MReadTree::DecEventNum(UInt_t dec)
|
---|
278 | {
|
---|
279 | //!
|
---|
280 | //! this function makes Process() read the event one (or more) before
|
---|
281 | //! the actual position (event) in the tree
|
---|
282 | //!
|
---|
283 | if (fNumEntry < dec/*+1*/)
|
---|
284 | {
|
---|
285 | *fLog << "MReadTree::SetPrevEvent: WARNING: " << fNumEntry/*-1*/ << "-" << dec << " out of Range." << endl;
|
---|
286 | return kFALSE;
|
---|
287 | }
|
---|
288 |
|
---|
289 | fNumEntry -= dec/*+1*/;
|
---|
290 | return kTRUE;
|
---|
291 | }
|
---|
292 |
|
---|
293 | // --------------------------------------------------------------------------
|
---|
294 | //
|
---|
295 | // Increase the number of the event which is read by Process() next
|
---|
296 | // by one or more
|
---|
297 | //
|
---|
298 | Bool_t MReadTree::IncEventNum(UInt_t inc)
|
---|
299 | {
|
---|
300 | //!
|
---|
301 | //! this function makes Process() read the next (or more) after
|
---|
302 | //! the actual position (event) in the tree
|
---|
303 | //! (Be careful: IncEventNum() or IncEventNum(1) does not chenge anything
|
---|
304 | //! in the standard behaviour of the task)
|
---|
305 | //!
|
---|
306 | if (fNumEntry+inc/*-1*/ >= fNumEntries)
|
---|
307 | {
|
---|
308 | *fLog << "MReadTree::SkipEvents: WARNING: " << fNumEntry/*-1*/ << "+" << inc << " out of Range." << endl;
|
---|
309 | return kFALSE;
|
---|
310 | }
|
---|
311 |
|
---|
312 | fNumEntry += inc/*-1*/;
|
---|
313 | return kTRUE;
|
---|
314 | }
|
---|
315 |
|
---|
316 | // --------------------------------------------------------------------------
|
---|
317 | //
|
---|
318 | // this function makes Process() read event number nr next
|
---|
319 | //
|
---|
320 | // Remark: You can use this function after instatiating you MReadTree-object
|
---|
321 | // to set the event number from which you want to start reading.
|
---|
322 | //
|
---|
323 | Bool_t MReadTree::SetEventNum(UInt_t nr)
|
---|
324 | {
|
---|
325 | if (nr>=fNumEntries)
|
---|
326 | {
|
---|
327 | *fLog << "MReadTree::SetEventNum: WARNING: " << nr << " out of Range." << endl;
|
---|
328 | return kFALSE;
|
---|
329 | }
|
---|
330 |
|
---|
331 | fNumEntry = nr;
|
---|
332 | return kTRUE;
|
---|
333 | }
|
---|
334 |
|
---|
335 | // --------------------------------------------------------------------------
|
---|
336 | //
|
---|
337 | // This function checks if the Branch Name is part of the Veto List.
|
---|
338 | // This means, that the preprocess doesn't enable the branch.
|
---|
339 | //
|
---|
340 | Bool_t MReadTree::HasVeto(const char *name) const
|
---|
341 | {
|
---|
342 | const size_t nlen = strlen(name)+1;
|
---|
343 |
|
---|
344 | char *pos = fVetoList->GetArray();
|
---|
345 | size_t len = fVetoList->GetSize();
|
---|
346 |
|
---|
347 | //
|
---|
348 | // now we compare the 'strings' in the list word by word
|
---|
349 | // (string or word in this context means a zero terminated
|
---|
350 | // array of chars
|
---|
351 | //
|
---|
352 |
|
---|
353 | for (;;)
|
---|
354 | {
|
---|
355 | //
|
---|
356 | // Search for the first byte of the name
|
---|
357 | //
|
---|
358 | char *c = (char*)memchr(pos, name[0], len);
|
---|
359 |
|
---|
360 | //
|
---|
361 | // if we don't find the first byte, the list cannot contain 'name'
|
---|
362 | //
|
---|
363 | if (!c)
|
---|
364 | return kFALSE;
|
---|
365 |
|
---|
366 | //
|
---|
367 | // calculate and check how many bytes remains in the list
|
---|
368 | //
|
---|
369 | len -= c-pos;
|
---|
370 |
|
---|
371 | //
|
---|
372 | // if there are not enough bytes to match the query we are done
|
---|
373 | //
|
---|
374 | if (len<nlen)
|
---|
375 | return kFALSE;
|
---|
376 |
|
---|
377 | //
|
---|
378 | // check if the next 'nlen' byte (including the trailing '\0'
|
---|
379 | // are matching
|
---|
380 | //
|
---|
381 | if (!memcmp(c, name, nlen))
|
---|
382 | return kTRUE;
|
---|
383 |
|
---|
384 | //
|
---|
385 | // we didn't find the string, goto the next 'string'
|
---|
386 | //
|
---|
387 | pos = (char*)memchr(c, '\0', len);
|
---|
388 |
|
---|
389 | //
|
---|
390 | // check if there is a 'next' string really
|
---|
391 | //
|
---|
392 | if (!pos)
|
---|
393 | return kFALSE;
|
---|
394 |
|
---|
395 | //
|
---|
396 | // calculate the remaining length
|
---|
397 | //
|
---|
398 | len -= pos-c;
|
---|
399 |
|
---|
400 | //
|
---|
401 | // if there are not enough bytes to match the query we are done
|
---|
402 | //
|
---|
403 | if (len<nlen)
|
---|
404 | return kFALSE;
|
---|
405 | }
|
---|
406 | }
|
---|
407 |
|
---|
408 |
|
---|
409 | // --------------------------------------------------------------------------
|
---|
410 | //
|
---|
411 | // If you don't want that a branch is enabled within the PreProcess you
|
---|
412 | // can set a veto for enabeling the branch. (This means also the
|
---|
413 | // corresponding object won't be created automatically)
|
---|
414 | //
|
---|
415 | // This functionality is for experienced users which don't want to
|
---|
416 | // read in branches which are not processed in the program (for
|
---|
417 | // speed reasons)
|
---|
418 | //
|
---|
419 | void MReadTree::VetoBranch(const char *name)
|
---|
420 | {
|
---|
421 | //
|
---|
422 | // Add this file as the last entry of the list
|
---|
423 | // (including the trailing '\0')
|
---|
424 | //
|
---|
425 | const int sz = fVetoList->GetSize();
|
---|
426 | const int tsz = strlen(name)+1;
|
---|
427 |
|
---|
428 | fVetoList->Set(sz+tsz);
|
---|
429 |
|
---|
430 | memcpy(fVetoList->GetArray()+sz, name, tsz);
|
---|
431 | }
|
---|