Fork me on GitHub

source: git/external/ExRootAnalysis/ExRootTreeReader.cc@ b315535

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b315535 was cab38f6, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

remove svn tags and fix formatting

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