source: HiSusy/trunk/Delphes-3.0.0/external/ExRootAnalysis/ExRootTreeReader.cc @ 1

Last change on this file since 1 was 1, checked in by zerwas, 11 years ago

first import of structure, PYTHIA8 and DELPHES

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