The Gaudi Framework  master (e98cfcff)
Loading...
Searching...
No Matches
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.
 
TBranch * getBranch (std::string_view section, std::string_view branch_name) override
 Access data branch by name: Get existing branch in read only mode.
 
int loadRefs (std::string_view section, std::string_view cnt, unsigned long entry, RootObjectRefs &refs) override
 Load references object from file.
 
void addParam (ParamMap &c, char *p)
 Helper function to read params table.
 
void addEntry (StringVec &c, char *val)
 Helper function to read string tables.
 
template<class C, class F>
StatusCode readBranch (TTree *t, const char *nam, C &v, F pmf)
 Helper function to read internal file tables.
 
bool get (const string &dsc, pair< string, ContainerSection > &e)
 Analyze the Sections table entries.
 
void analyzeMergeMap (StringVec &tmp)
 Build merge sections from the Sections table entries.
 
StatusCode readRefs () override
 Read reference tables.
 
string getEntry (const string &c)
 Helper function to convert string vectors to branch entries.
 
string getParam (const pair< string, string > &p)
 Helper function to convert parameter vectors to branch entries.
 
template<class C, class F>
StatusCode saveBranch (const char *nam, C &v, F pmf)
 Helper function to save internal tables.
 
StatusCode saveRefs () override
 Save/update reference tables.
 
- Public Member Functions inherited from Gaudi::RootDataConnection::Tool
TTree * refs () const
 
StringVecdbs () const
 
StringVecconts () const
 
StringVeclinks () const
 
ParamMapparams () const
 
MsgStreammsgSvc () const
 
IIncidentSvcincidentSvc () const
 
const std::string & name () const
 
Sectionssections () const
 
LinkSectionslinkSections () const
 
MergeSectionsmergeSections () const
 
virtual ~Tool ()=default
 Default destructor.
 

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.
 

Detailed Description

Description:

Concrete implementation to read objects from ROOT files.

Author
M.Frank
Version
1.0

Definition at line 28 of file RootTool.h.

Constructor & Destructor Documentation

◆ RootTool()

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

Standard constructor.

Definition at line 31 of file RootTool.h.

