Logo ROOT  
Reference Guide
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
13This class is used to 'drive' and hold a serie of TBranchProxy objects
14which represent and give access to the content of TTree object.
15This is intended to be used as part of a generate Selector class
16which will hold the directory and its associate
17*/
18
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
30namespace std {} using namespace std;
31
33
34namespace ROOT {
35namespace 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
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
52 fTree(tree),
53 fEntry(i)
54 {
55 // Simple constructor
56 }
57
59 // cint has a problem casting int to long long
60 fTree(tree),
61 fEntry(i)
62 {
63 // Simple constructor
64 }
65
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
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
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);
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
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
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
void Class()
Definition: Class.C:29
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
long long Long64_t
Definition: RtypesCore.h:69
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
R__EXTERN TEnv * gEnv
Definition: TEnv.h:171
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
@ kCanDelete
Definition: TObject.h:339
TRObject operator()(const T1 &t1) const
#define gPad
Definition: TVirtualPad.h:286
Base class for all the proxy object.
Definition: TBranchProxy.h:69
void Attach(Detail::TBranchProxy *p)
std::list< Detail::TBranchProxy * > fDirected
std::vector< TFriendProxy * > fFriends
TBranchProxyDirector(const TBranchProxyDirector &)
TH1F * CreateHistogram(const char *options)
void ResetReadEntry()
Refresh the cached read entry number from the original tree.
virtual Color_t GetFillColor() const
Return the fill area color.
Definition: TAttFill.h:30
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition: TAttFill.h:31
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition: TAttFill.h:39
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition: TAttLine.h:42
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:34
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition: TAttMarker.h:32
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition: TAttMarker.h:31
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition: TAttMarker.h:33
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition: TAttMarker.h:41
Double_t GetXmax() const
Definition: TAxis.h:134
Double_t GetXmin() const
Definition: TAxis.h:133
Int_t GetNbins() const
Definition: TAxis.h:121
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition: TEnv.cxx:491
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
The TH1 histogram class.
Definition: TH1.h:56
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
Definition: TH1.cxx:8381
@ kAllAxes
Definition: TH1.h:73
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:316
virtual UInt_t SetCanExtend(UInt_t extendBitMask)
Make the histogram axes extendable / not extendable according to the bit mask returns the previous bi...
Definition: TH1.cxx:6279
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition: TH1.cxx:8464
Iterator of linked list.
Definition: TList.h:200
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
Mother of all ROOT objects.
Definition: TObject.h:37
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
Basic string class.
Definition: TString.h:131
Ssiz_t Length() const
Definition: TString.h:405
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
A TTree represents a columnar dataset.
Definition: TTree.h:72
Double_t x[n]
Definition: legend1.C:17
void ResetReadEntry(TFriendProxy *fp)
void NotifyDirected(Detail::TBranchProxy *x)
VSD Structures.
Definition: StringConv.hxx:21
Definition: tree.py:1