Fork me on GitHub

source: git/external/ExRootAnalysis/ExRootTreeReader.cc@ 5862d1f

Last change on this file since 5862d1f was dc883b4, checked in by Pavel Demin <pavel.demin@…>, 4 years ago

add info to TreeWriter

  • Property mode set to 100644
File size: 3.8 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
[341014c]12#include "TBranchElement.h"
[d7d2da3]13#include "TCanvas.h"
[dc883b4]14#include "TParameter.h"
[d7d2da3]15#include "TClonesArray.h"
[341014c]16#include "TH2.h"
17#include "TStyle.h"
[d7d2da3]18
19#include <iostream>
20
21using namespace std;
22
23//------------------------------------------------------------------------------
24
25ExRootTreeReader::ExRootTreeReader(TTree *tree) :
26 fChain(tree), fCurrentTree(-1)
27{
28}
29
30//------------------------------------------------------------------------------
31
32ExRootTreeReader::~ExRootTreeReader()
33{
34 TBranchMap::iterator itBranchMap;
35
36 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
37 {
38 delete itBranchMap->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 Int_t treeEntry = fChain->LoadTree(entry);
50 if(treeEntry < 0) return kFALSE;
51
52 if(fChain->IsA() == TChain::Class())
53 {
[341014c]54 TChain *chain = static_cast<TChain *>(fChain);
[d7d2da3]55 if(chain->GetTreeNumber() != fCurrentTree)
56 {
57 fCurrentTree = chain->GetTreeNumber();
58 Notify();
59 }
60 }
61
62 TBranchMap::iterator itBranchMap;
63 TBranch *branch;
64
65 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
66 {
67 branch = itBranchMap->second.first;
68 if(branch)
69 {
70 branch->GetEntry(treeEntry);
71 }
72 }
73
74 return kTRUE;
75}
76
77//------------------------------------------------------------------------------
78
79TClonesArray *ExRootTreeReader::UseBranch(const char *branchName)
80{
81 TClonesArray *array = 0;
82
83 TBranchMap::iterator itBranchMap = fBranchMap.find(branchName);
84
85 if(itBranchMap != fBranchMap.end())
86 {
87 cout << "** WARNING: branch '" << branchName << "' is already in use" << endl;
88 array = itBranchMap->second.second;
89 }
90 else
91 {
92 TBranch *branch = fChain->GetBranch(branchName);
93 if(branch)
94 {
95 if(branch->IsA() == TBranchElement::Class())
96 {
[341014c]97 TBranchElement *element = static_cast<TBranchElement *>(branch);
[d7d2da3]98 const char *className = element->GetClonesName();
99 Int_t size = element->GetMaximum();
100 TClass *cl = gROOT->GetClass(className);
101 if(cl)
102 {
103 array = new TClonesArray(cl, size);
104 array->SetName(branchName);
105 fBranchMap.insert(make_pair(branchName, make_pair(branch, array)));
106 branch->SetAddress(&array);
107 }
108 }
109 }
110 }
111
112 if(!array)
113 {
114 cout << "** WARNING: cannot access branch '" << branchName << "', return NULL pointer" << endl;
115 }
116
117 return array;
118}
119
120//------------------------------------------------------------------------------
121
[dc883b4]122Double_t ExRootTreeReader::GetInfo(const char *name)
123{
124 TTree *tree = NULL;
125 TParameter<Double_t> *param = NULL;
126 Double_t result = -999.9;
127 if(fChain) tree = fChain->GetTree();
128 if(tree) param = static_cast<TParameter<Double_t> *>(tree->GetUserInfo()->FindObject(name));
129 if(param) result = param->GetVal();
130 return result;
131}
132
133//------------------------------------------------------------------------------
134
[d7d2da3]135Bool_t ExRootTreeReader::Notify()
136{
137 // Called when loading a new file.
138 // Get branch pointers.
139 if(!fChain) return kFALSE;
140
141 TBranchMap::iterator itBranchMap;
142 TBranch *branch;
143
144 for(itBranchMap = fBranchMap.begin(); itBranchMap != fBranchMap.end(); ++itBranchMap)
145 {
146 branch = fChain->GetBranch(itBranchMap->first);
147 if(branch)
148 {
149 itBranchMap->second.first = branch;
150 branch->SetAddress(&(itBranchMap->second.second));
151 }
152 else
153 {
154 cout << "** WARNING: cannot get branch '" << itBranchMap->first << "'" << endl;
155 }
156 }
157 return kTRUE;
158}
159
160//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.