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

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