31{ 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 114 of file RootTool.h.

114{ c.push_back( val ); }

◆ addParam()

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

Helper function to read params table.

Definition at line 106 of file RootTool.h.

106 {
107 char* q = strchr( p, '=' );
108 if ( q ) {
109 *q = 0;
110 c.emplace_back( p, ++q );
111 }
112 }

◆ analyzeMergeMap()

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

Build merge sections from the Sections table entries.

Definition at line 161 of file RootTool.h.

161 {
164 pair<string, ContainerSection> e;
165 MsgStream& msg = msgSvc();
166 RootRef r;
167 int cnt = 0;
168 ls.clear();
169 ms.clear();
170 msg << MSG::VERBOSE;
171 r.dbase = r.container = r.link = r.clid = r.svc = r.entry = 0;
172 for ( const auto& i : tmp ) {
173 if ( get( i, e ) ) {
174 msg << "Added Merge Section:" << e.first << endmsg;
175 ms[e.first].push_back( e.second );
176 if ( e.first == "Links" )
177 r.link = e.second.start;
178 else if ( e.first == "Containers" )
179 r.container = e.second.start;
180 else if ( e.first == "Databases" )
181 r.dbase = e.second.start;
182 else if ( e.first == "Params" )
183 r.svc = e.second.start;
184 } else if ( i == "[END-OF-SECTION]" ) {
185 r.entry = cnt;
186 if ( msg.isActive() ) {
187 msg << "Link Section [" << r.entry << "," << ls.size() << "] -> D:" << r.dbase << " C:" << r.container
188 << " L:" << r.link << " P:" << r.svc << endmsg;
189 }
190 ls.push_back( r );
191 cnt++;
192 }
193 }
194 }
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition MsgStream.h:198
LinkSections & linkSections() const
MergeSections & mergeSections() const
RootDataConnection::LinkSections LinkSections
RootDataConnection::MergeSections MergeSections
bool get(const string &dsc, pair< string, ContainerSection > &e)
Analyze the Sections table entries.
Definition RootTool.h:142
constexpr double ms
@ VERBOSE
Definition IMessageSvc.h:22

◆ get()

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

Analyze the Sections table entries.

Definition at line 142 of file RootTool.h.

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

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

33 {
34 std::string n{ branch_name };
35 n += '.';
36 for ( int i = 0, m = n.length() - 1; i < m; ++i )
37 if ( !isalnum( n[i] ) ) n[i] = '_';
38 TTree* t = c->getSection( section );
39 TBranch* b = t ? t->GetBranch( n.c_str() ) : nullptr;
40 if ( !b ) b = t ? t->GetBranch( std::string{ branch_name }.c_str() ) : nullptr;
41 if ( b ) b->SetAutoDelete( kFALSE );
42 return b;
43 }

◆ getEntry()

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

Helper function to convert string vectors to branch entries.

Definition at line 213 of file RootTool.h.

213{ return c; }

◆ getParam()

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

Helper function to convert parameter vectors to branch entries.

Definition at line 215 of file RootTool.h.

215{ 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 45 of file RootTool.h.

45 {
46 TBranch* b = getBranch( section, std::string{ cnt } + "#R" );
47 RootObjectRefs* prefs = &refs;
48 if ( b ) {
49 b->SetAddress( &prefs );
50 int nb = b->GetEntry( entry );
51 if ( nb >= 1 ) {
52 const MergeSections& ms = c->mergeSections();
53 if ( !ms.empty() ) {
54 MsgStream& msg = msgSvc();
56 pair<const RootRef*, const ContainerSection*> ls = c->getMergeSection( cnt, entry );
57 if ( ls.first ) {
58 if ( ls.first->dbase >= 0 ) {
59 // Now patch the references and links 'en block' to be efficient
60 // First the leafs from the TES
61 if ( msg.isActive() ) {
62 msg << "Refs: LS [" << entry << "] -> " << ls.first->dbase << "," << ls.first->container << ","
63 << ls.first->link << "," << ls.first->entry << endmsg;
64 }
65 for ( size_t j = 0, n = refs.refs.size(); j < n; ++j ) {
66 RootRef& r = refs.refs[j];
67 if ( r.entry >= 0 && r.dbase >= 0 ) {
68 int db = r.dbase + ls.first->dbase;
69 if ( c->getDb( db ) == c->fid() ) {
70 if ( r.dbase ) r.dbase += ls.first->dbase;
71 if ( r.container ) r.container += ls.first->container;
72 if ( r.link ) r.link += ls.first->link;
73 const string& rc = c->getCont( r.container );
74 auto k = ms.find( rc );
75 if ( k != ms.end() ) {
76 const auto& cs = ( *k ).second;
77 r.entry = ( ls.first->entry >= 0 && ls.first->entry < (int)cs.size() )
78 ? cs[ls.first->entry].start + r.entry
79 : -1;
80 if ( msg.isActive() ) {
81 msg << "Add link [" << r.entry << "," << ls.first->entry << "," << ls.first->container << ","
82 << r.container << "," << r.entry << "] to -> " << rc << endmsg;
83 }
84 } else {
85 msg << MSG::WARNING << c->fid() << " [" << c->getDb( db ) << "] Evt:" << entry
86 << " Invalid link to " << rc << endmsg;
88 }
89 }
90 }
91 }
93 std::for_each( std::begin( refs.links ), std::end( refs.links ),
94 [&]( int& i ) { i += ls.first->link; } );
95 }
96 return nb;
97 }
98 return -1;
99 }
100 return nb;
101 }
102 }
103 return -1;
104 }
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:33
@ WARNING
Definition IMessageSvc.h:22

◆ readBranch()

template<class C, class F>
StatusCode Gaudi::RootTool::readBranch ( TTree * t,
const char * nam,
C & v,
F pmf )
inline

Helper function to read internal file tables.

Definition at line 117 of file RootTool.h.

117 {
118 char text[2048];
119 TBranch* b = t->GetBranch( nam );
120 if ( b ) {
121 TLeaf* l = b->GetLeaf( nam );
122 if ( l ) {
123 StatusCode sc = StatusCode::SUCCESS;
124 b->SetAddress( text );
125 msgSvc() << MSG::VERBOSE;
126 for ( Long64_t i = 0, n = b->GetEntries(); i < n; ++i ) {
127 if ( b->GetEntry( i ) > 0 ) {
128 char* p = (char*)l->GetValuePointer();
129 msgSvc() << "Add Value[" << b->GetName() << "]:" << p << endmsg;
130 ( this->*pmf )( v, p );
131 } else {
133 }
134 }
135 return sc;
136 }
137 }
138 msgSvc() << MSG::ERROR << "Failed to read '" << nam << "' table." << endmsg;
140 }
constexpr static const auto SUCCESS
Definition StatusCode.h:99
@ ERROR
Definition IMessageSvc.h:22
dict l
Definition gaudirun.py:583

◆ readRefs()

StatusCode Gaudi::RootTool::readRefs ( )
inlineoverridevirtual

Read reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 196 of file RootTool.h.

196 {
197 TTree* t = (TTree*)c->file()->Get( "Sections" );
198 StatusCode sc = StatusCode::SUCCESS;
199 StringVec tmp;
200 if ( t && !( sc = readBranch( t, "Sections", tmp, &RootTool::addEntry ) ).isSuccess() )
201 return sc;
202 else if ( refs() ) {
203 analyzeMergeMap( tmp );
204 if ( !( sc = readBranch( refs(), "Databases", dbs(), &RootTool::addEntry ) ).isSuccess() ) return sc;
205 if ( !( sc = readBranch( refs(), "Containers", conts(), &RootTool::addEntry ) ).isSuccess() ) return sc;
206 if ( !( sc = readBranch( refs(), "Links", links(), &RootTool::addEntry ) ).isSuccess() ) return sc;
207 if ( !( sc = readBranch( refs(), "Params", params(), &RootTool::addParam ) ).isSuccess() ) return sc;
208 return sc;
209 }
210 return StatusCode::FAILURE;
211 }
RootDataConnection::StringVec StringVec
StatusCode readBranch(TTree *t, const char *nam, C &v, F pmf)
Helper function to read internal file tables.
Definition RootTool.h:117
void addParam(ParamMap &c, char *p)
Helper function to read params table.
Definition RootTool.h:106
void analyzeMergeMap(StringVec &tmp)
Build merge sections from the Sections table entries.
Definition RootTool.h:161
void addEntry(StringVec &c, char *val)
Helper function to read string tables.
Definition RootTool.h:114
constexpr static const auto FAILURE
Definition StatusCode.h:100

◆ saveBranch()

template<class C, class F>
StatusCode Gaudi::RootTool::saveBranch ( const char * nam,
C & v,
F pmf )
inline

Helper function to save internal tables.

Definition at line 218 of file RootTool.h.

218 {
219 Long64_t i, n;
220 string val, typ = nam;
221 StatusCode sc = StatusCode::SUCCESS;
222 TDirectory::TContext ctxt( c->file() );
223 TBranch* b = refs()->GetBranch( nam );
224 if ( !b ) b = refs()->Branch( nam, 0, ( typ + "/C" ).c_str() );
225 if ( b ) {
226 for ( i = b->GetEntries(), n = Long64_t( v.size() ); i < n; ++i ) {
227 val = ( this->*pmf )( v[size_t( i )] );
228 b->SetAddress( (char*)val.c_str() );
229 msgSvc() << MSG::VERBOSE << "Save Value[" << b->GetName() << "]:" << val << endmsg;
230 if ( b->Fill() < 0 ) sc = StatusCode::FAILURE;
231 }
232 return sc;
233 }
234 return StatusCode::FAILURE;
235 }

◆ saveRefs()

StatusCode Gaudi::RootTool::saveRefs ( )
inlineoverridevirtual

Save/update reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 237 of file RootTool.h.

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

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