|
Gaudi Framework, version v23r2 |
| Home | Generated: Thu Jun 28 2012 |
Description: More...
#include <src/RootTool.h>


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. | |
Description:
Concrete implementation to read objects from POOL files.
Concrete implementation to read objects from ROOT files.
Definition at line 16 of file RootTool.h.
| Gaudi::RootTool::RootTool | ( | RootDataConnection * | con ) | [inline] |
| void Gaudi::RootTool::addEntry | ( | StringVec & | c, |
| char * | val | ||
| ) | [inline] |
Helper function to read string tables.
Definition at line 108 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 100 of file RootTool.h.
| void Gaudi::RootTool::analyzeMergeMap | ( | StringVec & | tmp ) | [inline] |
Build merge sections from the Sections table entries.
Definition at line 152 of file RootTool.h.
{
StringVec::const_iterator i;
LinkSections& ls = linkSections();
MergeSections& ms = mergeSections();
pair<string,ContainerSection> e;
MsgStream& msg = msgSvc();
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 132 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;
}
Access data branch by name: Get existing branch in read only mode.
Definition at line 21 of file RootTool.h.
Helper function to convert string vectors to branch entries.
Definition at line 213 of file RootTool.h.
{ return c; }
Helper function to convert parameter vectors to branch entries.
Definition at line 215 of file RootTool.h.
{ return p.first+"="+p.second; }
| virtual int Gaudi::RootTool::loadRefs | ( | CSTR | section, |
| CSTR | cnt, | ||
| unsigned long | entry, | ||
| RootObjectRefs & | refs | ||
| ) | [inline, virtual] |
Load references object from file.
Link manager:
Definition at line 31 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() ) {
MsgStream& msg = msgSvc();
msgSvc() << MSG::VERBOSE;
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;
}
| StatusCode Gaudi::RootTool::readBranch | ( | TTree * | t, |
| const char * | nam, | ||
| C & | v, | ||
| F | pmf | ||
| ) | [inline] |
Helper function to read internal file tables.
Definition at line 110 of file RootTool.h.
{
char text[2048];
TBranch* b = t->GetBranch(nam);
if ( b ) {
TLeaf* l = b->GetLeaf(nam);
if ( l ) {
b->SetAddress(text);
msgSvc() << MSG::VERBOSE;
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);
}
}
return StatusCode::SUCCESS;
}
}
msgSvc() << MSG::ERROR << "Failed to read '" << nam << "' table." << endmsg;
return StatusCode::FAILURE;
}
| StatusCode Gaudi::RootTool::readRefs | ( | ) | [inline, virtual] |
Read reference tables.
Implements Gaudi::RootDataConnection::Tool.
Definition at line 193 of file RootTool.h.
{
TTree* t = (TTree*)c->file()->Get("Sections");
StringVec tmp;
if ( t && !readBranch(t, "Sections", tmp, &RootTool::addEntry).isSuccess() )
return StatusCode::FAILURE;
else if ( refs() ) {
analyzeMergeMap(tmp);
if ( !readBranch(refs(),"Databases", dbs(), &RootTool::addEntry).isSuccess() )
return StatusCode::FAILURE;
if ( !readBranch(refs(),"Containers",conts(), &RootTool::addEntry).isSuccess() )
return StatusCode::FAILURE;
if ( !readBranch(refs(),"Links", links(), &RootTool::addEntry).isSuccess() )
return StatusCode::FAILURE;
if ( !readBranch(refs(),"Params", params(), &RootTool::addParam).isSuccess() )
return StatusCode::FAILURE;
return StatusCode::SUCCESS;
}
return StatusCode::FAILURE;
}
| StatusCode Gaudi::RootTool::saveBranch | ( | const char * | nam, |
| C & | v, | ||
| F | pmf | ||
| ) | [inline] |
Helper function to save internal tables.
Definition at line 217 of file RootTool.h.
{
Long64_t i, n;
string val, typ = nam;
StatusCode sc = StatusCode::SUCCESS;
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() <= 1) sc = StatusCode::FAILURE;
}
return sc;
}
return StatusCode::FAILURE;
}
| StatusCode Gaudi::RootTool::saveRefs | ( | ) | [inline, virtual] |
Save/update reference tables.
Implements Gaudi::RootDataConnection::Tool.
Definition at line 236 of file RootTool.h.
{
if ( refs() ) {
if ( !saveBranch("Databases", dbs(), &RootTool::getEntry).isSuccess() )
return StatusCode::FAILURE;
if ( !saveBranch("Containers",conts(), &RootTool::getEntry).isSuccess() )
return StatusCode::FAILURE;
if ( !saveBranch("Links", links(), &RootTool::getEntry).isSuccess() )
return StatusCode::FAILURE;
if ( !saveBranch("Params", params(),&RootTool::getParam).isSuccess() )
return StatusCode::FAILURE;
return StatusCode::SUCCESS;
}
return StatusCode::FAILURE;
}