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

Last change on this file since 587 was 585, checked in by harald, 24 years ago
Adding some new code to start with the development of the usecase "Event Display".
File size: 5.8 KB
Line 
1/////////////////////////////////////////////////////////////////////////////
2// //
3// MReadTree //
4// //
5// This tasks opens all branches in a specified tree and creates the //
6// corresponding parameter containers if not already existing in the //
7// parameter list. //
8// //
9// The Process function reads one events from the tree. To go through the //
10// events of one tree make sure that the event number is increased from //
11// outside. It makes also possible to go back by decreasing the number. //
12// //
13/////////////////////////////////////////////////////////////////////////////
14
15#include "MReadTree.h"
16
17#include <iostream.h>
18#include <fstream.h>
19
20#include <TFile.h>
21#include <TTree.h>
22#include <TObjArray.h>
23
24#include "MTime.h"
25#include "MParList.h"
26
27ClassImp(MReadTree)
28
29MReadTree::MReadTree(const char *fname, const char *tname,
30 const char *name, const char *title)
31{
32 *fName = name ? name : "MReadTree";
33 *fTitle = title ? title : "Task to loop over all events in one single tree";
34
35 //
36 // open the input stream
37 //
38 fFileName = fname;
39 fTreeName = tname;
40}
41
42Bool_t MReadTree::PreProcess (MParList *pList)
43{
44 //
45 // open file
46 //
47 fFile = new TFile(fFileName, "READ");
48
49 if (!fFile->IsOpen())
50 {
51 cout << "MReadTree::PreProcess: ERROR: Cannot open file '";
52 cout << fFileName << "'" << endl;
53 return kFALSE;
54 }
55
56 fTree = (TTree*)fFile->Get(fTreeName);
57
58 if (!fTree)
59 {
60 cout << "MReadTree::PreProcess: ERROR: Cannot open tree '";
61 cout << fTreeName << "'" << endl;
62 return kFALSE;
63 }
64
65 fNumEntries = (UInt_t)fTree->GetEntries();
66 fNumEntry = 0;
67
68 cout << "File: '" << fFileName << "' Tree: '" << fTreeName;
69 cout << "' with " << fNumEntries << " Entries opened." << endl;
70
71 //
72 // Get all branches of this tree and
73 // create the Iterator to loop over all branches
74 //
75 TIter Next(fTree->GetListOfBranches());
76 TBranch *branch=NULL;
77
78 //
79 // loop over all tasks for processing
80 //
81 while ( (branch=(TBranch*)Next()) )
82 {
83 //
84 // Get Name of Branch
85 //
86 const char *name = branch->GetName();
87
88 //
89 // check if object is existing in the list
90 //
91 MParContainer *pcont = (MParContainer*)pList->FindObject(name);
92
93 if (!pcont)
94 {
95 //
96 // if object is not existing in the list try to create one
97 //
98 cout << "MReadTree::PreProcess - WARNING: '" << name << "' not found... creating." << endl;
99
100 //
101 // try to get class from root environment
102 //
103 TClass *cls = gROOT->GetClass(name);
104
105 if (!cls)
106 {
107 //
108 // if class is not existing in the root environment
109 // we cannot proceed reading this branch
110 //
111 cout << "MReadTree::PreProcess - Warning: Class '" << name << "' not existing in dictionary. Branch skipped." << endl;
112 continue;
113 }
114
115 //
116 // create the container and add it to the list
117 //
118 pcont = (MParContainer*)cls->New();
119 cout << pcont << endl;
120 pList->AddToList(pcont);
121 }
122
123 //
124 // here pcont is a pointer the to container in which the data from
125 // the actual branch should be stored - enable branch.
126 //
127 branch->SetAddress(&pcont);
128 }
129
130 return kTRUE;
131}
132
133Bool_t MReadTree::Process()
134{
135 //
136 // check for end of file
137 //
138 if (fNumEntry==fNumEntries)
139 return kFALSE;
140
141 //
142 // get entry
143 //
144 fTree->GetEntry(fNumEntry );
145
146 fNumEntry ++ ;
147
148 return kTRUE;
149}
150
151Bool_t MReadTree::PostProcess()
152{
153 //
154 // Close File
155 //
156 fFile->Close();
157
158 return kTRUE;
159}
160
161Bool_t MReadTree::GetEvent()
162{
163 // Get the Event with the current EventNumber fNumEntry
164 //
165 fTree->GetEntry(fNumEntry );
166
167 return kTRUE;
168}
169Bool_t MReadTree::DecEventNum(UInt_t dec)
170{
171 //
172 // Decrease the number of the event which is read by Process() next
173 // by one or more
174 //
175
176 //!
177 //! this function makes Process() read the event one (or more) before
178 //! the actual position (event) in the tree
179 //!
180 if (fNumEntry < dec/*+1*/)
181 {
182 cout << "MReadTree::SetPrevEvent: WARNING: " << fNumEntry/*-1*/ << "-" << dec << " out of Range." << endl;
183 return kFALSE;
184 }
185
186 fNumEntry -= dec/*+1*/;
187 return kTRUE;
188}
189
190Bool_t MReadTree::IncEventNum(UInt_t inc)
191{
192 //
193 // Increase the number of the event which is read by Process() next
194 // by one or more
195 //
196
197 //!
198 //! this function makes Process() read the next (or more) after
199 //! the actual position (event) in the tree
200 //! (Be careful: IncEventNum() or IncEventNum(1) does not chenge anything
201 //! in the standard behaviour of the task)
202 //!
203 if (fNumEntry+inc/*-1*/ >= fNumEntries)
204 {
205 cout << "MReadTree::SkipEvents: WARNING: " << fNumEntry/*-1*/ << "+" << inc << " out of Range." << endl;
206 return kFALSE;
207 }
208
209 fNumEntry += inc/*-1*/;
210 return kTRUE;
211}
212
213Bool_t MReadTree::SetEventNum(UInt_t nr)
214{
215 //
216 // this function makes Process() read event number nr next
217 //
218 if (nr>=fNumEntries)
219 {
220 cout << "MReadTree::SetEventNum: WARNING: " << nr << " out of Range." << endl;
221 return kFALSE;
222 }
223
224 fNumEntry = nr;
225 return kTRUE;
226}
227
Note: See TracBrowser for help on using the repository browser.