The Gaudi Framework  v33r0 (d5ea422b)
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 (std::string_view section, std::string_view branch_name) override
 Access data branch by name: Get existing branch in read only mode. More...
 
int loadRefs (std::string_view section, std::string_view 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 26 of file RootTool.h.

Constructor & Destructor Documentation

◆ RootTool()

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

Standard constructor.

Definition at line 29 of file RootTool.h.

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

Member Function Documentation

◆ addEntry()

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

Helper function to read string tables.

Definition at line 112 of file RootTool.h.

112 { c.push_back( val ); }
RootDataConnection * c
Pointer to containing data connection object.

◆ addParam()

void Gaudi::RootTool::addParam ( ParamMap c,
char *  p 
)
inline

Helper function to read params table.

Definition at line 104 of file RootTool.h.

104  {
105  char* q = strchr( p, '=' );
106  if ( q ) {
107  *q = 0;
108  c.emplace_back( p, ++q );
109  }
110  }
RootDataConnection * c
Pointer to containing data connection object.
T strchr(T... args)

◆ analyzeMergeMap()

void Gaudi::RootTool::analyzeMergeMap ( StringVec tmp)
inline

Build merge sections from the Sections table entries.

Definition at line 159 of file RootTool.h.

159  {
160  LinkSections& ls = linkSections();
162  pair<string, ContainerSection> e;
163  MsgStream& msg = msgSvc();
164  RootRef r;
165  int cnt = 0;
166  ls.clear();
167  ms.clear();
168  msg << MSG::VERBOSE;
169  r.dbase = r.container = r.link = r.clid = r.svc = r.entry = 0;
170  for ( const auto& i : tmp ) {
171  if ( get( i, e ) ) {
172  msg << "Added Merge Section:" << e.first << endmsg;
173  ms[e.first].push_back( e.second );
174  if ( e.first == "Links" )
175  r.link = e.second.start;
176  else if ( e.first == "Containers" )
177  r.container = e.second.start;
178  else if ( e.first == "Databases" )
179  r.dbase = e.second.start;
180  else if ( e.first == "Params" )
181  r.svc = e.second.start;
182  } else if ( i == "[END-OF-SECTION]" ) {
183  r.entry = cnt;
184  if ( msg.isActive() ) {
185  msg << "Link Section [" << r.entry << "," << ls.size() << "] -> D:" << r.dbase << " C:" << r.container
186  << " L:" << r.link << " P:" << r.svc << endmsg;
187  }
188  ls.push_back( r );
189  cnt++;
190  }
191  }
192  }
Definition of the MsgStream class used to transmit messages.
Definition: MsgStream.h:34
MergeSections & mergeSections() const
LinkSections & linkSections() const
bool get(const string &dsc, pair< string, ContainerSection > &e)
Analyze the Sections table entries.
Definition: RootTool.h:140
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:202

◆ get()

bool Gaudi::RootTool::get ( const string &  dsc,
pair< string, ContainerSection > &  e 
)
inline

Analyze the Sections table entries.

Definition at line 140 of file RootTool.h.

140  {
141  if ( dsc != "[END-OF-SECTION]" ) {
142  size_t id1 = dsc.find( "[CNT=" );
143  size_t id2 = dsc.find( "[START=" );
144  size_t id3 = dsc.find( "[LEN=" );
145  if ( id1 != string::npos && id2 != string::npos && id3 != string::npos ) {
146  string cnt = dsc.substr( id1 + 5, id2 - 1 - 5 );
147  int section_start = std::stoi( dsc.substr( id2 + 7, id3 - id2 - 8 ) );
148  int section_length = std::stoi( dsc.substr( id3 + 5, dsc.find( "]", id3 + 5 ) - id3 - 5 ) );
149  e.first = cnt;
150  e.second = ContainerSection( section_start, section_length );
151  return true;
152  }
153  }
154  e.first.clear();
155  e.second = ContainerSection( -1, -1 );
156  return false;
157  }
RootDataConnection::ContainerSection ContainerSection
T stoi(T... args)

◆ getBranch()

TBranch* Gaudi::RootTool::getBranch ( std::string_view  section,
std::string_view  branch_name 
)
inlineoverridevirtual

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

Implements Gaudi::RootDataConnection::Tool.

Definition at line 31 of file RootTool.h.

31  {
32  std::string n{branch_name};
33  n += '.';
34  for ( int i = 0, m = n.length() - 1; i < m; ++i )
35  if ( !isalnum( n[i] ) ) n[i] = '_';
36  TTree* t = c->getSection( section );
37  TBranch* b = t ? t->GetBranch( n.c_str() ) : nullptr;
38  if ( !b ) b = t ? t->GetBranch( std::string{branch_name}.c_str() ) : nullptr;
39  if ( b ) b->SetAutoDelete( kFALSE );
40  return b;
41  }
T isalnum(T... args)
STL class.
TTree * getSection(std::string_view sect, bool create=false)
Access TTree section from section name. The section is created if required.
constexpr double m
T c_str(T... args)
RootDataConnection * c
Pointer to containing data connection object.

◆ getEntry()

string Gaudi::RootTool::getEntry ( const string &  c)
inline

Helper function to convert string vectors to branch entries.

Definition at line 211 of file RootTool.h.

211 { return c; }
RootDataConnection * c
Pointer to containing data connection object.

◆ getParam()

string Gaudi::RootTool::getParam ( const pair< string, string > &  p)
inline

Helper function to convert parameter vectors to branch entries.

Definition at line 213 of file RootTool.h.

213 { return p.first + "=" + p.second; }

◆ loadRefs()

int Gaudi::RootTool::loadRefs ( std::string_view  section,
std::string_view  cnt,
unsigned long  entry,
RootObjectRefs refs 
)
inlineoverridevirtual

Load references object from file.

Link manager:

Implements Gaudi::RootDataConnection::Tool.

Definition at line 43 of file RootTool.h.

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

◆ readBranch()

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 115 of file RootTool.h.

115  {
116  char text[2048];
117  TBranch* b = t->GetBranch( nam );
118  if ( b ) {
119  TLeaf* l = b->GetLeaf( nam );
120  if ( l ) {
122  b->SetAddress( text );
123  msgSvc() << MSG::VERBOSE;
124  for ( Long64_t i = 0, n = b->GetEntries(); i < n; ++i ) {
125  if ( b->GetEntry( i ) > 0 ) {
126  char* p = (char*)l->GetValuePointer();
127  msgSvc() << "Add Value[" << b->GetName() << "]:" << p << endmsg;
128  ( this->*pmf )( v, p );
129  } else {
130  sc = RootDataConnection::Status::ROOT_READ_ERROR;
131  }
132  }
133  return sc;
134  }
135  }
136  msgSvc() << MSG::ERROR << "Failed to read '" << nam << "' table." << endmsg;
137  return RootDataConnection::Status::ROOT_READ_ERROR;
138  }
constexpr static const auto SUCCESS
Definition: StatusCode.h:96
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:61
dictionary l
Definition: gaudirun.py:543
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:202

◆ readRefs()

StatusCode Gaudi::RootTool::readRefs ( )
inlineoverridevirtual

Read reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 194 of file RootTool.h.

194  {
195  TTree* t = (TTree*)c->file()->Get( "Sections" );
196  StatusCode sc( StatusCode::SUCCESS, true );
197  StringVec tmp;
198  if ( t && !( sc = readBranch( t, "Sections", tmp, &RootTool::addEntry ) ).isSuccess() )
199  return sc;
200  else if ( refs() ) {
201  analyzeMergeMap( tmp );
202  if ( !( sc = readBranch( refs(), "Databases", dbs(), &RootTool::addEntry ) ).isSuccess() ) return sc;
203  if ( !( sc = readBranch( refs(), "Containers", conts(), &RootTool::addEntry ) ).isSuccess() ) return sc;
204  if ( !( sc = readBranch( refs(), "Links", links(), &RootTool::addEntry ) ).isSuccess() ) return sc;
205  if ( !( sc = readBranch( refs(), "Params", params(), &RootTool::addParam ) ).isSuccess() ) return sc;
206  return sc;
207  }
208  return StatusCode::FAILURE;
209  }
RootDataConnection::StringVec StringVec
constexpr static const auto SUCCESS
Definition: StatusCode.h:96
void analyzeMergeMap(StringVec &tmp)
Build merge sections from the Sections table entries.
Definition: RootTool.h:159
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:61
void addParam(ParamMap &c, char *p)
Helper function to read params table.
Definition: RootTool.h:104
void addEntry(StringVec &c, char *val)
Helper function to read string tables.
Definition: RootTool.h:112
StatusCode readBranch(TTree *t, const char *nam, C &v, F pmf)
Helper function to read internal file tables.
Definition: RootTool.h:115
constexpr static const auto FAILURE
Definition: StatusCode.h:97
RootDataConnection * c
Pointer to containing data connection object.
TFile * file() const
Direct access to TFile structure.

◆ saveBranch()

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 216 of file RootTool.h.

216  {
217  Long64_t i, n;
218  string val, typ = nam;
220  TDirectory::TContext ctxt( c->file() );
221  TBranch* b = refs()->GetBranch( nam );
222  if ( !b ) b = refs()->Branch( nam, 0, ( typ + "/C" ).c_str() );
223  if ( b ) {
224  for ( i = b->GetEntries(), n = Long64_t( v.size() ); i < n; ++i ) {
225  val = ( this->*pmf )( v[size_t( i )] );
226  b->SetAddress( (char*)val.c_str() );
227  msgSvc() << MSG::VERBOSE << "Save Value[" << b->GetName() << "]:" << val << endmsg;
228  if ( b->Fill() < 0 ) sc = StatusCode::FAILURE;
229  }
230  return sc;
231  }
232  return StatusCode::FAILURE;
233  }
constexpr static const auto SUCCESS
Definition: StatusCode.h:96
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:61
constexpr static const auto FAILURE
Definition: StatusCode.h:97
RootDataConnection * c
Pointer to containing data connection object.
TFile * file() const
Direct access to TFile structure.
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:202

◆ saveRefs()

StatusCode Gaudi::RootTool::saveRefs ( )
inlineoverridevirtual

Save/update reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 235 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  }
constexpr static const auto SUCCESS
Definition: StatusCode.h:96
string getEntry(const string &c)
Helper function to convert string vectors to branch entries.
Definition: RootTool.h:211
StatusCode saveBranch(const char *nam, C &v, F pmf)
Helper function to save internal tables.
Definition: RootTool.h:216
string getParam(const pair< string, string > &p)
Helper function to convert parameter vectors to branch entries.
Definition: RootTool.h:213
constexpr static const auto FAILURE
Definition: StatusCode.h:97

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