00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #include "TBranchProxyDirector.h"
00024 #include "TBranchProxy.h"
00025 #include "TFriendProxy.h"
00026 #include "TTree.h"
00027 #include "TEnv.h"
00028 #include "TH1F.h"
00029 #include "TPad.h"
00030 #include "TList.h"
00031
00032 #include <algorithm>
00033
00034 namespace std {} using namespace std;
00035
00036 ClassImp(ROOT::TBranchProxyDirector);
00037
00038 namespace ROOT {
00039
00040
00041 void Reset(TBranchProxy *x) { x->Reset(); }
00042
00043
00044 void ResetReadEntry(TFriendProxy *x) { x->ResetReadEntry(); }
00045
00046
00047 struct Update {
00048 Update(TTree *newtree) : fNewTree(newtree) {}
00049 TTree *fNewTree;
00050 void operator()(TFriendProxy *x) { x->Update(fNewTree); }
00051 };
00052
00053
00054 TBranchProxyDirector::TBranchProxyDirector(TTree* tree, Long64_t i) :
00055 fTree(tree),
00056 fEntry(i)
00057 {
00058
00059 }
00060
00061 TBranchProxyDirector::TBranchProxyDirector(TTree* tree, Int_t i) :
00062
00063 fTree(tree),
00064 fEntry(i)
00065 {
00066
00067 }
00068
00069 void TBranchProxyDirector::Attach(TBranchProxy* p) {
00070
00071
00072
00073
00074 fDirected.push_back(p);
00075 }
00076
00077 void TBranchProxyDirector::Attach(TFriendProxy* p) {
00078
00079
00080
00081
00082 fFriends.push_back(p);
00083 }
00084
00085 TH1F* TBranchProxyDirector::CreateHistogram(const char *options) {
00086
00087
00088 Int_t nbins = gEnv->GetValue("Hist.Binning.1D.x",100);
00089 Double_t vmin=0, vmax=0;
00090 Double_t xmin=0, xmax=0;
00091 Bool_t canRebin = kTRUE;
00092 TString opt( options );
00093 Bool_t optSame = opt.Contains("same");
00094 if (optSame) canRebin = kFALSE;
00095
00096 if (gPad && optSame) {
00097 TListIter np(gPad->GetListOfPrimitives());
00098 TObject *op;
00099 TH1 *oldhtemp = 0;
00100 while ((op = np()) && !oldhtemp) {
00101 if (op->InheritsFrom(TH1::Class())) oldhtemp = (TH1 *)op;
00102 }
00103 if (oldhtemp) {
00104 nbins = oldhtemp->GetXaxis()->GetNbins();
00105 vmin = oldhtemp->GetXaxis()->GetXmin();
00106 vmax = oldhtemp->GetXaxis()->GetXmax();
00107 } else {
00108 vmin = gPad->GetUxmin();
00109 vmax = gPad->GetUxmax();
00110 }
00111 } else {
00112 vmin = xmin;
00113 vmax = xmax;
00114 if (xmin < xmax) canRebin = kFALSE;
00115 }
00116 TH1F *hist = new TH1F("htemp","htemp",nbins,vmin,vmax);
00117 hist->SetLineColor(fTree->GetLineColor());
00118 hist->SetLineWidth(fTree->GetLineWidth());
00119 hist->SetLineStyle(fTree->GetLineStyle());
00120 hist->SetFillColor(fTree->GetFillColor());
00121 hist->SetFillStyle(fTree->GetFillStyle());
00122 hist->SetMarkerStyle(fTree->GetMarkerStyle());
00123 hist->SetMarkerColor(fTree->GetMarkerColor());
00124 hist->SetMarkerSize(fTree->GetMarkerSize());
00125 if (canRebin) hist->SetBit(TH1::kCanRebin);
00126 hist->GetXaxis()->SetTitle("var");
00127 hist->SetBit(kCanDelete);
00128 hist->SetDirectory(0);
00129
00130 if (opt.Length() && opt.Contains("e")) hist->Sumw2();
00131 return hist;
00132 }
00133
00134 void TBranchProxyDirector::SetReadEntry(Long64_t entry) {
00135
00136
00137 fEntry = entry;
00138 if (!fFriends.empty()) {
00139 for_each(fFriends.begin(),fFriends.end(),ResetReadEntry);
00140 }
00141 }
00142
00143 TTree* TBranchProxyDirector::SetTree(TTree *newtree) {
00144
00145
00146
00147
00148
00149 TTree* oldtree = fTree;
00150 fTree = newtree;
00151 fEntry = -1;
00152
00153
00154 for_each(fDirected.begin(),fDirected.end(),Reset);
00155 Update update(fTree);
00156 for_each(fFriends.begin(),fFriends.end(),update);
00157 return oldtree;
00158 }
00159
00160 }