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

Last change on this file since 655 was 654, checked in by tbretz, 24 years ago
*** empty log message ***
File size: 6.1 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 <fstream.h>
18
19#include <TFile.h>
20#include <TTree.h>
21#include <TObjArray.h>
22
23#include "MLog.h"
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 and check if file is really open
46 //
47 fFile = new TFile(fFileName, "READ");
48
49 if (!fFile->IsOpen())
50 {
51 *fLog << "MReadTree::PreProcess: ERROR: Cannot open file '";
52 *fLog << fFileName << "'" << endl;
53 return kFALSE;
54 }
55
56 //
57 // try to get the tree and check if it was found
58 //
59 fTree = (TTree*)fFile->Get(fTreeName);
60 if (!fTree)
61 {
62 *fLog << "MReadTree::PreProcess: ERROR: Cannot open tree '";
63 *fLog << fTreeName << "'" << endl;
64 return kFALSE;
65 }
66
67 //
68 // get number of events in this tree
69 //
70 fNumEntries = (UInt_t)fTree->GetEntries();
71
72 //
73 // set pointer to first event
74 //
75 fNumEntry = 0;
76
77 //
78 // output logging information
79 //
80 *fLog << "File: '" << fFileName << "' Tree: '" << fTreeName;
81 *fLog << "' with " << fNumEntries << " Entries opened." << endl;
82
83 //
84 // Get all branches of this tree and
85 // create the Iterator to loop over all branches
86 //
87 TIter Next(fTree->GetListOfBranches());
88 TBranch *branch=NULL;
89
90 //
91 // loop over all tasks for processing
92 //
93 while ( (branch=(TBranch*)Next()) )
94 {
95 //
96 // Get Name of Branch
97 //
98 const char *name = branch->GetName();
99
100 //
101 // check if object is existing in the list
102 //
103 MParContainer *pcont = (MParContainer*)pList->FindObject(name);
104
105 if (!pcont)
106 {
107 //
108 // if object is not existing in the list try to create one
109 //
110 *fLog << "MReadTree::PreProcess - WARNING: '" << name << "' not found... creating." << endl;
111
112 //
113 // try to get class from root environment
114 //
115 TClass *cls = gROOT->GetClass(name);
116
117 if (!cls)
118 {
119 //
120 // if class is not existing in the root environment
121 // we cannot proceed reading this branch
122 //
123 *fLog << "MReadTree::PreProcess - Warning: Class '" << name << "' not existing in dictionary. Branch skipped." << endl;
124 continue;
125 }
126
127 //
128 // create the container and add it to the list
129 //
130 pcont = (MParContainer*)cls->New();
131 *fLog << pcont << endl;
132 pList->AddToList(pcont);
133 }
134
135 //
136 // here pcont is a pointer the to container in which the data from
137 // the actual branch should be stored - enable branch.
138 //
139 branch->SetAddress(&pcont);
140 }
141
142 return kTRUE;
143}
144
145Bool_t MReadTree::Process()
146{
147 //
148 // check for end of file
149 //
150 if (fNumEntry==fNumEntries)
151 return kFALSE;
152
153 //
154 // get entry
155 //
156 fTree->GetEntry(fNumEntry );
157
158 fNumEntry ++ ;
159
160 return kTRUE;
161}
162
163Bool_t MReadTree::PostProcess()
164{
165 //
166 // Close File
167 //
168 fFile->Close();
169
170 return kTRUE;
171}
172
173Bool_t MReadTree::GetEvent()
174{
175 // Get the Event with the current EventNumber fNumEntry
176 //
177 fTree->GetEntry(fNumEntry );
178
179 return kTRUE;
180}
181Bool_t MReadTree::DecEventNum(UInt_t dec)
182{
183 //
184 // Decrease the number of the event which is read by Process() next
185 // by one or more
186 //
187
188 //!
189 //! this function makes Process() read the event one (or more) before
190 //! the actual position (event) in the tree
191 //!
192 if (fNumEntry < dec/*+1*/)
193 {
194 *fLog << "MReadTree::SetPrevEvent: WARNING: " << fNumEntry/*-1*/ << "-" << dec << " out of Range." << endl;
195 return kFALSE;
196 }
197
198 fNumEntry -= dec/*+1*/;
199 return kTRUE;
200}
201
202Bool_t MReadTree::IncEventNum(UInt_t inc)
203{
204 //
205 // Increase the number of the event which is read by Process() next
206 // by one or more
207 //
208
209 //!
210 //! this function makes Process() read the next (or more) after
211 //! the actual position (event) in the tree
212 //! (Be careful: IncEventNum() or IncEventNum(1) does not chenge anything
213 //! in the standard behaviour of the task)
214 //!
215 if (fNumEntry+inc/*-1*/ >= fNumEntries)
216 {
217 *fLog << "MReadTree::SkipEvents: WARNING: " << fNumEntry/*-1*/ << "+" << inc << " out of Range." << endl;
218 return kFALSE;
219 }
220
221 fNumEntry += inc/*-1*/;
222 return kTRUE;
223}
224
225Bool_t MReadTree::SetEventNum(UInt_t nr)
226{
227 //
228 // this function makes Process() read event number nr next
229 //
230 if (nr>=fNumEntries)
231 {
232 *fLog << "MReadTree::SetEventNum: WARNING: " << nr << " out of Range." << endl;
233 return kFALSE;
234 }
235
236 fNumEntry = nr;
237 return kTRUE;
238}
239
Note: See TracBrowser for help on using the repository browser.