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

Last change on this file since 1115 was 1114, checked in by tbretz, 23 years ago
*** empty log message ***
File size: 21.0 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 <mailto: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// 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.h>
54
55#include <TFile.h> // TFile::GetName
56#include <TChain.h>
57#include <TGProgressBar.h>
58#include <TChainElement.h>
59#include <TOrdCollection.h>
60
61#include "MLog.h"
62#include "MLogManip.h"
63
64#include "MTime.h"
65#include "MFilter.h"
66#include "MParList.h"
67#include "MTaskList.h"
68
69ClassImp(MReadTree);
70
71class MChain : public TChain
72{
73public:
74 MChain() : TChain() {}
75 MChain(const char *name, const char *title="") : TChain(name, title) {}
76
77 void ResetTree() { fTree = 0; }
78};
79
80// --------------------------------------------------------------------------
81//
82// Default constructor. It creates an TChain instance which represents the
83// the Tree you want to read and adds the given file (if you gave one).
84// More files can be added using MReadTree::AddFile.
85// Also an empty veto list is created. This list is used if you want to
86// veto (disable or "don't enable") a branch in the tree, it vetos also
87// the creation of the corresponding object.
88// An empty list of TObjects are also created. This objects are called
89// at any time the TChain starts to read from another file.
90//
91MReadTree::MReadTree(const char *tname, const char *fname,
92 const char *name, const char *title)
93 : fNumEntry(0), fBranchChoosing(kFALSE), fAutoEnable(kTRUE), fProgress(NULL)
94{
95 fName = name ? name : "MReadTree";
96 fTitle = title ? title : "Task to loop over all events in one single tree";
97
98 fVetoList = new TList;
99 fVetoList->SetOwner();
100
101 fNotify = new TList;
102
103 //
104 // open the input stream
105 //
106 fChain = new MChain(tname);
107
108 // root 3.02:
109 // In TChain::Addfile remove the limitation that the file name must contain
110 // the string ".root". ".root" is necessary only in case one wants to specify
111 // a Tree in a subdirectory of a Root file with eg, the format:
112
113 if (fname)
114 fChain->Add(fname);
115}
116
117// --------------------------------------------------------------------------
118//
119// Destructor. It deletes the TChain and veto list object
120//
121MReadTree::~MReadTree()
122{
123 //
124 // Delete all the pointers to pointers to the objects where the
125 // branche data gets stored.
126 //
127 TIter Next(fChain->GetStatus());
128
129 TChainElement *element = NULL;
130 while ((element=(TChainElement*)Next()))
131 delete (MParContainer**)element->GetBaddress();
132
133 //
134 // Delete the chain and the veto list
135 //
136 delete fChain;
137 delete fNotify;
138 delete fVetoList;
139}
140
141// --------------------------------------------------------------------------
142//
143// If the owner flag is set all TObjects which are scheduled via
144// AddNotify are deleted by the destructor of MReadTree
145//
146void MReadTree::SetOwner(Bool_t flag)
147{
148 flag ? fNotify->SetBit(kIsOwner) : fNotify->ResetBit(kIsOwner);
149}
150
151// --------------------------------------------------------------------------
152//
153// This function is called each time MReadTree changes the file to read
154// from. It calls all TObject::Notify() functions which are scheduled
155// via AddNotify.
156//
157Bool_t MReadTree::Notify()
158{
159 *fLog << inf << "MReadTree: Notify '" << fChain->GetName() << "' ";
160 *fLog << "(before processing event #" << GetEventNum()-1 << ")" << endl;
161
162 //fNotify->Notify();
163
164 return kTRUE;
165}
166
167// --------------------------------------------------------------------------
168//
169// If you want to read the given tree over several files you must add
170// the files here before PreProcess is called. Be careful: If the tree
171// doesn't have the same contents (branches) it may confuse your
172// program (trees which are are not existing in later files are not read
173// anymore, tree wich are not existing in the first file are never read)
174//
175// Name may use the wildcarding notation, eg "xxx*.root" means all files
176// starting with xxx in the current file system directory.
177//
178// AddFile returns the number of files added to the chain.
179//
180Int_t MReadTree::AddFile(const char *fname)
181{
182 //
183 // FIXME! A check is missing whether the file already exists or not.
184 //
185 //
186 // returns the number of file which were added
187 //
188 return fChain->Add(fname);
189}
190
191// --------------------------------------------------------------------------
192//
193// This function is called if Branch choosing method should get enabled.
194// Branch choosing means, that only the enabled branches are read into
195// memory. To use an enableing scheme we have to disable all branches first.
196// This is done, if this function is called the first time.
197//
198void MReadTree::EnableBranchChoosing()
199{
200 if (fBranchChoosing)
201 return;
202
203 *fLog << inf << "Branch choosing method enabled (only enabled branches are read)." << endl;
204 fChain->SetBranchStatus("*", kFALSE);
205 fBranchChoosing = kTRUE;
206}
207
208// --------------------------------------------------------------------------
209//
210// The first time this function is called all branches are disabled.
211// The given branch is enabled. By enabling only the branches you
212// are processing you can speed up your calculation many times (up to
213// a factor of 10 or 20)
214//
215void MReadTree::EnableBranch(const char *name)
216{
217 EnableBranchChoosing();
218 fChain->SetBranchStatus(name, kTRUE);
219}
220
221// --------------------------------------------------------------------------
222//
223// Set branch status of branch name
224//
225void MReadTree::SetBranchStatus(const char *name, Bool_t status)
226{
227 fChain->SetBranchStatus(name, status);
228
229 *fLog << inf << (status ? "Enabled" : "Disabled");
230 *fLog << " subbranch '" << name << "'." << endl;
231}
232
233// --------------------------------------------------------------------------
234//
235// Checks whether a branch with the given name exists in the chain
236// and sets the branch status of this branch corresponding to status.
237//
238void MReadTree::SetBranchStatus(TObject *branch, Bool_t status)
239{
240 //
241 // Get branch name
242 //
243 const char *name = branch->GetName();
244
245 //
246 // Check whether this branch really exists
247 //
248 if (fChain->GetBranch(name))
249 SetBranchStatus(name, status);
250
251 //
252 // Remove trailing '.' if one and try to enable the subbranch without
253 // the master barnch name. This is to be compatible with older mars
254 // and camera files.
255 //
256 const char *dot = strrchr(name, '.');
257 if (!dot)
258 return;
259
260 if (fChain->GetBranch(dot+1))
261 SetBranchStatus(dot+1, kTRUE);
262}
263
264// --------------------------------------------------------------------------
265//
266// Set the status of all branches in the list to status.
267//
268void MReadTree::SetBranchStatus(const TList *list, Bool_t status)
269{
270 //
271 // Loop over all subbranches in this master branch
272 //
273 TIter Next(list);
274
275 TObject *obj;
276 while ((obj=Next()))
277 SetBranchStatus(obj, status);
278}
279
280// --------------------------------------------------------------------------
281//
282// This is the implementation of the Auto Enabling Scheme.
283// For more information see MTask::AddBranchToList.
284// This function loops over all tasks and its filters in the tasklist
285// and enables all branches which are requested by the tasks and its
286// filters.
287//
288// To enable 'unknown' branches which are not in the branchlist of
289// the tasks you can call EnableBranch
290//
291void MReadTree::EnableBranches(MParList *plist)
292{
293 //
294 // check whether branch choosing must be switched on
295 //
296 //EnableBranchChoosing();
297
298 //
299 // request the tasklist from the parameter list.
300 // FIXME: Tasklist can have a different name
301 //
302 const MTaskList *tlist = (MTaskList*)plist->FindObject("MTaskList");
303 if (!tlist)
304 {
305 *fLog << warn << "Cannot use auto enabeling scheme for branches. 'MTaskList' not found." << endl;
306 return;
307 }
308
309 //
310 // Loop over all tasks iand its filters n the task list.
311 //
312 MTask *task;
313 TIter NextTask(tlist->GetList());
314 while ((task=(MTask*)NextTask()))
315 {
316 SetBranchStatus(task->GetListOfBranches(), kTRUE);
317 const MFilter *filter = task->GetFilter();
318
319 if (filter)
320 SetBranchStatus(filter->GetListOfBranches(), kTRUE);
321
322 }
323}
324
325// --------------------------------------------------------------------------
326//
327// The disables all subbranches of the given master branch.
328//
329void MReadTree::DisableSubBranches(TBranch *branch)
330{
331 //
332 // This is not necessary, it would work without. But the output
333 // may confuse the user...
334 //
335 if (fAutoEnable || fBranchChoosing)
336 return;
337
338 SetBranchStatus(branch->GetListOfBranches(), kFALSE);
339}
340
341// --------------------------------------------------------------------------
342//
343// The PreProcess loops (till now) over the branches in the given tree.
344// It checks if the corresponding containers (containers with the same
345// name than the branch name) are existing in the Parameter Container List.
346// If not, a container of objec type 'branch-name' is created (everything
347// after the last semicolon in the branch name is stripped). Only
348// branches which don't have a veto (see VetoBranch) are enabled If the
349// object isn't found in the root dictionary (a list of classes known by the
350// root environment) the branch is skipped and an error message is printed
351// out.
352//
353Bool_t MReadTree::PreProcess(MParList *pList)
354{
355 //
356 // Make sure, that all the following calls doesn't result in
357 // Notifications. This may be dangerous, because the notified
358 // tasks are not preprocessed.
359 //
360 fChain->SetNotify(NULL);
361
362 //
363 // get number of events in this tree
364 //
365 fNumEntries = (UInt_t)fChain->GetEntries();
366
367 if (!fNumEntries)
368 {
369 *fLog << warn << dbginf << "No entries found in file(s)." << endl;
370 return kFALSE;
371 }
372
373 //
374 // output logging information
375 //
376 *fLog << inf << fNumEntries << " entries found in file(s)." << endl;
377
378 *fLog << all << fChain->GetListOfBranches() << endl;
379
380 //
381 // Get all branches of this tree and
382 // create the Iterator to loop over all branches
383 //
384 TIter Next(fChain->GetListOfBranches());
385 TBranch *branch=NULL;
386
387 Int_t num=0;
388
389
390 *fLog << "Start..." << endl;
391 //
392 // loop over all tasks for processing
393 //
394 while ( (branch=(TBranch*)Next()) )
395 {
396 //
397 // Get Name of Branch and Object
398 //
399 const char *bname = branch->GetName();
400
401 *fLog << all << bname << endl;
402
403 TString oname(bname);
404
405 if (oname.EndsWith("."))
406 oname.Remove(oname.Length()-1);
407
408 //
409 // Check if enabeling the branch is allowed
410 //
411 if (fVetoList->FindObject(oname))
412 {
413 *fLog << inf << "Master branch " << bname << " has veto... skipped." << endl;
414 DisableSubBranches(branch);
415 continue;
416 }
417
418 MParContainer *obj = pList->FindCreateObj(oname);
419
420 if (!obj)
421 {
422 //
423 // if class is not existing in the (root) environment
424 // we cannot proceed reading this branch
425 //
426 *fLog << warn << dbginf << "Warning: Class '" << oname << "' not existing in dictionary. Branch skipped." << endl;
427 DisableSubBranches(branch);
428 continue;
429 }
430
431 //
432 // Create a pointer to the pointer to the object in which the
433 // branch data is stored. The pointers are stored in the TChain
434 // object and we get the pointers from there to delete it.
435 //
436 MParContainer **pcont= new MParContainer*;
437
438 //
439 // check if object is existing in the list
440 //
441 *pcont=pList->FindCreateObj(oname);
442
443 if (!*pcont)
444 {
445 //
446 // if class is not existing in the (root) environment
447 // we cannot proceed reading this branch
448 //
449 *fLog << warn << dbginf << "Warning: Class '" << oname << "' not existing in dictionary. Branch skipped." << endl;
450 DisableSubBranches(branch);
451 continue;
452 }
453
454 //
455 // Check whether a Pointer to a pointer already exists, if
456 // we created one already delete it.
457 //
458 TChainElement *element = (TChainElement*)fChain->GetStatus()->FindObject(bname);
459 if (element)
460 delete (MParContainer**)element->GetBaddress();
461
462 //
463 // here pcont is a pointer the to container in which the data from
464 // the actual branch should be stored - enable branch.
465 //
466 fChain->SetBranchAddress(bname, pcont);
467 *fLog << inf << "Master branch address " << bname << " setup for reading." << endl;
468
469 //*fLog << "Branch " << bname << " autodel: " << (int)branch->IsAutoDelete() << endl;
470 //branch->SetAutoDelete();
471
472 num++;
473 }
474
475 *fLog << inf << "MReadTree setup " << num << " master branches addresses." << endl;
476
477 //
478 // If auto enabling scheme isn't disabled, do auto enabling
479 //
480 if (fAutoEnable)
481 EnableBranches(pList);
482
483 //
484 // If a progress bar is given set its range.
485 //
486 if (fProgress)
487 fProgress->SetRange(0, fNumEntries);
488
489 //
490 // Now we can start notifying. Reset tree makes sure, that TChain thinks
491 // that the correct file is not yet initialized and reinitilizes it
492 // as soon as the first event is read. This is necessary to call
493 // the notifiers when the first event is read, but after the
494 // PreProcess-function.
495 //
496 fChain->ResetTree();
497 fChain->SetNotify(this);
498
499 return kTRUE;
500}
501
502// --------------------------------------------------------------------------
503//
504// Set the ready to save flag of all containers which branchaddresses are
505// set for. This is necessary to copy data.
506//
507void MReadTree::SetReadyToSave(Bool_t flag)
508{
509 TIter Next(fChain->GetStatus());
510
511 TChainElement *element = NULL;
512 while ((element=(TChainElement*)Next()))
513 {
514 //
515 // Check whether the branch is enabled
516 //
517 if (!element->GetStatus())
518 continue;
519
520 //
521 // Get the pointer to the pointer of the corresponding container
522 //
523 MParContainer **pcont = (MParContainer**)element->GetBaddress();
524
525 //
526 // Check whether the pointer is not NULL
527 //
528 if (!pcont || !*pcont)
529 continue;
530
531 //
532 // Set the ready to save status of the container.
533 //
534 (*pcont)->SetReadyToSave(flag);
535 }
536
537 //
538 // Set the ready to save status of this task (used?), too
539 //
540 MTask::SetReadyToSave(flag);
541}
542
543// --------------------------------------------------------------------------
544//
545// The Process-function reads one event from the tree (this contains all
546// enabled branches) and increases the position in the file by one event.
547// (Remark: The position can also be set by some member functions
548// If the end of the file is reached the Eventloop is stopped.
549//
550Bool_t MReadTree::Process()
551{
552 //
553 // This is necessary due to a bug in TCHain::LoadTree in root.
554 //
555 if (fNumEntry >= fNumEntries)
556 return kFALSE;
557
558 Bool_t rc = fChain->GetEntry(fNumEntry++) != 0;
559
560 if (rc)
561 SetReadyToSave();
562
563 return rc;
564}
565
566// --------------------------------------------------------------------------
567//
568// Get the Event with the current EventNumber fNumEntry
569//
570Bool_t MReadTree::GetEvent()
571{
572 Bool_t rc = fChain->GetEntry(fNumEntry) != 0;
573
574 if (rc)
575 SetReadyToSave();
576
577 return rc;
578}
579
580// --------------------------------------------------------------------------
581//
582// Decrease the number of the event which is read by Process() next
583// by one or more
584//
585Bool_t MReadTree::DecEventNum(UInt_t dec)
586{
587 if (fNumEntry-dec >= fNumEntries)
588 {
589 *fLog << warn << "MReadTree::DecEventNum: WARNING - Event " << fNumEntry << "-";
590 *fLog << dec << "=" << (Int_t)fNumEntry-dec << " out of Range." << endl;
591 return kFALSE;
592 }
593
594 fNumEntry -= dec;
595 return kTRUE;
596}
597
598// --------------------------------------------------------------------------
599//
600// Increase the number of the event which is read by Process() next
601// by one or more
602//
603Bool_t MReadTree::IncEventNum(UInt_t inc)
604{
605 if (fNumEntry+inc >= fNumEntries)
606 {
607 *fLog << warn << "MReadTree::IncEventNum: WARNING - Event " << fNumEntry << "+";
608 *fLog << inc << "=" << (Int_t)fNumEntry+inc << " out of Range." << endl;
609 return kFALSE;
610 }
611
612 fNumEntry += inc;
613 return kTRUE;
614}
615
616// --------------------------------------------------------------------------
617//
618// This function makes Process() read event number nr next
619//
620// Remark: You can use this function after instatiating you MReadTree-object
621// to set the event number from which you want to start reading.
622//
623Bool_t MReadTree::SetEventNum(UInt_t nr)
624{
625 if (nr >= fNumEntries)
626 {
627 *fLog << warn << "MReadTree::SetEventNum: WARNING - " << nr << " out of Range." << endl;
628 return kFALSE;
629 }
630
631 fNumEntry = nr;
632 return kTRUE;
633}
634
635// --------------------------------------------------------------------------
636//
637// For the branch with the given name:
638// 1) no object is automatically created
639// 2) the branch address for this branch is not set
640// (because we lack the object, see 1)
641// 3) The whole branch (exactly: all its subbranches) are disabled
642// this means are not read in memory by TTree:GetEntry
643//
644void MReadTree::VetoBranch(const char *name)
645{
646 TNamed *branch = new TNamed(name, "");
647 fVetoList->Add(branch);
648}
649
650// --------------------------------------------------------------------------
651//
652// Return the name of the file we are actually reading from.
653//
654const char *MReadTree::GetFileName() const
655{
656 return fChain->GetFile()->GetName();
657}
658
659// --------------------------------------------------------------------------
660//
661// This schedules a TObject which Notify(9 function is called in case
662// of MReadTree (TChain) switches from one file in the chain to another
663// one.
664//
665void MReadTree::AddNotify(TObject *obj)
666{
667 fNotify->Add(obj);
668}
669
670void MReadTree::Print(Option_t *o) const
671{
672 *fLog << all << GetDescriptor() << dec << endl;
673 *fLog << setfill('-') << setw(strlen(GetDescriptor())) << "" << endl;
674 *fLog << " Files:" << endl;
675
676 int i = 0;
677 TIter Next(fChain->GetListOfFiles());
678 TObject *obj = NULL;
679 while ((obj=Next()))
680 *fLog << " " << i++ << ") " << obj->GetName() << endl;
681
682 *fLog << " Entries: " << fNumEntries << endl;
683 *fLog << " Next Entry: " << fNumEntry << endl;
684}
Note: See TracBrowser for help on using the repository browser.