All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
RootTool.h
Go to the documentation of this file.
1 #include <algorithm>
2 /*
3 * Gaudi namespace declaration
4 */
5 namespace Gaudi {
6 
16  class RootTool : virtual public RootDataConnection::Tool {
17  public:
19  RootTool(RootDataConnection* con) { c = con; }
21  virtual TBranch* getBranch(CSTR section, CSTR branch_name) {
22  std::string n = branch_name+".";
23  for(int i=0, m=n.length()-1; i<m; ++i) if ( !isalnum(n[i]) ) n[i]='_';
24  TTree* t = c->getSection(section);
25  TBranch* b = t ? t->GetBranch(n.c_str()) : 0;
26  if ( !b ) b = t ? t->GetBranch(branch_name.c_str()) : 0;
27  if ( b ) b->SetAutoDelete(kFALSE);
28  return b;
29  }
31  virtual int loadRefs(CSTR section, CSTR cnt, unsigned long entry, RootObjectRefs& refs) {
32  TBranch* b = getBranch(section,cnt+"#R");
33  RootObjectRefs* prefs = &refs;
34  if ( b ) {
35  b->SetAddress(&prefs);
36  int nb = b->GetEntry(entry);
37  if ( nb >= 1 ) {
38  const MergeSections& ms = c->mergeSections();
39  if ( !ms.empty() ) {
40  MsgStream& msg = msgSvc();
41  msgSvc() << MSG::VERBOSE;
42  pair<const RootRef*,const ContainerSection*> ls = c->getMergeSection(cnt,entry);
43  if ( ls.first ) {
44  if ( ls.first->dbase >= 0 ) {
45  // Now patch the references and links 'en block' to be efficient
46  // First the leafs from the TES
47  if ( msg.isActive() ) {
48  msg << "Refs: LS [" << entry << "] -> "
49  << ls.first->dbase << "," << ls.first->container
50  << "," << ls.first->link
51  << "," << ls.first->entry
52  << endmsg;
53  }
54  for(size_t j=0, n=refs.refs.size(); j<n; ++j) {
55  RootRef& r = refs.refs[j];
56  if ( r.entry>= 0 && r.dbase >= 0 ) {
57  int db = r.dbase + ls.first->dbase;
58  if ( c->getDb(db) == c->fid() ) {
59  if ( r.dbase ) r.dbase += ls.first->dbase;
60  if ( r.container ) r.container += ls.first->container;
61  if ( r.link ) r.link += ls.first->link;
62  const string& rc = c->getCont(r.container);
63  auto k = ms.find(rc);
64  if ( k != ms.end() ) {
65  const auto& cs = (*k).second;
66  r.entry = ( ls.first->entry >= 0 && ls.first->entry < (int)cs.size() )
67  ? cs[ls.first->entry].start + r.entry : -1;
68  if ( msg.isActive() ) {
69  msg << "Add link [" << r.entry
70  << "," << ls.first->entry
71  << "," << ls.first->container
72  << "," << r.container
73  << "," << r.entry
74  << "] to -> " << rc << endmsg;
75  }
76  }
77  else {
78  msg << MSG::WARNING << c->fid() << " [" << c->getDb(db)
79  << "] Evt:" << entry << " Invalid link to " << rc << endmsg;
80  msg << MSG::VERBOSE;
81  }
82  }
83  }
84  }
86  std::for_each(std::begin(refs.links),std::end(refs.links),
87  [&](int& i) { i += ls.first->link; } );
88  }
89  return nb;
90  }
91  return -1;
92  }
93  return nb;
94  }
95  }
96  return -1;
97  }
99  void addParam(ParamMap& c, char* p) {
100  char* q = strchr(p,'=');
101  if ( q ) {
102  *q = 0;
103  c.emplace_back(p,++q);
104  }
105  }
107  void addEntry(StringVec& c, char* val) { c.push_back(val); }
109  template <class C, class F> StatusCode readBranch(TTree* t, const char* nam, C& v, F pmf) {
110  char text[2048];
111  TBranch* b = t->GetBranch(nam);
112  if ( b ) {
113  TLeaf* l = b->GetLeaf(nam);
114  if ( l ) {
116  b->SetAddress(text);
117  msgSvc() << MSG::VERBOSE;
118  for(Long64_t i=0, n=b->GetEntries(); i<n; ++i) {
119  if ( b->GetEntry(i)>0 ) {
120  char* p = (char*)l->GetValuePointer();
121  msgSvc() << "Add Value[" << b->GetName() << "]:" << p << endmsg;
122  (this->*pmf)(v,p);
123  }
124  else {
126  }
127  }
128  return sc;
129  }
130  }
131  msgSvc() << MSG::ERROR << "Failed to read '" << nam << "' table." << endmsg;
133  }
135  bool get(const string& dsc, pair<string,ContainerSection>& e) {
136  if ( dsc != "[END-OF-SECTION]" ) {
137  size_t id1 = dsc.find("[CNT=");
138  size_t id2 = dsc.find("[START=");
139  size_t id3 = dsc.find("[LEN=");
140  if ( id1 != string::npos && id2 != string::npos && id3 != string::npos ) {
141  string cnt = dsc.substr(id1+5, id2-1-5);
142  int section_start = std::stoi(dsc.substr(id2+7,id3-id2-8));
143  int section_length = std::stoi(dsc.substr(id3+5,dsc.find("]",id3+5)-id3-5));
144  e.first = cnt;
145  e.second = ContainerSection(section_start,section_length);
146  return true;
147  }
148  }
149  e.first.clear();
150  e.second = ContainerSection(-1,-1);
151  return false;
152  }
155  LinkSections& ls = linkSections();
157  pair<string,ContainerSection> e;
158  MsgStream& msg = msgSvc();
159  RootRef r;
160  int cnt = 0;
161  ls.clear();
162  ms.clear();
163  msg << MSG::VERBOSE;
164  r.dbase = r.container = r.link = r.clid = r.svc = r.entry = 0;
165  for(const auto& i:tmp) {
166  if ( get(i,e) ) {
167  msg << "Added Merge Section:" << e.first << endmsg;
168  ms[e.first].push_back(e.second);
169  if ( e.first == "Links" )
170  r.link = e.second.start;
171  else if ( e.first == "Containers" )
172  r.container = e.second.start;
173  else if ( e.first == "Databases" )
174  r.dbase = e.second.start;
175  else if ( e.first == "Params" )
176  r.svc = e.second.start;
177  }
178  else if ( i == "[END-OF-SECTION]" ) {
179  r.entry = cnt;
180  if ( msg.isActive() ) {
181  msg << "Link Section [" << r.entry << "," << ls.size()
182  << "] -> D:" << r.dbase
183  << " C:" << r.container
184  << " L:" << r.link
185  << " P:" << r.svc
186  << endmsg;
187  }
188  ls.push_back(r);
189  cnt++;
190  }
191  }
192  }
195  TTree* t = (TTree*)c->file()->Get("Sections");
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() )
203  return sc;
204  if ( !(sc=readBranch(refs(),"Containers",conts(), &RootTool::addEntry)).isSuccess() )
205  return sc;
206  if ( !(sc=readBranch(refs(),"Links", links(), &RootTool::addEntry)).isSuccess() )
207  return sc;
208  if ( !(sc=readBranch(refs(),"Params", params(), &RootTool::addParam)).isSuccess() )
209  return sc;
210  return sc;
211  }
212  return StatusCode::FAILURE;
213  }
215  string getEntry(const string& c) { return c; }
217  string getParam(const pair<string,string>& p) { return p.first+"="+p.second; }
219  template <class C, class F> StatusCode saveBranch(const char* nam, C& v, F pmf) {
220  Long64_t i, n;
221  string val, typ = nam;
223  TDirectory::TContext ctxt(c->file());
224  TBranch* b = refs()->GetBranch(nam);
225  if ( !b ) b = refs()->Branch(nam,0,(typ+"/C").c_str());
226  if ( b ) {
227  for(i=b->GetEntries(), n=Long64_t(v.size()); i<n; ++i) {
228  val = (this->*pmf)(v[size_t(i)]);
229  b->SetAddress((char*)val.c_str());
230  msgSvc() << MSG::VERBOSE << "Save Value[" << b->GetName() << "]:" << val << endmsg;
231  if ( b->Fill() < 0 ) sc = StatusCode::FAILURE;
232  }
233  return sc;
234  }
235  return StatusCode::FAILURE;
236  }
239  if ( refs() ) {
240  if ( !saveBranch("Databases", dbs(), &RootTool::getEntry).isSuccess() )
241  return StatusCode::FAILURE;
242  if ( !saveBranch("Containers",conts(), &RootTool::getEntry).isSuccess() )
243  return StatusCode::FAILURE;
244  if ( !saveBranch("Links", links(), &RootTool::getEntry).isSuccess() )
245  return StatusCode::FAILURE;
246  if ( !saveBranch("Params", params(),&RootTool::getParam).isSuccess() )
247  return StatusCode::FAILURE;
248  return StatusCode::SUCCESS;
249  }
250  return StatusCode::FAILURE;
251  }
252  };
253 }
Definition of the MsgStream class used to transmit messages.
Definition: MsgStream.h:24
const MergeSections & mergeSections() const
Access merged data section inventory.
virtual int loadRefs(CSTR section, CSTR cnt, unsigned long entry, RootObjectRefs &refs)
Load references object from file.
Definition: RootTool.h:31
const std::string & fid() const
Access file id.
TFile * file() const
Direct access to TFile structure.
const std::string & getDb(int which) const
Access database/file name from saved index.
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:244
auto begin(reverse_wrapper< T > &w)
Definition: reverse.h:45
RootDataConnection::StringVec StringVec
bool isActive() const
Accessor: is MsgStream active.
Definition: MsgStream.h:128
string getEntry(const string &c)
Helper function to convert string vectors to branch entries.
Definition: RootTool.h:215
std::vector< int > links
The links of the link manager.
Definition: RootRefs.h:73
Persistent reference object containing all leafs and links corresponding to a Gaudi DataObject...
Definition: RootRefs.h:71
int dbase
Data members to define object location in the persistent world.
Definition: RootRefs.h:32
MergeSections & mergeSections() const
Persistent reference object.
Definition: RootRefs.h:30
void analyzeMergeMap(StringVec &tmp)
Build merge sections from the Sections table entries.
Definition: RootTool.h:154
Description:
Definition: RootTool.h:16
StatusCode saveBranch(const char *nam, C &v, F pmf)
Helper function to save internal tables.
Definition: RootTool.h:219
auto end(reverse_wrapper< T > &w)
Definition: reverse.h:47
Helper class to facilitate an abstraction layer for reading POOL style files with this package...
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:26
constexpr double m
Definition: SystemOfUnits.h:93
RootDataConnection::ParamMap ParamMap
std::pair< const RootRef *, const ContainerSection * > getMergeSection(const std::string &container, int entry) const
Access link section for single container and entry.
const std::string CSTR
RootDataConnection::MergeSections MergeSections
string getParam(const pair< string, string > &p)
Helper function to convert parameter vectors to branch entries.
Definition: RootTool.h:217
tuple rc
Definition: IOTest.py:92
const std::string & getCont(int which) const
Access container name from saved index.
StatusCode saveRefs()
Save/update reference tables.
Definition: RootTool.h:238
dictionary l
Definition: gaudirun.py:422
void addParam(ParamMap &c, char *p)
Helper function to read params table.
Definition: RootTool.h:99
void addEntry(StringVec &c, char *val)
Helper function to read string tables.
Definition: RootTool.h:107
StatusCode readBranch(TTree *t, const char *nam, C &v, F pmf)
Helper function to read internal file tables.
Definition: RootTool.h:109
StatusCode readRefs()
Read reference tables.
Definition: RootTool.h:194
RootDataConnection::ContainerSection ContainerSection
RootDataConnection::LinkSections LinkSections
RootTool(RootDataConnection *con)
Standard constructor.
Definition: RootTool.h:19
RootDataConnection * c
Pointer to containing data connection object.
list i
Definition: ana.py:128
Concrete implementation of the IDataConnection interface to access ROOT files.
std::vector< RootRef > refs
The references corresponding to the next layer of items in the data store.
Definition: RootRefs.h:75
constexpr double ms
Helper functions to set/get the application return code.
Definition: __init__.py:1
virtual TBranch * getBranch(CSTR section, CSTR branch_name)
Access data branch by name: Get existing branch in read only mode.
Definition: RootTool.h:21
TTree * getSection(const std::string &sect, bool create=false)
Access TTree section from section name. The section is created if required.
LinkSections & linkSections() const