Gaudi Framework, version v23r6

Home   Generated: Wed Jan 30 2013
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Public Member Functions | List of all members
Gaudi::RootTool Class Reference

Description: More...

#include <src/RootTool.h>

Inheritance diagram for Gaudi::RootTool:
Inheritance graph
[legend]
Collaboration diagram for Gaudi::RootTool:
Collaboration graph
[legend]

Public Member Functions

 RootTool (RootDataConnection *con)
 Standard constructor.
 
virtual TBranch * getBranch (CSTR section, CSTR branch_name)
 Access data branch by name: Get existing branch in read only mode.
 
virtual int loadRefs (CSTR section, CSTR cnt, unsigned long entry, RootObjectRefs &refs)
 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 ()
 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 ()
 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
 
const std::stringname () const
 
Sectionssections () const
 
LinkSectionslinkSections () const
 
MergeSectionsmergeSections () const
 
virtual ~Tool ()
 Default destructor.
 
virtual void release ()
 Use releasePtr() helper to delete object.
 
virtual TBranch * getBranch (const std::string &section, const std::string &n)=0
 Access data branch by name: Get existing branch in read only mode.
 
virtual RootRef poolRef (size_t) const
 Internal overload to facilitate the access to POOL files.
 
virtual int loadRefs (const std::string &section, const std::string &cnt, unsigned long entry, RootObjectRefs &refs)=0
 Load references object.
 

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

Constructor & Destructor Documentation

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

Standard constructor.

Definition at line 18 of file RootTool.h.

{ c = con; }

Member Function Documentation

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

Helper function to read string tables.

Definition at line 107 of file RootTool.h.

{ c.push_back(val); }
void Gaudi::RootTool::addParam ( ParamMap c,
char *  p 
)
inline

Helper function to read params table.

Definition at line 99 of file RootTool.h.

{
char* q = strchr(p,'=');
if ( q ) {
*q = 0;
c.push_back(make_pair(p,++q));
}
}
void Gaudi::RootTool::analyzeMergeMap ( StringVec tmp)
inline

Build merge sections from the Sections table entries.

Definition at line 155 of file RootTool.h.

{
StringVec::const_iterator i;
pair<string,ContainerSection> e;
RootRef r;
int cnt = 0;
ls.clear();
ms.clear();
msg << MSG::VERBOSE;
r.dbase = r.container = r.link = r.clid = r.svc = r.entry = 0;
for(i=tmp.begin(); i!=tmp.end();++i) {
if ( get(*i,e) ) {
msg << "Added Merge Section:" << e.first << endmsg;
ms[e.first].push_back(e.second);
if ( e.first == "Links" )
r.link = e.second.start;
else if ( e.first == "Containers" )
r.container = e.second.start;
else if ( e.first == "Databases" )
r.dbase = e.second.start;
else if ( e.first == "Params" )
r.svc = e.second.start;
}
else if ( (*i) == "[END-OF-SECTION]" ) {
r.entry = cnt;
if ( msg.isActive() ) {
msg << "Link Section [" << r.entry << "," << ls.size()
<< "] -> D:" << r.dbase
<< " C:" << r.container
<< " L:" << r.link
<< " P:" << r.svc
<< endmsg;
}
ls.push_back(r);
cnt++;
}
}
}
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.

{
if ( dsc != "[END-OF-SECTION]" ) {
size_t id1 = dsc.find("[CNT=");
size_t id2 = dsc.find("[START=");
size_t id3 = dsc.find("[LEN=");
if ( id1 != string::npos && id2 != string::npos && id3 != string::npos ) {
string tmp;
string cnt = dsc.substr(id1+5, id2-1-5);
int section_start = ::atoi((tmp=dsc.substr(id2+7,id3-id2-8)).c_str());
int section_length = ::atoi((tmp=dsc.substr(id3+5,dsc.find("]",id3+5)-id3-5)).c_str());
e.first = cnt;
e.second = ContainerSection(section_start,section_length);
return true;
}
}
e.first = "";
e.second = ContainerSection(-1,-1);
return false;
}
virtual TBranch* Gaudi::RootTool::getBranch ( CSTR  section,
CSTR  branch_name 
)
inlinevirtual

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

Definition at line 20 of file RootTool.h.

{
std::string n = branch_name+".";
for(int i=0, m=n.length()-1; i<m; ++i) if ( !isalnum(n[i]) ) n[i]='_';
TTree* t = c->getSection(section);
TBranch* b = t ? t->GetBranch(n.c_str()) : 0;
if ( !b ) b = t ? t->GetBranch(branch_name.c_str()) : 0;
if ( b ) b->SetAutoDelete(kFALSE);
return b;
}
string Gaudi::RootTool::getEntry ( const string c)
inline

Helper function to convert string vectors to branch entries.

Definition at line 217 of file RootTool.h.

{ return c; }
string Gaudi::RootTool::getParam ( const pair< string, string > &  p)
inline

Helper function to convert parameter vectors to branch entries.

Definition at line 219 of file RootTool.h.

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

Load references object from file.

Link manager:

Definition at line 30 of file RootTool.h.

