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

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