Fork me on GitHub

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

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

first commit

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