{
TBranch* b = getBranch(section,cnt+"#R");
RootObjectRefs* prefs = &refs;
if ( b ) {
b->SetAddress(&prefs);
int nb = b->GetEntry(entry);
if ( nb >= 1 ) {
const MergeSections& ms = c->mergeSections();
if ( !ms.empty() ) {
pair<const RootRef*,const ContainerSection*> ls = c->getMergeSection(cnt,entry);
if ( ls.first ) {
if ( ls.first->dbase >= 0 ) {
// Now patch the references and links 'en block' to be efficient
// First the leafs from the TES
if ( msg.isActive() ) {
msg << "Refs: LS [" << entry << "] -> "
<< ls.first->dbase << "," << ls.first->container
<< "," << ls.first->link
<< "," << ls.first->entry
<< endmsg;
}
for(size_t j=0, n=refs.refs.size(); j<n; ++j) {
RootRef& r = refs.refs[j];
if ( r.entry>= 0 && r.dbase >= 0 ) {
int db = r.dbase + ls.first->dbase;
if ( c->getDb(db) == c->fid() ) {
if ( r.dbase ) r.dbase += ls.first->dbase;
if ( r.container ) r.container += ls.first->container;
if ( r.link ) r.link += ls.first->link;
const string& rc = c->getCont(r.container);
MergeSections::const_iterator k = ms.find(rc);
if ( k != ms.end() ) {
const ContainerSections& cs = (*k).second;
r.entry = ( ls.first->entry >= 0 && ls.first->entry < (int)cs.size() )
? cs[ls.first->entry].start + r.entry : -1;
if ( msg.isActive() ) {
msg << "Add link [" << r.entry
<< "," << ls.first->entry
<< "," << ls.first->container
<< "," << r.container
<< "," << r.entry
<< "] to -> " << rc << endmsg;
}
}
else {
msg << MSG::WARNING << c->fid() << " [" << c->getDb(db)
<< "] Evt:" << entry << " Invalid link to " << rc << endmsg;
msg << MSG::VERBOSE;
}
}
}
}
for(vector<int>::iterator i=refs.links.begin(); i!=refs.links.end();++i) {
(*i) += ls.first->link;
}
}
return nb;
}
return -1;
}
return nb;
}
}
return -1;
}
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.

{
char text[2048];
TBranch* b = t->GetBranch(nam);
if ( b ) {
TLeaf* l = b->GetLeaf(nam);
if ( l ) {
b->SetAddress(text);
for(Long64_t i=0, n=b->GetEntries(); i<n; ++i) {
if ( b->GetEntry(i)>0 ) {
char* p = (char*)l->GetValuePointer();
msgSvc() << "Add Value[" << b->GetName() << "]:" << p << endmsg;
(this->*pmf)(v,p);
}
else {
sc = RootDataConnection::ROOT_READ_ERROR;
}
}
return sc;
}
}
msgSvc() << MSG::ERROR << "Failed to read '" << nam << "' table." << endmsg;
return RootDataConnection::ROOT_READ_ERROR;
}
StatusCode Gaudi::RootTool::readRefs ( )
inlinevirtual

Read reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 196 of file RootTool.h.

{
TTree* t = (TTree*)c->file()->Get("Sections");
StringVec tmp;
if ( t && !(sc=readBranch(t, "Sections", tmp, &RootTool::addEntry)).isSuccess() )
return sc;
else if ( refs() ) {
if ( !(sc=readBranch(refs(),"Databases", dbs(), &RootTool::addEntry)).isSuccess() )
return sc;
if ( !(sc=readBranch(refs(),"Containers",conts(), &RootTool::addEntry)).isSuccess() )
return sc;
if ( !(sc=readBranch(refs(),"Links", links(), &RootTool::addEntry)).isSuccess() )
return sc;
if ( !(sc=readBranch(refs(),"Params", params(), &RootTool::addParam)).isSuccess() )
return sc;
return sc;
}
}
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 221 of file RootTool.h.

{
Long64_t i, n;
string val, typ = nam;
TDirectory::TContext ctxt(c->file());
TBranch* b = refs()->GetBranch(nam);
if ( !b ) b = refs()->Branch(nam,0,(typ+"/C").c_str());
if ( b ) {
for(i=b->GetEntries(), n=Long64_t(v.size()); i<n; ++i) {
val = (this->*pmf)(v[size_t(i)]);
b->SetAddress((char*)val.c_str());
msgSvc() << MSG::VERBOSE << "Save Value[" << b->GetName() << "]:" << val << endmsg;
if ( b->Fill() < 0 ) sc = StatusCode::FAILURE;
}
return sc;
}
}
StatusCode Gaudi::RootTool::saveRefs ( )
inlinevirtual

Save/update reference tables.

Implements Gaudi::RootDataConnection::Tool.

Definition at line 240 of file RootTool.h.

{
if ( refs() ) {
if ( !saveBranch("Databases", dbs(), &RootTool::getEntry).isSuccess() )
if ( !saveBranch("Containers",conts(), &RootTool::getEntry).isSuccess() )
if ( !saveBranch("Links", links(), &RootTool::getEntry).isSuccess() )
if ( !saveBranch("Params", params(),&RootTool::getParam).isSuccess() )
}
}

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

Generated at Wed Jan 30 2013 17:13:50 for Gaudi Framework, version v23r6 by Doxygen version 1.8.2 written by Dimitri van Heesch, © 1997-2004