Fork me on GitHub

source: svn/trunk/Utilities/ExRootAnalysis/src/ExRootTreeReader.cc@ 220

Last change on this file since 220 was 220, checked in by Xavier Rouby, 16 years ago

include statements have been cleaned

File size: 3.5 KB
RevLine 
[3]1
2/** \class ExRootTreeReader
3 *
4 * Class simplifying access to ROOT tree branches
5 *
[220]6 * $Date: 2009-02-02 11:32:01 $
7 * $Revision: 1.2 $
[3]8 *
9 *
10 * \author S. Ovyn - UCL, Louvain-la-Neuve
11 *
12 */
13
[220]14#include "ExRootTreeReader.h"
[3]15#include "TBranchElement.h"
16#include "TClonesArray.h"
[220]17#include "TChain.h"
[3]18#include <iostream>
19
20using namespace std;
21
22//------------------------------------------------------------------------------
23
24ExRootTreeReader::ExRootTreeReader(TTree *tree) :
25 fChain(tree), fCurrentTree(-1), fIsInitDone(kFALSE)
26{
27
28}
29
30//------------------------------------------------------------------------------
31
32ExRootTreeReader::~ExRootTreeReader()
33{
34 TBranchMap::iterator it_map;
35
36 for(it_map = fBranchMap.begin(); it_map != fBranchMap.end(); ++it_map)
37 {
38 delete it_map->second.second;
39 }
40}
41
42//------------------------------------------------------------------------------
43
44Bool_t ExRootTreeReader::ReadEntry(Long64_t entry)
45{
46 // Read contents of entry.
47 if(!fChain) return kFALSE;
48
49 if(!fIsInitDone) Init();
50
51 Int_t treeEntry = fChain->LoadTree(entry);
52 if(treeEntry < 0) return kFALSE;
53
54 if(fChain->IsA() == TChain::Class())
55 {
56 TChain *chain = static_cast<TChain*>(fChain);
57 if(chain->GetTreeNumber() != fCurrentTree)
58 {
59 fCurrentTree = chain->GetTreeNumber();
60 Notify();
61 }
62 }
63
64 TBranchMap::iterator it_map;
65 TBranch *branch;
66
67 for(it_map = fBranchMap.begin(); it_map != fBranchMap.end(); ++it_map)
68 {
69 branch = it_map->second.first;
70 if(branch)
71 {
72 branch->GetEntry(treeEntry);
73 }
74 }
75
76 return kTRUE;
77}
78
79//------------------------------------------------------------------------------
80
81TClonesArray *ExRootTreeReader::UseBranch(const TString &branchName)
82{
83 fIsInitDone = kFALSE;
84
85 TClonesArray *array = 0;
86
87 TBranchMap::iterator it_map = fBranchMap.find(branchName);
88
89 if(it_map != fBranchMap.end())
90 {
91 array = it_map->second.second;
92 }
93 else
94 {
95 TBranch *branch = fChain->GetBranch(branchName);
96 if(branch)
97 {
98 if(branch->IsA() == TBranchElement::Class())
99 {
100 TBranchElement *element = static_cast<TBranchElement*>(branch);
101 const char *className = element->GetClonesName();
102 Int_t size = element->GetMaximum();
103 TClass *cl = gROOT->GetClass(className);
104 if(cl)
105 {
106 array = new TClonesArray(cl, size);
107 fBranchMap.insert(make_pair(branchName, make_pair(branch, array)));
108 }
109 }
110 }
111 }
112
113 if(!array)
114 {
115 cout << "** WARNING: cannot access branch '" << branchName << "', return NULL pointer" << endl;
116 }
117
118 return array;
119}
120
121//------------------------------------------------------------------------------
122
123void ExRootTreeReader::Init()
124{
125 TBranchMap::iterator it_map;
126
127 for(it_map = fBranchMap.begin(); it_map != fBranchMap.end(); ++it_map)
128 {
129 fChain->SetBranchAddress(it_map->first, &(it_map->second.second));
130 }
131
132 Notify();
133
134 fIsInitDone = kTRUE;
135}
136
137//------------------------------------------------------------------------------
138
139void ExRootTreeReader::Notify()
140{
141 // Called when loading a new file.
142 // Get branch pointers.
143 if(!fChain) return;
144
145 TBranchMap::iterator it_map;
146
147 for(it_map = fBranchMap.begin(); it_map != fBranchMap.end(); ++it_map)
148 {
149 it_map->second.first = fChain->GetBranch(it_map->first);
150 if(!it_map->second.first)
151 {
152 cout << "** WARNING: cannot get branch '" << it_map->first << "'" << endl;
153 }
154 }
155}
156
157//------------------------------------------------------------------------------
158
159
Note: See TracBrowser for help on using the repository browser.