Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TBranchProxyDirector.h
Go to the documentation of this file.
1// @(#)root/base:$Id$
2// Author: Philippe Canal 13/05/2003
3
4/*************************************************************************
5 * Copyright (C) 1995-2004, 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#ifndef ROOT_TBranchProxyDirector
13#define ROOT_TBranchProxyDirector
14
15#include "RtypesCore.h"
16#include <vector>
17#include <list>
18#include <algorithm>
19
20class TH1F;
21class TTree;
22
23namespace ROOT {
24namespace Detail {
25 class TBranchProxy;
26 class TFriendProxy;
27}
28
29namespace Internal{
30 class TFriendProxy;
31
32 /// Helper function to call SetReadEntry on all TFriendProxy
34
36
37 //This class could actually be the selector itself.
38 TTree *fTree; ///< TTree we are currently looking at.
39 Long64_t fEntry; ///< Entry currently being read (in the local TTree rather than the TChain)
40
41 std::list<Detail::TBranchProxy*> fDirected;
42 std::vector<TFriendProxy*> fFriends;
43
46
47 public:
48
50 TBranchProxyDirector(TTree* tree, Int_t i); // cint has (had?) a problem casting int to long long
51
53 void Attach(TFriendProxy* f);
54 TH1F* CreateHistogram(const char *options);
55
56 /// Return the current 'local' entry number; i.e. in the 'local' TTree rather than the TChain.
57 /// This value will be passed directly to TBranch::GetEntry.
58 Long64_t GetReadEntry() const { return fEntry; }
59
60 TTree* GetTree() const { return fTree; };
61 // void Print();
62
63 /// Move to a new entry to read
64 /// entry is the 'local' entry number; i.e. in the 'local' TTree rather than the TChain.
65 /// This value will be passed directly to TBranch::GetEntry.
66 void SetReadEntry(Long64_t entry) {
67 fEntry = entry;
68 if (!fFriends.empty()) {
69 std::for_each(fFriends.begin(), fFriends.end(), ResetReadEntry);
70 }
71 }
72 TTree* SetTree(TTree *newtree);
73 bool Notify();
74
75 };
76
77} // namespace Internal
78} // namespace ROOT
79
80#endif
#define f(i)
Definition RSha256.hxx:104
int Int_t
Definition RtypesCore.h:45
long long Long64_t
Definition RtypesCore.h:80
winID h TVirtualViewer3D TVirtualGLPainter p
Base class for all the proxy object.
Long64_t GetReadEntry() const
Return the current 'local' entry number; i.e.
Long64_t fEntry
Entry currently being read (in the local TTree rather than the TChain)
TTree * fTree
TTree we are currently looking at.
std::list< Detail::TBranchProxy * > fDirected
void SetReadEntry(Long64_t entry)
Move to a new entry to read entry is the 'local' entry number; i.e.
std::vector< TFriendProxy * > fFriends
TTree * SetTree(TTree *newtree)
Set the BranchProxy to be looking at a new tree.
TBranchProxyDirector(const TBranchProxyDirector &)
TH1F * CreateHistogram(const char *options)
Create a temporary 1D histogram.
void Attach(Detail::TBranchProxy *p)
Attach a TBranchProxy object to this director.
TBranchProxyDirector & operator=(const TBranchProxyDirector &)
Concrete implementation of the proxy around a Friend Tree.
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
A TTree represents a columnar dataset.
Definition TTree.h:79
void ResetReadEntry(TFriendProxy *fp)
Helper function to call SetReadEntry on all TFriendProxy.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...