The Gaudi Framework  v29r0 (ff2e7097)
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 (CSTR section, CSTR branch_name) override
 Access data branch by name: Get existing branch in read only mode. More...
 
int loadRefs (CSTR section, CSTR 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 TBranch * getBranch (const std::string &section, const std::string &n)=0
 Access data branch by name: Get existing branch in read only mode. More...
 
virtual RootRef poolRef (size_t) const
 Internal overload to facilitate the access to POOL files. More...
 
virtual int loadRefs (const std::string &section, const std::string &cnt, unsigned long entry, RootObjectRefs &refs)=0
 Load references object. 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 17 of file RootTool.h.

Constructor & Destructor Documentation

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

Standard constructor.

Definition at line 21 of file RootTool.h.

21 { 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 106 of file RootTool.h.

106 { 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 97 of file RootTool.h.

98  {
99  char* q = strchr( p, '=' );
100  if ( q ) {
101  *q = 0;
102  c.emplace_back( p, ++q );
103  }
104  }
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 155 of file RootTool.h.

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

Analyze the Sections table entries.

Definition at line 135 of file RootTool.h.

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

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

Definition at line 23 of file RootTool.h.

24  {
25  std::string n = branch_name + ".";
26  for ( int i = 0, m = n.length() - 1; i < m; ++i )
27  if ( !isalnum( n[i] ) ) n[i] = '_';
28  TTree* t = c->getSection( section );
29  TBranch* b = t ? t->GetBranch( n.c_str() ) : 0;
30  if ( !b ) b = t ? t->GetBranch( branch_name.c_str() ) : 0;
31  if ( b ) b->SetAutoDelete( kFALSE );
32  return b;
33  }
T isalnum(T...args)
TTree * getSection(const std::string &sect, bool create=false)
Access TTree section from section name. The section is created if required.
STL class.
constexpr double m
Definition: SystemOfUnits.h:94
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 209 of file RootTool.h.

209 { 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 211 of file RootTool.h.

211 { return p.first + "=" + p.second; }
int Gaudi::RootTool::loadRefs ( CSTR  section,
CSTR  cnt,
unsigned long  entry,
RootObjectRefs refs 
)
inlineoverride

Load references object from file.

Link manager:

Definition at line 35 of file RootTool.h.

36  {
37  TBranch* b = getBranch( section, cnt + "#R" );
38  RootObjectRefs* prefs = &refs;
39  if ( b ) {
40  b->SetAddress( &prefs );
41  int nb = b->GetEntry( entry );
42  if ( nb >= 1 ) {
43  const MergeSections& ms = c->mergeSections();
44  if ( !ms.empty() ) {
45  MsgStream& msg = msgSvc();
46  msgSvc() << MSG::VERBOSE;
47  pair<const RootRef*, const ContainerSection*> ls = c->getMergeSection( cnt, entry );
48  if ( ls.first ) {
49  if ( ls.first->dbase >= 0 ) {
50  // Now patch the references and links 'en block' to be efficient
51  // First the leafs from the TES
52  if ( msg.isActive() ) {
53  msg << "Refs: LS [" << entry << "] -> " << ls.first->dbase << "," << ls.first->container << ","
54  << ls.first->link << "," << ls.first->entry << endmsg;
55  }
56  for ( size_t j = 0, n = refs.refs.size(); j < n; ++j ) {
57  RootRef& r = refs.refs[j];
58  if ( r.entry >= 0 && r.dbase >= 0 ) {
59  int db = r.dbase + ls.first->dbase;
60  if ( c->getDb( db ) == c->fid() ) {
61  if ( r.dbase ) r.dbase += ls.first->dbase;
62  if ( r.container ) r.container += ls.first->container;
63  if ( r.link ) r.link += ls.first->link;
64  const string& rc = c->getCont( r.container );
65  auto k = ms.find( rc );
66  if ( k != ms.end() ) {
67  const auto& cs = ( *k ).second;
68  r.entry = ( ls.first->entry >= 0 && ls.first->entry < (int)cs.size() )
69  ? cs[ls.first->entry].start + r.entry
70  : -1;
71  if ( msg.isActive() ) {
72  msg << "Add link [" << r.entry << "," << ls.first->entry << "," << ls.first->container << ","
73  << r.container << "," << r.entry << "] to -> " << rc << endmsg;
74  }
75  } else {
76  msg << MSG::WARNING << c->fid() << " [" << c->getDb( db ) << "] Evt:" << entry
77  << " Invalid link to " << rc << endmsg;
78  msg << MSG::VERBOSE;
79  }
80  }
81  }
82  }
84  std::for_each( std::begin( refs.links ), std::end( refs.links ),
85  [&]( int& i ) { i += ls.first->link; } );
86  }
87  return nb;
88  }
89  return -1;
90  }
91  return nb;
92  }
93  }
94  return -1;
95  }
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.
TBranch * getBranch(CSTR section, CSTR branch_name) override
Access data branch by name: Get existing branch in read only mode.
Definition: RootTool.h:23
const std::string & getDb(int which) const
Access database/file name from saved index.
bool isActive() const
Accessor: is MsgStream active.
Definition: MsgStream.h:116
T end(T...args)
std::pair< const RootRef *, const ContainerSection * > getMergeSection(const std::string &container, int entry) const
Access link section for single container and entry.
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.
T for_each(T...args)
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:209
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 109 of file RootTool.h.

110  {
111  char text[2048];
112  TBranch* b = t->GetBranch( nam );
113  if ( b ) {
114  TLeaf* l = b->GetLeaf( nam );
115  if ( l ) {
117  b->SetAddress( text );
118  msgSvc() << MSG::VERBOSE;
119  for ( Long64_t i = 0, n = b->GetEntries(); i < n; ++i ) {
120  if ( b->GetEntry( i ) > 0 ) {
121  char* p = (char*)l->GetValuePointer();
122  msgSvc() << "Add Value[" << b->GetName() << "]:" << p << endmsg;
123  ( this->*pmf )( v, p );
124  } else {
126  }
127  }
128  return sc;
129  }
130  }
131  msgSvc() << MSG::ERROR << "Failed to read '" << nam << "' table." << endmsg;
133  }
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
dictionary l
Definition: gaudirun.py:440
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:209
StatusCode Gaudi::RootTool::readRefs ( )
inlineoverridevirtual

Read reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 191 of file RootTool.h.

192  {
193  TTree* t = (TTree*)c->file()->Get( "Sections" );
194  StatusCode sc( StatusCode::SUCCESS, true );
195  StringVec tmp;
196  if ( t && !( sc = readBranch( t, "Sections", tmp, &RootTool::addEntry ) ).isSuccess() )
197  return sc;
198  else if ( refs() ) {
199  analyzeMergeMap( tmp );
200  if ( !( sc = readBranch( refs(), "Databases", dbs(), &RootTool::addEntry ) ).isSuccess() ) return sc;
201  if ( !( sc = readBranch( refs(), "Containers", conts(), &RootTool::addEntry ) ).isSuccess() ) return sc;
202  if ( !( sc = readBranch( refs(), "Links", links(), &RootTool::addEntry ) ).isSuccess() ) return sc;
203  if ( !( sc = readBranch( refs(), "Params", params(), &RootTool::addParam ) ).isSuccess() ) return sc;
204  return sc;
205  }
206  return StatusCode::FAILURE;
207  }
TFile * file() const
Direct access to TFile structure.
RootDataConnection::StringVec StringVec
void analyzeMergeMap(StringVec &tmp)
Build merge sections from the Sections table entries.
Definition: RootTool.h:155
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
void addParam(ParamMap &c, char *p)
Helper function to read params table.
Definition: RootTool.h:97
void addEntry(StringVec &c, char *val)
Helper function to read string tables.
Definition: RootTool.h:106
StatusCode readBranch(TTree *t, const char *nam, C &v, F pmf)
Helper function to read internal file tables.
Definition: RootTool.h:109
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 214 of file RootTool.h.

215  {
216  Long64_t i, n;
217  string val, typ = nam;
219  TDirectory::TContext ctxt( c->file() );
220  TBranch* b = refs()->GetBranch( nam );
221  if ( !b ) b = refs()->Branch( nam, 0, ( typ + "/C" ).c_str() );
222  if ( b ) {
223  for ( i = b->GetEntries(), n = Long64_t( v.size() ); i < n; ++i ) {
224  val = ( this->*pmf )( v[size_t( i )] );
225  b->SetAddress( (char*)val.c_str() );
226  msgSvc() << MSG::VERBOSE << "Save Value[" << b->GetName() << "]:" << val << endmsg;
227  if ( b->Fill() < 0 ) sc = StatusCode::FAILURE;
228  }
229  return sc;
230  }
231  return StatusCode::FAILURE;
232  }
TFile * file() const
Direct access to TFile structure.
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
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:209
StatusCode Gaudi::RootTool::saveRefs ( )
inlineoverridevirtual

Save/update reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 234 of file RootTool.h.

235  {
236  if ( refs() ) {
237  if ( !saveBranch( "Databases", dbs(), &RootTool::getEntry ).isSuccess() ) return StatusCode::FAILURE;
238  if ( !saveBranch( "Containers", conts(), &RootTool::getEntry ).isSuccess() ) return StatusCode::FAILURE;
239  if ( !saveBranch( "Links", links(), &RootTool::getEntry ).isSuccess() ) return StatusCode::FAILURE;
240  if ( !saveBranch( "Params", params(), &RootTool::getParam ).isSuccess() ) return StatusCode::FAILURE;
241  return StatusCode::SUCCESS;
242  }
243  return StatusCode::FAILURE;
244  }
string getEntry(const string &c)
Helper function to convert string vectors to branch entries.
Definition: RootTool.h:209
StatusCode saveBranch(const char *nam, C &v, F pmf)
Helper function to save internal tables.
Definition: RootTool.h:214
string getParam(const pair< string, string > &p)
Helper function to convert parameter vectors to branch entries.
Definition: RootTool.h:211

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