The Gaudi Framework  v31r0 (aeb156f0)
Gaudi::RootTool Class Reference

Description: More...

#include <src/RootTool.h>

Inheritance diagram for Gaudi::RootTool:
Collaboration diagram for Gaudi::RootTool:

Public Member Functions

 RootTool (RootDataConnection *con)
 Standard constructor. More...
 
TBranch * getBranch (boost::string_ref section, boost::string_ref branch_name) override
 Access data branch by name: Get existing branch in read only mode. More...
 
int loadRefs (boost::string_ref section, boost::string_ref cnt, unsigned long entry, RootObjectRefs &refs) override
 Load references object from file. More...
 
void addParam (ParamMap &c, char *p)
 Helper function to read params table. More...
 
void addEntry (StringVec &c, char *val)
 Helper function to read string tables. More...
 
template<class C , class F >
StatusCode readBranch (TTree *t, const char *nam, C &v, F pmf)
 Helper function to read internal file tables. More...
 
bool get (const string &dsc, pair< string, ContainerSection > &e)
 Analyze the Sections table entries. More...
 
void analyzeMergeMap (StringVec &tmp)
 Build merge sections from the Sections table entries. More...
 
StatusCode readRefs () override
 Read reference tables. More...
 
string getEntry (const string &c)
 Helper function to convert string vectors to branch entries. More...
 
string getParam (const pair< string, string > &p)
 Helper function to convert parameter vectors to branch entries. More...
 
template<class C , class F >
StatusCode saveBranch (const char *nam, C &v, F pmf)
 Helper function to save internal tables. More...
 
StatusCode saveRefs () override
 Save/update reference tables. More...
 
- Public Member Functions inherited from Gaudi::RootDataConnection::Tool
TTree * refs () const
 
StringVecdbs () const
 
StringVecconts () const
 
StringVeclinks () const
 
ParamMapparams () const
 
MsgStreammsgSvc () const
 
const std::stringname () const
 
Sectionssections () const
 
LinkSectionslinkSections () const
 
MergeSectionsmergeSections () const
 
virtual ~Tool ()=default
 Default destructor. More...
 
virtual RootRef poolRef (size_t) const
 Internal overload to facilitate the access to POOL files. More...
 

Additional Inherited Members

- Protected Types inherited from Gaudi::RootDataConnection::Tool
typedef RootDataConnection::StringVec StringVec
 
typedef RootDataConnection::ParamMap ParamMap
 
typedef RootDataConnection::Sections Sections
 
typedef RootDataConnection::MergeSections MergeSections
 
typedef RootDataConnection::LinkSections LinkSections
 
typedef RootDataConnection::ContainerSection ContainerSection
 
typedef RootDataConnection::ContainerSections ContainerSections
 
- Protected Attributes inherited from Gaudi::RootDataConnection::Tool
RootDataConnectionc
 Pointer to containing data connection object. More...
 

Detailed Description

Description:

Concrete implementation to read objects from POOL files.

Author
M.Frank
Version
1.0

Concrete implementation to read objects from ROOT files.

Author
M.Frank
Version
1.0

Definition at line 16 of file RootTool.h.

Constructor & Destructor Documentation

Gaudi::RootTool::RootTool ( RootDataConnection con)
inline

Standard constructor.

Definition at line 19 of file RootTool.h.

19 { c = con; }
RootDataConnection * c
Pointer to containing data connection object.

Member Function Documentation

void Gaudi::RootTool::addEntry ( StringVec c,
char *  val 
)
inline

Helper function to read string tables.

Definition at line 102 of file RootTool.h.

102 { c.push_back( val ); }
RootDataConnection * c
Pointer to containing data connection object.
void Gaudi::RootTool::addParam ( ParamMap c,
char *  p 
)
inline

Helper function to read params table.

Definition at line 94 of file RootTool.h.

94  {
95  char* q = strchr( p, '=' );
96  if ( q ) {
97  *q = 0;
98  c.emplace_back( p, ++q );
99  }
100  }
RootDataConnection * c
Pointer to containing data connection object.
T strchr(T...args)
void Gaudi::RootTool::analyzeMergeMap ( StringVec tmp)
inline

Build merge sections from the Sections table entries.

Definition at line 149 of file RootTool.h.

