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 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2003
|
---|
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 | // If the chain switches from one file to another file all //
|
---|
45 | // TObject::Notify() functions are called of TObject objects which were //
|
---|
46 | // added to the Notifier list view MReadTree::AddNotify. If MReadTree //
|
---|
47 | // is the owner (viw MReadTree::SetOwner) all this objects are deleted //
|
---|
48 | // by the destructor of MReadTree //
|
---|
49 | // //
|
---|
50 | /////////////////////////////////////////////////////////////////////////////
|
---|
51 | #include "MReadTree.h"
|
---|
52 |
|
---|
53 | #include <fstream>
|
---|
54 |
|
---|
55 | #include <TFile.h> // TFile::GetName
|
---|
56 | #include <TSystem.h> // gSystem->ExpandPath
|
---|
57 | #include <TChainElement.h>
|
---|
58 | #include <TOrdCollection.h>
|
---|
59 |
|
---|
60 | #include "MChain.h"
|
---|
61 | #include "MFilter.h"
|
---|
62 | #include "MParList.h"
|
---|
63 | #include "MTaskList.h"
|
---|
64 |
|
---|
65 | #include "MLog.h"
|
---|
66 | #include "MLogManip.h"
|
---|
67 |
|
---|
68 |
|
---|
69 | ClassImp(MReadTree);
|
---|
70 |
|
---|
71 | using namespace std;
|
---|
72 |
|
---|
73 | // --------------------------------------------------------------------------
|
---|
74 | //
|
---|
75 | // Default constructor. Don't use it.
|
---|
76 | //
|
---|
77 | MReadTree::MReadTree()
|
---|
78 | : fNumEntry(0), fNumEntries(0), fBranchChoosing(kFALSE), fAutoEnable(kTRUE)
|
---|
79 | {
|
---|
80 | fName = "MRead";
|
---|
81 | fTitle = "Task to loop over all events in one single tree";
|
---|
82 |
|
---|
83 | fVetoList = NULL;
|
---|
84 | fNotify = NULL;
|
---|
85 |
|
---|
86 | fChain = NULL;
|
---|
87 | }
|
---|
88 |
|
---|
89 | // --------------------------------------------------------------------------
|
---|
90 | //
|
---|
91 | // Constructor. It creates an TChain instance which represents the
|
---|
92 | // the Tree you want to read and adds the given file (if you gave one).
|
---|
93 | // More files can be added using MReadTree::AddFile.
|
---|
94 | // Also an empty veto list is created. This list is used if you want to
|
---|
95 | // veto (disable or "don't enable") a branch in the tree, it vetos also
|
---|
96 | // the creation of the corresponding object.
|
---|
97 | // An empty list of TObjects are also created. This objects are called
|
---|
98 | // at any time the TChain starts to read from another file.
|
---|
99 | //
|
---|
100 | MReadTree::MReadTree(const char *tname, const char *fname,
|
---|
101 | const char *name, const char *title)
|
---|
102 | : fNumEntry(0), fNumEntries(0), fBranchChoosing(kFALSE), fAutoEnable(kTRUE)
|
---|
103 | {
|
---|
104 | fName = name ? name : "MRead";
|
---|
105 | fTitle = title ? title : "Task to loop over all events in one single tree";
|
---|
106 |
|
---|
107 | fVetoList = new TList;
|
---|
108 | fVetoList->SetOwner();
|
---|
109 |
|
---|
110 | fNotify = new TList;
|
---|
111 |
|
---|
112 | //
|
---|
113 | // open the input stream
|
---|
114 | //
|
---|
115 | fChain = new MChain(tname);
|
---|
116 |
|
---|
117 | // root 3.02:
|
---|
118 | // In TChain::Addfile remove the limitation that the file name must contain
|
---|
119 | // the string ".root". ".root" is necessary only in case one wants to specify
|
---|
120 | // a Tree in a subdirectory of a Root file with eg, the format:
|
---|
121 |
|
---|
122 | if (fname)
|
---|
123 | AddFile(fname);
|
---|
124 | // if (fChain->Add(fname)>0)
|
---|
125 | // SetBit(kChainWasChanged);
|
---|
126 | }
|
---|
127 |
|
---|
128 | // --------------------------------------------------------------------------
|
---|
129 | //
|
---|
130 | // Destructor. It deletes the TChain and veto list object
|
---|
131 | //
|
---|
132 | MReadTree::~MReadTree()
|
---|
133 | {
|
---|
134 | //
|
---|
135 | // Delete all the pointers to pointers to the objects where the
|
---|
136 | // branche data gets stored. FIXME: When PreProcessed twice this
|
---|
137 | // creates a memory leak!
|
---|
138 | //
|
---|
139 | TIter Next(fChain->GetStatus());
|
---|
140 |
|
---|
141 | TChainElement *element = NULL;
|
---|
142 | while ((element=(TChainElement*)Next()))
|
---|
143 | if (element->GetBaddress())
|
---|
144 | delete (MParContainer**)element->GetBaddress();
|
---|
145 |
|
---|
146 | //
|
---|
147 | // Delete the chain and the veto list
|
---|
148 | //
|
---|
149 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,03,00)
|
---|
150 | if (fChain->GetFile())
|
---|
151 | delete fChain->GetFile();
|
---|
152 | #endif
|
---|
153 | delete fChain;
|
---|
154 |
|
---|
155 | delete fNotify;
|
---|
156 | delete fVetoList;
|
---|
157 | }
|
---|
158 |
|
---|
159 | // --------------------------------------------------------------------------
|
---|
160 | //
|
---|
161 | // This check whether all branches in the tree have the same size. If
|
---|
162 | // this is not the case the events in the different branches cannot
|
---|
163 | // be ordered correctly.
|
---|
164 | //
|
---|
165 | Bool_t MReadTree::CheckBranchSize()
|
---|
166 | {
|
---|
167 | TArrayI entries(fChain->GetStatus()->GetSize());
|
---|
168 | Int_t num=0;
|
---|
169 |
|
---|
170 | // Loop over all branches which have a corresponding container
|
---|
171 | TIter Next(fChain->GetStatus());
|
---|
172 |
|
---|
173 | TChainElement *element = NULL;
|
---|
174 | while ((element=(TChainElement*)Next()))
|
---|
175 | {
|
---|
176 | // Get branch name and find pointer to corresponding branch
|
---|
177 | const TString name = element->GetName();
|
---|
178 | const TBranch *b = fChain->FindBranch(name);
|
---|
179 |
|
---|
180 | // Skip element without corresponding branches (like "*")
|
---|
181 | if (!b)
|
---|
182 | continue;
|
---|
183 |
|
---|
184 | entries[num++] = (Int_t)b->GetEntries();
|
---|
185 | }
|
---|
186 |
|
---|
187 | // Check the number of entries of the branches pair-wise
|
---|
188 | for (int i=0; i<num; i++)
|
---|
189 | for (int j=i; j<num; j++)
|
---|
190 | {
|
---|
191 | if (entries[i]==entries[j])
|
---|
192 | continue;
|
---|
193 |
|
---|
194 | *fLog << err << "ERROR - File corruption detected:" << endl;
|
---|
195 | *fLog << " Due to several circumstances (such at a bug in MReadTree or wrong" << endl;
|
---|
196 | *fLog << " usage of the file UPDATE mode) you may have produced a file in which" << endl;
|
---|
197 | *fLog << " at least two branches in the same tree (" << fChain->GetName() << ") have different" << endl;
|
---|
198 | *fLog << " number of entries. Sorry, but this file (" << GetFileName() << ")" << endl;
|
---|
199 | *fLog << " is unusable." << endl;
|
---|
200 | return kFALSE;
|
---|
201 | }
|
---|
202 |
|
---|
203 | return kTRUE;
|
---|
204 | }
|
---|
205 |
|
---|
206 | // --------------------------------------------------------------------------
|
---|
207 | //
|
---|
208 | // If the owner flag is set all TObjects which are scheduled via
|
---|
209 | // AddNotify are deleted by the destructor of MReadTree
|
---|
210 | //
|
---|
211 | void MReadTree::SetOwner(Bool_t flag)
|
---|
212 | {
|
---|
213 | flag ? fNotify->SetBit(kIsOwner) : fNotify->ResetBit(kIsOwner);
|
---|
214 | }
|
---|
215 |
|
---|
216 | // --------------------------------------------------------------------------
|
---|
217 | //
|
---|
218 | // This function is called each time MReadTree changes the file to read
|
---|
219 | // from. It calls all TObject::Notify() functions which are scheduled
|
---|
220 | // via AddNotify.
|
---|
221 | //
|
---|
222 | Bool_t MReadTree::Notify()
|
---|
223 | {
|
---|
224 | //
|
---|
225 | // Do a consistency check for all branches
|
---|
226 | //
|
---|
227 | if (!CheckBranchSize())
|
---|
228 | return kFALSE;
|
---|
229 |
|
---|
230 | *fLog << inf << GetDescriptor() << ": Switching to #" << GetFileIndex();
|
---|
231 | *fLog << " '" << GetFileName() << "' (before event #";
|
---|
232 | *fLog << GetNumEntry()-1 << ")" << endl;
|
---|
233 |
|
---|
234 | if (!fNotify)
|
---|
235 | return kTRUE;
|
---|
236 |
|
---|
237 | TIter Next(fNotify);
|
---|
238 | TObject *o=NULL;
|
---|
239 | while ((o=Next()))
|
---|
240 | if (!o->Notify())
|
---|
241 | {
|
---|
242 | *fLog << err << "Calling Notify() for object " << o->GetName() << " failed... abort." << endl;
|
---|
243 | return kFALSE;
|
---|
244 | }
|
---|
245 |
|
---|
246 | return kTRUE;
|
---|
247 | }
|
---|
248 |
|
---|
249 | // --------------------------------------------------------------------------
|
---|
250 | //
|
---|
251 | // If you want to read the given tree over several files you must add
|
---|
252 | // the files here before PreProcess is called. Be careful: If the tree
|
---|
253 | // doesn't have the same contents (branches) it may confuse your
|
---|
254 | // program (trees which are are not existing in later files are not read
|
---|
255 | // anymore, tree wich are not existing in the first file are never read)
|
---|
256 | //
|
---|
257 | // Name may use the wildcarding notation, eg "xxx*.root" means all files
|
---|
258 | // starting with xxx in the current file system directory.
|
---|
259 | //
|
---|
260 | // AddFile returns the number of files added to the chain.
|
---|
261 | //
|
---|
262 | // For more information see TChain::Add
|
---|
263 | //
|
---|
264 | Int_t MReadTree::AddFile(const char *fname, Int_t entries)
|
---|
265 | {
|
---|
266 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,03,01)
|
---|
267 | //
|
---|
268 | // This is a workaround to get rid of crashed if the file doesn't
|
---|
269 | // exist due to non initialized TFile::fProcessIDs
|
---|
270 | //
|
---|
271 | // (Code taken from TFile::TFile
|
---|
272 | //
|
---|
273 | TString newname; // char-array must overcome comming block
|
---|
274 |
|
---|
275 | if (strrchr(fname, '?') || strrchr(fname, '*'))
|
---|
276 | {
|
---|
277 | *fLog << warn;
|
---|
278 | *fLog<< "WARNING: Using widcards with older root versions:" << endl;
|
---|
279 | *fLog << " You may encounter crashes closing the files..." << endl;
|
---|
280 | }
|
---|
281 | else
|
---|
282 | {
|
---|
283 | const char *name;
|
---|
284 |
|
---|
285 | if ((name = gSystem->ExpandPathName(fname)))
|
---|
286 | {
|
---|
287 | newname = name;
|
---|
288 | delete [] name;
|
---|
289 | }
|
---|
290 |
|
---|
291 | if (newname.IsNull())
|
---|
292 | {
|
---|
293 | *fLog << err << dbginf << "Error expanding path " << fname << "." << endl;
|
---|
294 | return 0;
|
---|
295 | }
|
---|
296 |
|
---|
297 | if (gSystem->AccessPathName(newname, kFileExists))
|
---|
298 | {
|
---|
299 | *fLog << err << "ERROR - File '" << fname << "' does not exist." << endl;
|
---|
300 | return 0;
|
---|
301 | }
|
---|
302 |
|
---|
303 | fname = newname.Data();
|
---|
304 | }
|
---|
305 | #endif
|
---|
306 |
|
---|
307 | //
|
---|
308 | // FIXME! A check is missing whether the file already exists or not.
|
---|
309 | //
|
---|
310 | const Int_t numfiles = fChain->Add(fname, entries<0?TChain::kBigNumber:entries);
|
---|
311 |
|
---|
312 | if (numfiles>0)
|
---|
313 | SetBit(kChainWasChanged);
|
---|
314 | else
|
---|
315 | *fLog << warn << "WARNING: '" << fname << "' not added to " << GetDescriptor() << endl;
|
---|
316 |
|
---|
317 | return numfiles;
|
---|
318 | }
|
---|
319 |
|
---|
320 | /*
|
---|
321 | // --------------------------------------------------------------------------
|
---|
322 | //
|
---|
323 | //
|
---|
324 | Int_t MReadTree::AddFile(const TChainElement &obj)
|
---|
325 | {
|
---|
326 | return AddFile(obj.GetTitle(), obj.GetEntries());
|
---|
327 | }
|
---|
328 | */
|
---|
329 |
|
---|
330 | // --------------------------------------------------------------------------
|
---|
331 | //
|
---|
332 | // Adds all files from another MReadTree to this instance
|
---|
333 | //
|
---|
334 | // Returns the number of file which were added
|
---|
335 | //
|
---|
336 | Int_t MReadTree::AddFiles(const MReadTree &read)
|
---|
337 | {
|
---|
338 | const Int_t rc = fChain->Add(read.fChain);
|
---|
339 |
|
---|
340 | if (rc>0)
|
---|
341 | SetBit(kChainWasChanged);
|
---|
342 |
|
---|
343 | /*
|
---|
344 | Int_t rc = 0;
|
---|
345 |
|
---|
346 | TIter Next(read.fChain->GetListOfFiles());
|
---|
347 | TObject *obj = NULL;
|
---|
348 | while ((obj=Next()))
|
---|
349 | rc += AddFile(*(TChainElement*)obj);
|
---|
350 | */
|
---|
351 |
|
---|
352 | return rc;
|
---|
353 | }
|
---|
354 |
|
---|
355 | // --------------------------------------------------------------------------
|
---|
356 | //
|
---|
357 | // This function is called if Branch choosing method should get enabled.
|
---|
358 | // Branch choosing means, that only the enabled branches are read into
|
---|
359 | // memory. To use an enableing scheme we have to disable all branches first.
|
---|
360 | // This is done, if this function is called the first time.
|
---|
361 | //
|
---|
362 | void MReadTree::EnableBranchChoosing()
|
---|
363 | {
|
---|
364 | if (fBranchChoosing)
|
---|
365 | return;
|
---|
366 |
|
---|
367 | *fLog << inf << GetDescriptor() << ": Branch choosing enabled (only enabled branches are read)." << endl;
|
---|
368 | fChain->SetBranchStatus("*", kFALSE);
|
---|
369 | fBranchChoosing = kTRUE;
|
---|
370 | }
|
---|
371 |
|
---|
372 | // --------------------------------------------------------------------------
|
---|
373 | //
|
---|
374 | // The first time this function is called all branches are disabled.
|
---|
375 | // The given branch is enabled. By enabling only the branches you
|
---|
376 | // are processing you can speed up your calculation many times (up to
|
---|
377 | // a factor of 10 or 20)
|
---|
378 | //
|
---|
379 | void MReadTree::EnableBranch(const char *name)
|
---|
380 | {
|
---|
381 | if (fChain->GetListOfFiles()->GetEntries()==0)
|
---|
382 | {
|
---|
383 | *fLog << err << "Chain contains no file... Branch '";
|
---|
384 | *fLog << name << "' ignored." << endl;
|
---|
385 | return;
|
---|
386 | }
|
---|
387 |
|
---|
388 | EnableBranchChoosing();
|
---|
389 |
|
---|
390 | TNamed branch(name, "");
|
---|
391 | SetBranchStatus(&branch, kTRUE);
|
---|
392 | }
|
---|
393 |
|
---|
394 | // --------------------------------------------------------------------------
|
---|
395 | //
|
---|
396 | // Set branch status of branch name
|
---|
397 | //
|
---|
398 | void MReadTree::SetBranchStatus(const char *name, Bool_t status)
|
---|
399 | {
|
---|
400 | fChain->SetBranchStatus(name, status);
|
---|
401 |
|
---|
402 | *fLog << inf << (status ? "Enabled" : "Disabled");
|
---|
403 | *fLog << " subbranch '" << name << "'." << endl;
|
---|
404 | }
|
---|
405 |
|
---|
406 | // --------------------------------------------------------------------------
|
---|
407 | //
|
---|
408 | // Checks whether a branch with the given name exists in the chain
|
---|
409 | // and sets the branch status of this branch corresponding to status.
|
---|
410 | //
|
---|
411 | void MReadTree::SetBranchStatus(TObject *branch, Bool_t status)
|
---|
412 | {
|
---|
413 | //
|
---|
414 | // Get branch name
|
---|
415 | //
|
---|
416 | const char *name = branch->GetName();
|
---|
417 |
|
---|
418 | //
|
---|
419 | // Check whether this branch really exists
|
---|
420 | //
|
---|
421 | TString bn(name);
|
---|
422 | if (bn.EndsWith("*"))
|
---|
423 | bn.Remove(bn.Length()-1);
|
---|
424 |
|
---|
425 | if (fChain->GetBranch(bn))
|
---|
426 | SetBranchStatus(name, status);
|
---|
427 |
|
---|
428 | //
|
---|
429 | // Remove trailing '.' if one and try to enable the subbranch without
|
---|
430 | // the master branch name. This is to be compatible with older mars
|
---|
431 | // and camera files.
|
---|
432 | //
|
---|
433 | const char *dot = strrchr(name, '.');
|
---|
434 | if (!dot)
|
---|
435 | return;
|
---|
436 |
|
---|
437 | if (fChain->GetBranch(dot+1))
|
---|
438 | SetBranchStatus(dot+1, status);
|
---|
439 | }
|
---|
440 |
|
---|
441 | // --------------------------------------------------------------------------
|
---|
442 | //
|
---|
443 | // Set the status of all branches in the list to status.
|
---|
444 | //
|
---|
445 | void MReadTree::SetBranchStatus(const TList *list, Bool_t status)
|
---|
446 | {
|
---|
447 | //
|
---|
448 | // Loop over all subbranches in this master branch
|
---|
449 | //
|
---|
450 | TIter Next(list);
|
---|
451 |
|
---|
452 | TObject *obj;
|
---|
453 | while ((obj=Next()))
|
---|
454 | SetBranchStatus(obj, status);
|
---|
455 | }
|
---|
456 |
|
---|
457 | // --------------------------------------------------------------------------
|
---|
458 | //
|
---|
459 | // This is the implementation of the Auto Enabling Scheme.
|
---|
460 | // For more information see MTask::AddBranchToList.
|
---|
461 | // This function loops over all tasks and its filters in the tasklist
|
---|
462 | // and enables all branches which are requested by the tasks and its
|
---|
463 | // filters.
|
---|
464 | //
|
---|
465 | // To enable 'unknown' branches which are not in the branchlist of
|
---|
466 | // the tasks you can call EnableBranch
|
---|
467 | //
|
---|
468 | void MReadTree::EnableBranches(MParList *plist)
|
---|
469 | {
|
---|
470 | //
|
---|
471 | // check whether branch choosing must be switched on
|
---|
472 | //
|
---|
473 | EnableBranchChoosing();
|
---|
474 |
|
---|
475 | //
|
---|
476 | // request the tasklist from the parameter list.
|
---|
477 | // FIXME: Tasklist can have a different name
|
---|
478 | //
|
---|
479 | const MTaskList *tlist = (MTaskList*)plist->FindObject("MTaskList");
|
---|
480 | if (!tlist)
|
---|
481 | {
|
---|
482 | *fLog << warn << GetDescriptor() << "Cannot use auto enabeling scheme for branches. 'MTaskList' not found." << endl;
|
---|
483 | return;
|
---|
484 | }
|
---|
485 |
|
---|
486 | //
|
---|
487 | // This loop is not necessary. We could do it like in the commented
|
---|
488 | // loop below. But this loop makes sure, that we don't try to enable
|
---|
489 | // one branch several times. This would not harm, but we would get
|
---|
490 | // an output for each attempt. To have several outputs for one subbranch
|
---|
491 | // may confuse the user, this we don't want.
|
---|
492 | // This loop creates a new list of subbranches and for each branch
|
---|
493 | // which is added we check before whether it already exists or not.
|
---|
494 | //
|
---|
495 | TList list;
|
---|
496 |
|
---|
497 | MTask *task;
|
---|
498 | TIter NextTask(tlist->GetList());
|
---|
499 | while ((task=(MTask*)NextTask()))
|
---|
500 | {
|
---|
501 | TObject *obj;
|
---|
502 |
|
---|
503 | TIter NextTBranch(task->GetListOfBranches());
|
---|
504 | while ((obj=NextTBranch()))
|
---|
505 | if (!list.FindObject(obj->GetName()))
|
---|
506 | list.Add(obj);
|
---|
507 |
|
---|
508 | const MFilter *filter = task->GetFilter();
|
---|
509 |
|
---|
510 | if (!filter)
|
---|
511 | continue;
|
---|
512 |
|
---|
513 | TIter NextFBranch(filter->GetListOfBranches());
|
---|
514 | while ((obj=NextFBranch()))
|
---|
515 | if (!list.FindObject(obj->GetName()))
|
---|
516 | list.Add(obj);
|
---|
517 | }
|
---|
518 |
|
---|
519 | SetBranchStatus(&list, kTRUE);
|
---|
520 | /*
|
---|
521 | //
|
---|
522 | // Loop over all tasks iand its filters n the task list.
|
---|
523 | //
|
---|
524 | MTask *task;
|
---|
525 | TIter NextTask(tlist->GetList());
|
---|
526 | while ((task=(MTask*)NextTask()))
|
---|
527 | {
|
---|
528 | SetBranchStatus(task->GetListOfBranches(), kTRUE);
|
---|
529 |
|
---|
530 | const MFilter *filter = task->GetFilter();
|
---|
531 | if (!filter)
|
---|
532 | continue;
|
---|
533 |
|
---|
534 | SetBranchStatus(filter->GetListOfBranches(), kTRUE);
|
---|
535 |
|
---|
536 | }
|
---|
537 | */
|
---|
538 | }
|
---|
539 |
|
---|
540 | // --------------------------------------------------------------------------
|
---|
541 | //
|
---|
542 | // If the chain has been changed (by calling AddFile or using a file
|
---|
543 | // in the constructors argument) the number of entries is newly
|
---|
544 | // calculated from the files in the chain - this may take a while.
|
---|
545 | // The number of entries is returned.
|
---|
546 | //
|
---|
547 | UInt_t MReadTree::GetEntries()
|
---|
548 | {
|
---|
549 | if (TestBit(kChainWasChanged))
|
---|
550 | {
|
---|
551 | *fLog << inf << "Scanning chain... " << flush;
|
---|
552 | fNumEntries = (UInt_t)fChain->GetEntries();
|
---|
553 | *fLog << fNumEntries << " events found." << endl;
|
---|
554 | ResetBit(kChainWasChanged);
|
---|
555 | }
|
---|
556 | return fNumEntries;
|
---|
557 | }
|
---|
558 |
|
---|
559 | // --------------------------------------------------------------------------
|
---|
560 | //
|
---|
561 | // The disables all subbranches of the given master branch.
|
---|
562 | //
|
---|
563 | void MReadTree::DisableSubBranches(TBranch *branch)
|
---|
564 | {
|
---|
565 | //
|
---|
566 | // This is not necessary, it would work without. But the output
|
---|
567 | // may confuse the user...
|
---|
568 | //
|
---|
569 | if (fAutoEnable || fBranchChoosing)
|
---|
570 | return;
|
---|
571 |
|
---|
572 | SetBranchStatus(branch->GetListOfBranches(), kFALSE);
|
---|
573 | }
|
---|
574 |
|
---|
575 | // --------------------------------------------------------------------------
|
---|
576 | //
|
---|
577 | // The PreProcess loops (till now) over the branches in the given tree.
|
---|
578 | // It checks if the corresponding containers (containers with the same
|
---|
579 | // name than the branch name) are existing in the Parameter Container List.
|
---|
580 | // If not, a container of objec type 'branch-name' is created (everything
|
---|
581 | // after the last semicolon in the branch name is stripped). Only
|
---|
582 | // branches which don't have a veto (see VetoBranch) are enabled If the
|
---|
583 | // object isn't found in the root dictionary (a list of classes known by the
|
---|
584 | // root environment) the branch is skipped and an error message is printed
|
---|
585 | // out.
|
---|
586 | // If a selector is specified it is preprocessed after the
|
---|
587 | // MReadTree::PreProcess
|
---|
588 | //
|
---|
589 | Int_t MReadTree::PreProcess(MParList *pList)
|
---|
590 | {
|
---|
591 | //
|
---|
592 | // Make sure, that all the following calls doesn't result in
|
---|
593 | // Notifications. This may be dangerous, because the notified
|
---|
594 | // tasks are not preprocessed.
|
---|
595 | //
|
---|
596 | fChain->SetNotify(NULL);
|
---|
597 |
|
---|
598 | //
|
---|
599 | // check for files and for the tree!
|
---|
600 | //
|
---|
601 | if (!fChain->GetFile())
|
---|
602 | {
|
---|
603 | *fLog << err << GetDescriptor() << ": No file or no tree with name " << fChain->GetName() << " in file." << endl;
|
---|
604 | return kFALSE;
|
---|
605 | }
|
---|
606 |
|
---|
607 | //
|
---|
608 | // get number of events in this tree
|
---|
609 | //
|
---|
610 | if (!GetEntries())
|
---|
611 | {
|
---|
612 | *fLog << err << GetDescriptor() << ": No entries found in file(s)" << endl;
|
---|
613 | return kFALSE;
|
---|
614 | }
|
---|
615 |
|
---|
616 | //
|
---|
617 | // output logging information
|
---|
618 | //
|
---|
619 | *fLog << inf << fNumEntries << " entries found in file(s)." << endl;
|
---|
620 |
|
---|
621 | //
|
---|
622 | // Get all branches of this tree and
|
---|
623 | // create the Iterator to loop over all branches
|
---|
624 | //
|
---|
625 | TIter Next(fChain->GetListOfBranches());
|
---|
626 | TBranch *branch=NULL;
|
---|
627 |
|
---|
628 | Int_t num=0;
|
---|
629 |
|
---|
630 | //
|
---|
631 | // loop over all tasks for processing
|
---|
632 | //
|
---|
633 | while ( (branch=(TBranch*)Next()) )
|
---|
634 | {
|
---|
635 | //
|
---|
636 | // Get Name of Branch and Object
|
---|
637 | //
|
---|
638 | const char *bname = branch->GetName();
|
---|
639 |
|
---|
640 | TString oname(bname);
|
---|
641 | if (oname.EndsWith("."))
|
---|
642 | oname.Remove(oname.Length()-1);
|
---|
643 |
|
---|
644 | //
|
---|
645 | // Check if enabeling the branch is allowed
|
---|
646 | //
|
---|
647 | if (fVetoList->FindObject(oname))
|
---|
648 | {
|
---|
649 | *fLog << inf << "Master branch " << bname << " has veto... skipped." << endl;
|
---|
650 | DisableSubBranches(branch);
|
---|
651 | continue;
|
---|
652 | }
|
---|
653 |
|
---|
654 | //
|
---|
655 | // Create a pointer to the pointer to the object in which the
|
---|
656 | // branch data is stored. The pointers are stored in the TChain
|
---|
657 | // object and we get the pointers from there to delete it.
|
---|
658 | //
|
---|
659 | MParContainer **pcont= new MParContainer*;
|
---|
660 |
|
---|
661 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,02,06)
|
---|
662 | const char *classname = oname;
|
---|
663 | #else
|
---|
664 | const char *classname = branch->GetClassName();
|
---|
665 | #endif
|
---|
666 |
|
---|
667 | //
|
---|
668 | // check if object is existing in the list
|
---|
669 | //
|
---|
670 | *pcont=pList->FindCreateObj(classname, oname);
|
---|
671 |
|
---|
672 | if (!*pcont)
|
---|
673 | {
|
---|
674 | //
|
---|
675 | // if class is not existing in the (root) environment
|
---|
676 | // we cannot proceed reading this branch
|
---|
677 | //
|
---|
678 | *fLog << warn << dbginf << "Warning: Class '" << classname;
|
---|
679 | *fLog << "' for " << oname << " not existing in dictionary. Branch skipped." << endl;
|
---|
680 | DisableSubBranches(branch);
|
---|
681 | continue;
|
---|
682 | }
|
---|
683 |
|
---|
684 | //
|
---|
685 | // Check whether a Pointer to a pointer already exists.
|
---|
686 | // If we created one already, delete it.
|
---|
687 | //
|
---|
688 | TChainElement *element = (TChainElement*)fChain->GetStatus()->FindObject(bname);
|
---|
689 | if (element)
|
---|
690 | delete (MParContainer**)element->GetBaddress();
|
---|
691 |
|
---|
692 | //
|
---|
693 | // here pcont is a pointer the to container in which the data from
|
---|
694 | // the actual branch should be stored - enable branch.
|
---|
695 | //
|
---|
696 | fChain->SetBranchAddress(bname, pcont);
|
---|
697 |
|
---|
698 | *fLog << inf << "Master branch address " << bname << " [";
|
---|
699 | *fLog << classname << "] setup for reading." << endl;
|
---|
700 |
|
---|
701 | //*fLog << "Branch " << bname << " autodel: " << (int)branch->IsAutoDelete() << endl;
|
---|
702 | //branch->SetAutoDelete();
|
---|
703 |
|
---|
704 | num++;
|
---|
705 | }
|
---|
706 |
|
---|
707 | *fLog << inf << GetDescriptor() << " setup " << num << " master branches addresses." << endl;
|
---|
708 |
|
---|
709 | //
|
---|
710 | // If auto enabling scheme isn't disabled, do auto enabling
|
---|
711 | //
|
---|
712 | if (fAutoEnable)
|
---|
713 | EnableBranches(pList);
|
---|
714 |
|
---|
715 | //
|
---|
716 | // Now we can start notifying. Reset tree makes sure, that TChain thinks
|
---|
717 | // that the correct file is not yet initialized and reinitilizes it
|
---|
718 | // as soon as the first event is read. This is necessary to call
|
---|
719 | // the notifiers when the first event is read, but after the
|
---|
720 | // PreProcess-function.
|
---|
721 | //
|
---|
722 | fChain->ResetTree();
|
---|
723 | fChain->SetNotify(this);
|
---|
724 |
|
---|
725 | return GetSelector() ? GetSelector()->CallPreProcess(pList) : kTRUE;
|
---|
726 | }
|
---|
727 |
|
---|
728 | // --------------------------------------------------------------------------
|
---|
729 | //
|
---|
730 | // Set the ready to save flag of all containers which branchaddresses are
|
---|
731 | // set for. This is necessary to copy data.
|
---|
732 | //
|
---|
733 | void MReadTree::SetReadyToSave(Bool_t flag)
|
---|
734 | {
|
---|
735 | TIter Next(fChain->GetStatus());
|
---|
736 |
|
---|
737 | TChainElement *element = NULL;
|
---|
738 | while ((element=(TChainElement*)Next()))
|
---|
739 | {
|
---|
740 | //
|
---|
741 | // Check whether the branch is enabled
|
---|
742 | //
|
---|
743 | if (!element->GetStatus())
|
---|
744 | continue;
|
---|
745 |
|
---|
746 | //
|
---|
747 | // Get the pointer to the pointer of the corresponding container
|
---|
748 | //
|
---|
749 | MParContainer **pcont = (MParContainer**)element->GetBaddress();
|
---|
750 |
|
---|
751 | //
|
---|
752 | // Check whether the pointer is not NULL
|
---|
753 | //
|
---|
754 | if (!pcont || !*pcont)
|
---|
755 | continue;
|
---|
756 |
|
---|
757 | //
|
---|
758 | // Set the ready to save status of the container.
|
---|
759 | //
|
---|
760 | (*pcont)->SetReadyToSave(flag);
|
---|
761 | }
|
---|
762 |
|
---|
763 | //
|
---|
764 | // Set the ready to save status of this task (used?), too
|
---|
765 | //
|
---|
766 | MTask::SetReadyToSave(flag);
|
---|
767 | }
|
---|
768 |
|
---|
769 | // --------------------------------------------------------------------------
|
---|
770 | //
|
---|
771 | // The Process-function reads one event from the tree (this contains all
|
---|
772 | // enabled branches) and increases the position in the file by one event.
|
---|
773 | // (Remark: The position can also be set by some member functions
|
---|
774 | // If the end of the file is reached the Eventloop is stopped.
|
---|
775 | // In case an event selector is given its value is checked before
|
---|
776 | // reading the event. If it returns kAFLSE the event is skipped.
|
---|
777 | //
|
---|
778 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,02,06)
|
---|
779 | #include "MRawEvtData.h"
|
---|
780 | #endif
|
---|
781 | Int_t MReadTree::Process()
|
---|
782 | {
|
---|
783 | //
|
---|
784 | // This is necessary due to a bug in TChain::LoadTree in root.
|
---|
785 | // will be fixed in 3.03
|
---|
786 | //
|
---|
787 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,03,01)
|
---|
788 | if (fNumEntry >= GetEntries())
|
---|
789 | return kFALSE;
|
---|
790 | #endif
|
---|
791 |
|
---|
792 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,02,06)
|
---|
793 | //
|
---|
794 | // This fixes 99.9% of a memory leak using a root version prior
|
---|
795 | // to 3.02/??
|
---|
796 | //
|
---|
797 | TChainElement *element=NULL;
|
---|
798 | TIter Next(fChain->GetStatus());
|
---|
799 | while ((element=(TChainElement*)Next()))
|
---|
800 | {
|
---|
801 | MParContainer **c = (MParContainer**)element->GetBaddress();
|
---|
802 | if (!c) continue;
|
---|
803 | if ((*c)->InheritsFrom(MRawEvtData::Class()))
|
---|
804 | static_cast<MRawEvtData*>(*c)->DeletePixels(kFALSE);
|
---|
805 |
|
---|
806 | }
|
---|
807 | #endif
|
---|
808 |
|
---|
809 | if (GetSelector())
|
---|
810 | {
|
---|
811 | //
|
---|
812 | // Make sure selector is processed
|
---|
813 | //
|
---|
814 | if (!GetSelector()->CallProcess())
|
---|
815 | {
|
---|
816 | *fLog << err << dbginf << "Processing Selector failed." << endl;
|
---|
817 | return kFALSE;
|
---|
818 | }
|
---|
819 |
|
---|
820 | //
|
---|
821 | // Skip Event
|
---|
822 | //
|
---|
823 | if (!GetSelector()->IsConditionTrue())
|
---|
824 | {
|
---|
825 | fNumEntry++;
|
---|
826 | return kCONTINUE;
|
---|
827 | }
|
---|
828 | }
|
---|
829 |
|
---|
830 | const Bool_t rc = fChain->GetEntry(fNumEntry++) != 0;
|
---|
831 |
|
---|
832 | if (rc)
|
---|
833 | SetReadyToSave();
|
---|
834 |
|
---|
835 | return rc;
|
---|
836 | }
|
---|
837 |
|
---|
838 | // --------------------------------------------------------------------------
|
---|
839 | //
|
---|
840 | // If a selector is given the selector is post processed
|
---|
841 | //
|
---|
842 | Int_t MReadTree::PostProcess()
|
---|
843 | {
|
---|
844 | return GetSelector() ? GetSelector()->CallPostProcess() : kTRUE;
|
---|
845 | }
|
---|
846 |
|
---|
847 | // --------------------------------------------------------------------------
|
---|
848 | //
|
---|
849 | // Get the Event with the current EventNumber fNumEntry
|
---|
850 | //
|
---|
851 | Bool_t MReadTree::GetEvent()
|
---|
852 | {
|
---|
853 | Bool_t rc = fChain->GetEntry(fNumEntry) != 0;
|
---|
854 |
|
---|
855 | if (rc)
|
---|
856 | SetReadyToSave();
|
---|
857 |
|
---|
858 | return rc;
|
---|
859 | }
|
---|
860 |
|
---|
861 | // --------------------------------------------------------------------------
|
---|
862 | //
|
---|
863 | // Decrease the number of the event which is read by Process() next
|
---|
864 | // by one or more
|
---|
865 | //
|
---|
866 | Bool_t MReadTree::DecEventNum(UInt_t dec)
|
---|
867 | {
|
---|
868 | if (fNumEntry-dec >= GetEntries())
|
---|
869 | {
|
---|
870 | *fLog << warn << GetDescriptor() << ": DecEventNum, WARNING - Event " << fNumEntry << "-";
|
---|
871 | *fLog << dec << "=" << (Int_t)fNumEntry-dec << " out of Range." << endl;
|
---|
872 | return kFALSE;
|
---|
873 | }
|
---|
874 |
|
---|
875 | fNumEntry -= dec;
|
---|
876 | return kTRUE;
|
---|
877 | }
|
---|
878 |
|
---|
879 | // --------------------------------------------------------------------------
|
---|
880 | //
|
---|
881 | // Increase the number of the event which is read by Process() next
|
---|
882 | // by one or more
|
---|
883 | //
|
---|
884 | Bool_t MReadTree::IncEventNum(UInt_t inc)
|
---|
885 | {
|
---|
886 | if (fNumEntry+inc >= GetEntries())
|
---|
887 | {
|
---|
888 | *fLog << warn << GetDescriptor() << ": IncEventNum, WARNING - Event " << fNumEntry << "+";
|
---|
889 | *fLog << inc << "=" << (Int_t)fNumEntry+inc << " out of Range." << endl;
|
---|
890 | return kFALSE;
|
---|
891 | }
|
---|
892 |
|
---|
893 | fNumEntry += inc;
|
---|
894 | return kTRUE;
|
---|
895 | }
|
---|
896 |
|
---|
897 | // --------------------------------------------------------------------------
|
---|
898 | //
|
---|
899 | // This function makes Process() read event number nr next
|
---|
900 | //
|
---|
901 | // Remark: You can use this function after instatiating you MReadTree-object
|
---|
902 | // to set the event number from which you want to start reading.
|
---|
903 | //
|
---|
904 | Bool_t MReadTree::SetEventNum(UInt_t nr)
|
---|
905 | {
|
---|
906 | if (nr >= GetEntries())
|
---|
907 | {
|
---|
908 | *fLog << warn << GetDescriptor() << ": SetEventNum, WARNING - " << nr << " out of Range." << endl;
|
---|
909 | return kFALSE;
|
---|
910 | }
|
---|
911 |
|
---|
912 | fNumEntry = nr;
|
---|
913 | return kTRUE;
|
---|
914 | }
|
---|
915 |
|
---|
916 | // --------------------------------------------------------------------------
|
---|
917 | //
|
---|
918 | // For the branch with the given name:
|
---|
919 | // 1) no object is automatically created
|
---|
920 | // 2) the branch address for this branch is not set
|
---|
921 | // (because we lack the object, see 1)
|
---|
922 | // 3) The whole branch (exactly: all its subbranches) are disabled
|
---|
923 | // this means are not read in memory by TTree:GetEntry
|
---|
924 | //
|
---|
925 | void MReadTree::VetoBranch(const char *name)
|
---|
926 | {
|
---|
927 | fVetoList->Add(new TNamed(name, ""));
|
---|
928 | }
|
---|
929 |
|
---|
930 | // --------------------------------------------------------------------------
|
---|
931 | //
|
---|
932 | // Return the name of the file we are actually reading from.
|
---|
933 | //
|
---|
934 | TString MReadTree::GetFileName() const
|
---|
935 | {
|
---|
936 | const TFile *file = fChain->GetFile();
|
---|
937 |
|
---|
938 | if (!file)
|
---|
939 | return TString("<unknown>");
|
---|
940 |
|
---|
941 | TString name(file->GetName());
|
---|
942 | name.Remove(0, name.Last('/')+1);
|
---|
943 | return name;
|
---|
944 | }
|
---|
945 |
|
---|
946 | // --------------------------------------------------------------------------
|
---|
947 | //
|
---|
948 | // Return the number of the file in the chain, -1 in case of an error
|
---|
949 | //
|
---|
950 | Int_t MReadTree::GetFileIndex() const
|
---|
951 | {
|
---|
952 | return fChain->GetTreeNumber();
|
---|
953 | /*
|
---|
954 | const TString filename = fChain->GetFile()->GetName();
|
---|
955 |
|
---|
956 | int i=0;
|
---|
957 | TObject *file = NULL;
|
---|
958 |
|
---|
959 | TIter Next(fChain->GetListOfFiles());
|
---|
960 | while ((file=Next()))
|
---|
961 | {
|
---|
962 | if (filename==gSystem->ExpandPathName(file->GetTitle()))
|
---|
963 | return i;
|
---|
964 | i++;
|
---|
965 | }
|
---|
966 | return -1;
|
---|
967 | */
|
---|
968 | }
|
---|
969 |
|
---|
970 | // --------------------------------------------------------------------------
|
---|
971 | //
|
---|
972 | // This schedules a TObject which Notify(9 function is called in case
|
---|
973 | // of MReadTree (TChain) switches from one file in the chain to another
|
---|
974 | // one.
|
---|
975 | //
|
---|
976 | void MReadTree::AddNotify(TObject *obj)
|
---|
977 | {
|
---|
978 | fNotify->Add(obj);
|
---|
979 | }
|
---|
980 |
|
---|
981 | void MReadTree::Print(Option_t *o) const
|
---|
982 | {
|
---|
983 | *fLog << all << underline << GetDescriptor() << ":" << endl << dec;
|
---|
984 | *fLog << " Files [Tree]:" << endl;
|
---|
985 |
|
---|
986 | int i = 0;
|
---|
987 | TIter Next(fChain->GetListOfFiles());
|
---|
988 | TObject *obj = NULL;
|
---|
989 | while ((obj=Next()))
|
---|
990 | *fLog << " " << i++ << ") " << obj->GetTitle() << " [" << obj->GetName() << "]" << endl;
|
---|
991 |
|
---|
992 | *fLog << " Total Number of Entries: " << fNumEntries << endl;
|
---|
993 | *fLog << " Next Entry to read: " << fNumEntry << endl;
|
---|
994 | }
|
---|
995 |
|
---|
996 | // --------------------------------------------------------------------------
|
---|
997 | //
|
---|
998 | // Implementation of SavePrimitive. Used to write the call to a constructor
|
---|
999 | // to a macro. In the original root implementation it is used to write
|
---|
1000 | // gui elements to a macro-file.
|
---|
1001 | //
|
---|
1002 | void MReadTree::StreamPrimitive(ofstream &out) const
|
---|
1003 | {
|
---|
1004 | out << " " << ClassName() << " " << GetUniqueName() << "(\"";
|
---|
1005 | out << fChain->GetName() << "\", \"" << fName << "\", \"" << fTitle << "\");" << endl;
|
---|
1006 |
|
---|
1007 | TIter Next(fChain->GetListOfFiles());
|
---|
1008 | TObject *obj = NULL;
|
---|
1009 | while ((obj=Next()))
|
---|
1010 | out << " " << GetUniqueName() << ".AddFile(\"" << obj->GetTitle() << "\");" << endl;
|
---|
1011 |
|
---|
1012 | if (!fAutoEnable)
|
---|
1013 | out << " " << GetUniqueName() << ".DisableAutoScheme();" << endl;
|
---|
1014 |
|
---|
1015 | if (fNumEntry!=0)
|
---|
1016 | out << " " << GetUniqueName() << ".SetEventNum(" << fNumEntry << ");" << endl;
|
---|
1017 |
|
---|
1018 |
|
---|
1019 | }
|
---|