Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TBranchProxyDirector.cxx
Go to the documentation of this file.
1 // @(#)root/base:$Id$
2 // Author: Philippe Canal 13/05/2003
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun, Fons Rademakers and al. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** TBranchProxyDirector
13 This class is used to 'drive' and hold a serie of TBranchProxy objects
14 which represent and give access to the content of TTree object.
15 This is intended to be used as part of a generate Selector class
16 which will hold the directory and its associate
17 */
18 
19 #include "TBranchProxyDirector.h"
20 #include "TBranchProxy.h"
21 #include "TFriendProxy.h"
22 #include "TTree.h"
23 #include "TEnv.h"
24 #include "TH1F.h"
25 #include "TPad.h"
26 #include "TList.h"
27 
28 #include <algorithm>
29 
30 namespace std {} using namespace std;
31 
32 ClassImp(ROOT::Internal::TBranchProxyDirector);
33 
34 namespace ROOT {
35 namespace Internal {
36 
37  // Helper function to call Reset on each TBranchProxy
38  void NotifyDirected(Detail::TBranchProxy *x) { x->Notify(); }
39 
40  // Helper function to call SetReadEntry on all TFriendProxy
41  void ResetReadEntry(TFriendProxy *fp) { fp->ResetReadEntry(); }
42 
43  // Helper class to call Update on all TFriendProxy
44  struct Update {
45  Update(TTree *newtree) : fNewTree(newtree) {}
46  TTree *fNewTree;
47  void operator()(TFriendProxy *x) { x->Update(fNewTree); }
48  };
49 
50 
51  TBranchProxyDirector::TBranchProxyDirector(TTree* tree, Long64_t i) :
52  fTree(tree),
53  fEntry(i)
54  {
55  // Simple constructor
56  }
57 
58  TBranchProxyDirector::TBranchProxyDirector(TTree* tree, Int_t i) :
59  // cint has a problem casting int to long long
60  fTree(tree),
61  fEntry(i)
62  {
63  // Simple constructor
64  }
65 
66  void TBranchProxyDirector::Attach(Detail::TBranchProxy* p) {
67 
68  // Attach a TBranchProxy object to this director. The director just
69  // 'remembers' this BranchProxy and does not own it. It will be use
70  // to apply Tree wide operation (like reseting).
71  fDirected.push_back(p);
72  }
73 
74  void TBranchProxyDirector::Attach(TFriendProxy* p) {
75 
76  // Attach a TFriendProxy object to this director. The director just
77  // 'remembers' this BranchProxy and does not own it. It will be use
78  // to apply Tree wide operation (like reseting).
79  fFriends.push_back(p);
80  }
81 
82  TH1F* TBranchProxyDirector::CreateHistogram(const char *options) {
83  // Create a temporary 1D histogram.
84 
85  Int_t nbins = gEnv->GetValue("Hist.Binning.1D.x",100);
86  Double_t vmin=0, vmax=0;
87  Double_t xmin=0, xmax=0;
88  Bool_t canExtend = kTRUE;
89  TString opt( options );
90  Bool_t optSame = opt.Contains("same");
91  if (optSame) canExtend = kFALSE;
92 
93  if (gPad && optSame) {
94  TListIter np(gPad->GetListOfPrimitives());
95  TObject *op;
96  TH1 *oldhtemp = 0;
97  while ((op = np()) && !oldhtemp) {
98  if (op->InheritsFrom(TH1::Class())) oldhtemp = (TH1 *)op;
99  }
100  if (oldhtemp) {
101  nbins = oldhtemp->GetXaxis()->GetNbins();
102  vmin = oldhtemp->GetXaxis()->GetXmin();
103  vmax = oldhtemp->GetXaxis()->GetXmax();
104  } else {
105  vmin = gPad->GetUxmin();
106  vmax = gPad->GetUxmax();
107  }
108  } else {
109  vmin = xmin;
110  vmax = xmax;
111  if (xmin < xmax) canExtend = kFALSE;
112  }
113  TH1F *hist = new TH1F("htemp","htemp",nbins,vmin,vmax);
114  hist->SetLineColor(fTree->GetLineColor());
115  hist->SetLineWidth(fTree->GetLineWidth());
116  hist->SetLineStyle(fTree->GetLineStyle());
117  hist->SetFillColor(fTree->GetFillColor());
118  hist->SetFillStyle(fTree->GetFillStyle());
119  hist->SetMarkerStyle(fTree->GetMarkerStyle());
120  hist->SetMarkerColor(fTree->GetMarkerColor());
121  hist->SetMarkerSize(fTree->GetMarkerSize());
122  if (canExtend) hist->SetCanExtend(TH1::kAllAxes);
123  hist->GetXaxis()->SetTitle("var");
124  hist->SetBit(kCanDelete);
125  hist->SetDirectory(0);
126 
127  if (opt.Length() && opt.Contains("e")) hist->Sumw2();
128  return hist;
129  }
130 
131  TTree* TBranchProxyDirector::SetTree(TTree *newtree) {
132 
133  // Set the BranchProxy to be looking at a new tree.
134  // Reset all.
135  // Return the old tree.
136 
137  TTree* oldtree = fTree;
138  fTree = newtree;
139  if(!Notify()) return nullptr;
140  return oldtree;
141  }
142 
143  Bool_t TBranchProxyDirector::Notify() {
144  fEntry = -1;
145  bool retVal = true;
146  for_each(fDirected.begin(),fDirected.end(),NotifyDirected);
147  for (auto brProxy : fDirected) {
148  retVal = retVal && brProxy->Notify();
149  }
150  Update update(fTree);
151  for_each(fFriends.begin(),fFriends.end(),update);
152  return retVal;
153  }
154 
155 } // namespace Internal
156 } // namespace ROOT