149  {
150  LinkSections& ls = linkSections();
152  pair<string, ContainerSection> e;
153  MsgStream& msg = msgSvc();
154  RootRef r;
155  int cnt = 0;
156  ls.clear();
157  ms.clear();
158  msg << MSG::VERBOSE;
159  r.dbase = r.container = r.link = r.clid = r.svc = r.entry = 0;
160  for ( const auto& i : tmp ) {
161  if ( get( i, e ) ) {
162  msg << "Added Merge Section:" << e.first << endmsg;
163  ms[e.first].push_back( e.second );
164  if ( e.first == "Links" )
165  r.link = e.second.start;
166  else if ( e.first == "Containers" )
167  r.container = e.second.start;
168  else if ( e.first == "Databases" )
169  r.dbase = e.second.start;
170  else if ( e.first == "Params" )
171  r.svc = e.second.start;
172  } else if ( i == "[END-OF-SECTION]" ) {
173  r.entry = cnt;
174  if ( msg.isActive() ) {
175  msg << "Link Section [" << r.entry << "," << ls.size() << "] -> D:" << r.dbase << " C:" << r.container
176  << " L:" << r.link << " P:" << r.svc << endmsg;
177  }
178  ls.push_back( r );
179  cnt++;
180  }
181  }
182  }
Definition of the MsgStream class used to transmit messages.
Definition: MsgStream.h:24
bool isActive() const
Accessor: is MsgStream active.
Definition: MsgStream.h:111
MergeSections & mergeSections() const
RootDataConnection::MergeSections MergeSections
RootDataConnection::LinkSections LinkSections
constexpr double ms
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:192
LinkSections & linkSections() const
bool Gaudi::RootTool::get ( const string &  dsc,
pair< string, ContainerSection > &  e 
)
inline

Analyze the Sections table entries.

Definition at line 130 of file RootTool.h.

130  {
131  if ( dsc != "[END-OF-SECTION]" ) {
132  size_t id1 = dsc.find( "[CNT=" );
133  size_t id2 = dsc.find( "[START=" );
134  size_t id3 = dsc.find( "[LEN=" );
135  if ( id1 != string::npos && id2 != string::npos && id3 != string::npos ) {
136  string cnt = dsc.substr( id1 + 5, id2 - 1 - 5 );
137  int section_start = std::stoi( dsc.substr( id2 + 7, id3 - id2 - 8 ) );
138  int section_length = std::stoi( dsc.substr( id3 + 5, dsc.find( "]", id3 + 5 ) - id3 - 5 ) );
139  e.first = cnt;
140  e.second = ContainerSection( section_start, section_length );
141  return true;
142  }
143  }
144  e.first.clear();
145  e.second = ContainerSection( -1, -1 );
146  return false;
147  }
RootDataConnection::ContainerSection ContainerSection
T stoi(T...args)
TBranch* Gaudi::RootTool::getBranch ( boost::string_ref  section,
boost::string_ref  branch_name 
)
inlineoverridevirtual

Access data branch by name: Get existing branch in read only mode.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 21 of file RootTool.h.

21  {
22  std::string n = branch_name.to_string() + ".";
23  for ( int i = 0, m = n.length() - 1; i < m; ++i )
24  if ( !isalnum( n[i] ) ) n[i] = '_';
25  TTree* t = c->getSection( section );
26  TBranch* b = t ? t->GetBranch( n.c_str() ) : nullptr;
27  if ( !b ) b = t ? t->GetBranch( branch_name.to_string().c_str() ) : nullptr;
28  if ( b ) b->SetAutoDelete( kFALSE );
29  return b;
30  }
T isalnum(T...args)
STL class.
constexpr double m
Definition: SystemOfUnits.h:92
TTree * getSection(boost::string_ref sect, bool create=false)
Access TTree section from section name. The section is created if required.
T length(T...args)
T c_str(T...args)
RootDataConnection * c
Pointer to containing data connection object.
string Gaudi::RootTool::getEntry ( const string &  c)
inline

Helper function to convert string vectors to branch entries.

Definition at line 201 of file RootTool.h.

201 { return c; }
RootDataConnection * c
Pointer to containing data connection object.
string Gaudi::RootTool::getParam ( const pair< string, string > &  p)
inline

Helper function to convert parameter vectors to branch entries.

Definition at line 203 of file RootTool.h.

203 { return p.first + "=" + p.second; }
int Gaudi::RootTool::loadRefs ( boost::string_ref  section,
boost::string_ref  cnt,
unsigned long  entry,
RootObjectRefs refs 
)
inlineoverridevirtual

Load references object from file.

Link manager:

Implements Gaudi::RootDataConnection::Tool.

Definition at line 32 of file RootTool.h.

