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

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