source: trunk/MagicSoft/Mars/mbase/MReadTree.cc@ 992

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