33  {
34  TBranch* b = getBranch( section, cnt.to_string() + "#R" );
35  RootObjectRefs* prefs = &refs;
36  if ( b ) {
37  b->SetAddress( &prefs );
38  int nb = b->GetEntry( entry );
39  if ( nb >= 1 ) {
40  const MergeSections& ms = c->mergeSections();
41  if ( !ms.empty() ) {
42  MsgStream& msg = msgSvc();
43  msgSvc() << MSG::VERBOSE;
44  pair<const RootRef*, const ContainerSection*> ls = c->getMergeSection( cnt, entry );
45  if ( ls.first ) {
46  if ( ls.first->dbase >= 0 ) {
47  // Now patch the references and links 'en block' to be efficient
48  // First the leafs from the TES
49  if ( msg.isActive() ) {
50  msg << "Refs: LS [" << entry << "] -> " << ls.first->dbase << "," << ls.first->container << ","
51  << ls.first->link << "," << ls.first->entry << endmsg;
52  }
53  for ( size_t j = 0, n = refs.refs.size(); j < n; ++j ) {
54  RootRef& r = refs.refs[j];
55  if ( r.entry >= 0 && r.dbase >= 0 ) {
56  int db = r.dbase + ls.first->dbase;
57  if ( c->getDb( db ) == c->fid() ) {
58  if ( r.dbase ) r.dbase += ls.first->dbase;
59  if ( r.container ) r.container += ls.first->container;
60  if ( r.link ) r.link += ls.first->link;
61  const string& rc = c->getCont( r.container );
62  auto k = ms.find( rc );
63  if ( k != ms.end() ) {
64  const auto& cs = ( *k ).second;
65  r.entry = ( ls.first->entry >= 0 && ls.first->entry < (int)cs.size() )
66  ? cs[ls.first->entry].start + r.entry
67  : -1;
68  if ( msg.isActive() ) {
69  msg << "Add link [" << r.entry << "," << ls.first->entry << "," << ls.first->container << ","
70  << r.container << "," << r.entry << "] to -> " << rc << endmsg;
71  }
72  } else {
73  msg << MSG::WARNING << c->fid() << " [" << c->getDb( db ) << "] Evt:" << entry
74  << " Invalid link to " << rc << endmsg;
75  msg << MSG::VERBOSE;
76  }
77  }
78  }
79  }
81  std::for_each( std::begin( refs.links ), std::end( refs.links ),
82  [&]( int& i ) { i += ls.first->link; } );
83  }
84  return nb;
85  }
86  return -1;
87  }
88  return nb;
89  }
90  }
91  return -1;
92  }
Definition of the MsgStream class used to transmit messages.
Definition: MsgStream.h:24
const MergeSections & mergeSections() const
Access merged data section inventory.
const std::string & fid() const
Access file id.
const std::string & getDb(int which) const
Access database/file name from saved index.
bool isActive() const
Accessor: is MsgStream active.
Definition: MsgStream.h:111
T end(T...args)
TBranch * getBranch(boost::string_ref section, boost::string_ref branch_name) override
Access data branch by name: Get existing branch in read only mode.
Definition: RootTool.h:21
RootDataConnection::MergeSections MergeSections
const std::string & getCont(int which) const
Access container name from saved index.
T find(T...args)
T begin(T...args)
RootDataConnection * c
Pointer to containing data connection object.
std::pair< const RootRef *, const ContainerSection * > getMergeSection(boost::string_ref container, int entry) const
Access link section for single container and entry.
T for_each(T...args)
constexpr double ms
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:192
template<class C , class F >
StatusCode Gaudi::RootTool::readBranch ( TTree *  t,
const char *  nam,
C &  v,
pmf 
)
inline

Helper function to read internal file tables.

Definition at line 105 of file RootTool.h.

105  {
106  char text[2048];
107  TBranch* b = t->GetBranch( nam );
108  if ( b ) {
109  TLeaf* l = b->GetLeaf( nam );
110  if ( l ) {
112  b->SetAddress( text );
113  msgSvc() << MSG::VERBOSE;
114  for ( Long64_t i = 0, n = b->GetEntries(); i < n; ++i ) {
115  if ( b->GetEntry( i ) > 0 ) {
116  char* p = (char*)l->GetValuePointer();
117  msgSvc() << "Add Value[" << b->GetName() << "]:" << p << endmsg;
118  ( this->*pmf )( v, p );
119  } else {
120  sc = RootDataConnection::Status::ROOT_READ_ERROR;
121  }
122  }
123  return sc;
124  }
125  }
126  msgSvc() << MSG::ERROR << "Failed to read '" << nam << "' table." << endmsg;
127  return RootDataConnection::Status::ROOT_READ_ERROR;
128  }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
dictionary l
Definition: gaudirun.py:517
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:192
StatusCode Gaudi::RootTool::readRefs ( )
inlineoverridevirtual

Read reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 184 of file RootTool.h.

184  {
185  TTree* t = (TTree*)c->file()->Get( "Sections" );
186  StatusCode sc( StatusCode::SUCCESS, true );
187  StringVec tmp;
188  if ( t && !( sc = readBranch( t, "Sections", tmp, &RootTool::addEntry ) ).isSuccess() )
189  return sc;
190  else if ( refs() ) {
191  analyzeMergeMap( tmp );
192  if ( !( sc = readBranch( refs(), "Databases", dbs(), &RootTool::addEntry ) ).isSuccess() ) return sc;
193  if ( !( sc = readBranch( refs(), "Containers", conts(), &RootTool::addEntry ) ).isSuccess() ) return sc;
194  if ( !( sc = readBranch( refs(), "Links", links(), &RootTool::addEntry ) ).isSuccess() ) return sc;
195  if ( !( sc = readBranch( refs(), "Params", params(), &RootTool::addParam ) ).isSuccess() ) return sc;
196  return sc;
197  }
198  return StatusCode::FAILURE;
199  }
TFile * file() const
Direct access to TFile structure.
RootDataConnection::StringVec StringVec
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
void analyzeMergeMap(StringVec &tmp)
Build merge sections from the Sections table entries.
Definition: RootTool.h:149
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
void addParam(ParamMap &c, char *p)
Helper function to read params table.
Definition: RootTool.h:94
void addEntry(StringVec &c, char *val)
Helper function to read string tables.
Definition: RootTool.h:102
StatusCode readBranch(TTree *t, const char *nam, C &v, F pmf)
Helper function to read internal file tables.
Definition: RootTool.h:105
constexpr static const auto FAILURE
Definition: StatusCode.h:86
RootDataConnection * c
Pointer to containing data connection object.
template<class C , class F >
StatusCode Gaudi::RootTool::saveBranch ( const char *  nam,
C &  v,
pmf 
)
inline

Helper function to save internal tables.

Definition at line 206 of file RootTool.h.

206  {
207  Long64_t i, n;
208  string val, typ = nam;
210  TDirectory::TContext ctxt( c->file() );
211  TBranch* b = refs()->GetBranch( nam );
212  if ( !b ) b = refs()->Branch( nam, 0, ( typ + "/C" ).c_str() );
213  if ( b ) {
214  for ( i = b->GetEntries(), n = Long64_t( v.size() ); i < n; ++i ) {
215  val = ( this->*pmf )( v[size_t( i )] );
216  b->SetAddress( (char*)val.c_str() );
217  msgSvc() << MSG::VERBOSE << "Save Value[" << b->GetName() << "]:" << val << endmsg;
218  if ( b->Fill() < 0 ) sc = StatusCode::FAILURE;
219  }
220  return sc;
221  }
222  return StatusCode::FAILURE;
223  }
TFile * file() const
Direct access to TFile structure.
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
constexpr static const auto FAILURE
Definition: StatusCode.h:86
RootDataConnection * c
Pointer to containing data connection object.
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:192
StatusCode Gaudi::RootTool::saveRefs ( )
inlineoverridevirtual

Save/update reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 225 of file RootTool.h.

225  {
226  if ( refs() ) {
227  if ( !saveBranch( "Databases", dbs(), &RootTool::getEntry ).isSuccess() ) return StatusCode::FAILURE;
228  if ( !saveBranch( "Containers", conts(), &RootTool::getEntry ).isSuccess() ) return StatusCode::FAILURE;
229  if ( !saveBranch( "Links", links(), &RootTool::getEntry ).isSuccess() ) return StatusCode::FAILURE;
230  if ( !saveBranch( "Params", params(), &RootTool::getParam ).isSuccess() ) return StatusCode::FAILURE;
231  return StatusCode::SUCCESS;
232  }
233  return StatusCode::FAILURE;
234  }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
string getEntry(const string &c)
Helper function to convert string vectors to branch entries.
Definition: RootTool.h:201
StatusCode saveBranch(const char *nam, C &v, F pmf)
Helper function to save internal tables.
Definition: RootTool.h:206
string getParam(const pair< string, string > &p)
Helper function to convert parameter vectors to branch entries.
Definition: RootTool.h:203
constexpr static const auto FAILURE
Definition: StatusCode.h:86

The documentation for this class was generated from the following file: