The Gaudi Framework  v30r3 (a5ef0a68)
RootEvtSelector.cpp
Go to the documentation of this file.
1 //====================================================================
2 // RootEvtSelector.cpp
3 //--------------------------------------------------------------------
4 //
5 // Package : RootCnv
6 //
7 // Author : M.Frank
8 //====================================================================
9 #ifndef GAUDIROOTCNV_ROOTEVTSELECTORCONTEXT_H
10 #define GAUDIROOTCNV_ROOTEVTSELECTORCONTEXT_H
11 
12 // Include files
14 #include <vector>
15 
16 // Forward declarations
17 class TBranch;
18 
19 /*
20 * Gaudi namespace declaration
21 */
22 namespace Gaudi
23 {
24 
34  {
35  public:
38 
39  private:
43  Files m_files;
45  Files::const_iterator m_fiter;
47  long m_entry;
49  TBranch* m_branch;
52 
53  public:
55  RootEvtSelectorContext( const RootEvtSelector* s ) : m_sel( s ), m_entry( -1 ), m_branch( nullptr ) {}
57  const Files& files() const { return m_files; }
59  void setFiles( const Files& f )
60  {
61  m_files = f;
62  m_fiter = m_files.begin();
63  }
65  void* identifier() const override { return const_cast<RootEvtSelector*>( m_sel ); }
67  Files::const_iterator fileIterator() const { return m_fiter; }
69  void setFileIterator( Files::const_iterator i ) { m_fiter = i; }
71  long entry() const { return m_entry; }
73  void setEntry( long e ) { m_entry = e; }
75  void setFID( const std::string& fid ) { m_fid = fid; }
77  const std::string& fid() const { return m_fid; }
79  TBranch* branch() const { return m_branch; }
81  void setBranch( TBranch* b ) { m_branch = b; }
82  };
83 }
84 #endif // GAUDIROOTCNV_ROOTEVTSELECTORCONTEXT_H
85 
86 // Include files
88 #include "GaudiKernel/ClassID.h"
92 #include "GaudiKernel/MsgStream.h"
94 #include "RootCnv/RootCnvSvc.h"
96 #include "TBranch.h"
97 
98 using namespace Gaudi;
99 using namespace std;
100 
101 // Helper method to issue error messages
102 StatusCode RootEvtSelector::error( const string& msg ) const
103 {
104  MsgStream log( msgSvc(), name() );
105  log << MSG::ERROR << msg << endmsg;
106  return StatusCode::FAILURE;
107 }
108 
109 // IService implementation: Db event selector override
111 {
112  // Initialize base class
113  StatusCode status = Service::initialize();
114  if ( !status.isSuccess() ) {
115  return error( "Error initializing base class Service!" );
116  }
117 
118  auto ipers = serviceLocator()->service<IPersistencySvc>( m_persName );
119  if ( !ipers ) {
120  return error( "Unable to locate IPersistencySvc interface of " + m_persName );
121  }
122  IConversionSvc* cnvSvc = nullptr;
123  Gaudi::Utils::TypeNameString itm( m_cnvSvcName );
124  status = ipers->getService( itm.name(), cnvSvc );
125  if ( !status.isSuccess() ) {
126  status = ipers->getService( itm.type(), cnvSvc );
127  if ( !status.isSuccess() ) {
128  return error( "Unable to locate IConversionSvc interface of database type " + m_cnvSvcName );
129  }
130  }
131  m_dbMgr = dynamic_cast<RootCnvSvc*>( cnvSvc );
132  if ( !m_dbMgr ) {
133  cnvSvc->release();
134  return error( "Unable to localize service:" + m_cnvSvcName );
135  }
136  m_dbMgr->addRef();
137 
138  // Get DataSvc
139  auto eds = serviceLocator()->service<IDataManagerSvc>( "EventDataSvc" );
140  if ( !eds ) {
141  return error( "Unable to localize service EventDataSvc" );
142  }
143  m_rootCLID = eds->rootCLID();
144  m_rootName = eds->rootName();
145  MsgStream log( msgSvc(), name() );
146  log << MSG::DEBUG << "Selection root:" << m_rootName << " CLID:" << m_rootCLID << endmsg;
147  return status;
148 }
149 
150 // IService implementation: Service finalization
152 {
153  // Initialize base class
154  if ( m_dbMgr ) m_dbMgr->release();
155  m_dbMgr = nullptr; // release
156  return Service::finalize();
157 }
158 
159 // Create a new event loop context
160 StatusCode RootEvtSelector::createContext( Context*& refpCtxt ) const
161 {
162  refpCtxt = new RootEvtSelectorContext( this );
163  return StatusCode::SUCCESS;
164 }
165 
166 // Access last item in the iteration
167 StatusCode RootEvtSelector::last( Context& /*refContext*/ ) const { return StatusCode::FAILURE; }
168 
169 // Get next iteration item from the event loop context
170 StatusCode RootEvtSelector::next( Context& ctxt ) const
171 {
172  RootEvtSelectorContext* pCtxt = dynamic_cast<RootEvtSelectorContext*>( &ctxt );
173  if ( pCtxt ) {
174  TBranch* b = pCtxt->branch();
175  if ( !b ) {
176  auto fileit = pCtxt->fileIterator();
177  pCtxt->setBranch( nullptr );
178  pCtxt->setEntry( -1 );
179  if ( fileit != pCtxt->files().end() ) {
180  RootDataConnection* con = nullptr;
181  string in = *fileit;
182  StatusCode sc = m_dbMgr->connectDatabase( in, IDataConnection::READ, &con );
183  if ( sc.isSuccess() ) {
184  string section = m_rootName[0] == '/' ? m_rootName.substr( 1 ) : m_rootName;
185  b = con->getBranch( section, m_rootName );
186  if ( b ) {
187  pCtxt->setFID( con->fid() );
188  pCtxt->setBranch( b );
189  return next( ctxt );
190  }
191  }
192  m_dbMgr->disconnect( in ).ignore();
193  pCtxt->setFileIterator( ++fileit );
194  return next( ctxt );
195  }
196  return StatusCode::FAILURE;
197  }
198  long ent = pCtxt->entry();
199  Long64_t nent = b->GetEntries();
200  if ( nent > ( ent + 1 ) ) {
201  pCtxt->setEntry( ++ent );
202  return StatusCode::SUCCESS;
203  }
204  auto fit = pCtxt->fileIterator();
205  pCtxt->setFileIterator( ++fit );
206  pCtxt->setEntry( -1 );
207  pCtxt->setBranch( nullptr );
208  pCtxt->setFID( "" );
209  return next( ctxt );
210  }
211  return StatusCode::FAILURE;
212 }
213 
214 // Get next iteration item from the event loop context
215 StatusCode RootEvtSelector::next( Context& ctxt, int jump ) const
216 {
217  if ( jump > 0 ) {
218  for ( int i = 0; i < jump; ++i ) {
219  StatusCode status = next( ctxt );
220  if ( !status.isSuccess() ) {
221  return status;
222  }
223  }
224  return StatusCode::SUCCESS;
225  }
226  return StatusCode::FAILURE;
227 }
228 
229 // Get previous iteration item from the event loop context
230 StatusCode RootEvtSelector::previous( Context& /* ctxt */ ) const
231 {
232  return error( "EventSelector Iterator, operator -- not supported " );
233 }
234 
235 // Get previous iteration item from the event loop context
236 StatusCode RootEvtSelector::previous( Context& ctxt, int jump ) const
237 {
238  if ( jump > 0 ) {
239  for ( int i = 0; i < jump; ++i ) {
240  StatusCode status = previous( ctxt );
241  if ( !status.isSuccess() ) {
242  return status;
243  }
244  }
245  return StatusCode::SUCCESS;
246  }
247  return StatusCode::FAILURE;
248 }
249 
250 // Rewind the dataset
251 StatusCode RootEvtSelector::rewind( Context& ctxt ) const
252 {
253  RootEvtSelectorContext* pCtxt = dynamic_cast<RootEvtSelectorContext*>( &ctxt );
254  if ( pCtxt ) {
255  auto fileit = pCtxt->fileIterator();
256  if ( fileit != pCtxt->files().end() ) {
257  string input = *fileit;
258  m_dbMgr->disconnect( input ).ignore();
259  }
260  pCtxt->setFID( "" );
261  pCtxt->setEntry( -1 );
262  pCtxt->setBranch( nullptr );
263  pCtxt->setFileIterator( pCtxt->files().begin() );
264  return StatusCode::SUCCESS;
265  }
266  return StatusCode::FAILURE;
267 }
268 
269 // Create new Opaque address corresponding to the current record
270 StatusCode RootEvtSelector::createAddress( const Context& ctxt, IOpaqueAddress*& pAddr ) const
271 {
272  const RootEvtSelectorContext* pctxt = dynamic_cast<const RootEvtSelectorContext*>( &ctxt );
273  if ( pctxt ) {
274  long ent = pctxt->entry();
275  if ( ent >= 0 ) {
276  auto fileit = pctxt->fileIterator();
277  if ( fileit != pctxt->files().end() ) {
278  const string par[2] = {pctxt->fid(), m_rootName};
279  const unsigned long ipar[2] = {0, (unsigned long)ent};
280  return m_dbMgr->createAddress( m_dbMgr->repSvcType(), m_rootCLID, &par[0], &ipar[0], pAddr );
281  }
282  }
283  }
284  pAddr = nullptr;
285  return StatusCode::FAILURE;
286 }
287 
288 // Release existing event iteration context
290 {
291  RootEvtSelectorContext* pCtxt = dynamic_cast<RootEvtSelectorContext*>( ctxt );
292  if ( pCtxt ) {
293  delete pCtxt;
294  return StatusCode::SUCCESS;
295  }
296  return StatusCode::FAILURE;
297 }
298 
299 // Will set a new criteria for the selection of the next list of events and will change
300 // the state of the context in a way to point to the new list.
301 StatusCode RootEvtSelector::resetCriteria( const string& criteria, Context& context ) const
302 {
303  MsgStream log( msgSvc(), name() );
304  RootEvtSelectorContext* ctxt = dynamic_cast<RootEvtSelectorContext*>( &context );
305  string db, typ, item, sel, stmt, aut, addr;
306  if ( ctxt ) {
307  if ( criteria.compare( 0, 5, "FILE " ) == 0 ) {
308  // The format for the criteria is:
309  // FILE filename1, filename2 ...
310  db = criteria.substr( 5 );
311  } else {
312  using Parser = Gaudi::Utils::AttribStringParser;
313  for ( auto attrib : Parser( criteria ) ) {
314  string tmp = attrib.tag.substr( 0, 3 );
315  if ( tmp == "DAT" ) {
316  db = std::move( attrib.value );
317  } else if ( tmp == "OPT" ) {
318  if ( attrib.value.compare( 0, 3, "REA" ) != 0 ) {
319  log << MSG::ERROR << "Option:\"" << attrib.value << "\" not valid" << endmsg;
320  return StatusCode::FAILURE;
321  }
322  } else if ( tmp == "TYP" ) {
323  typ = std::move( attrib.value );
324  } else if ( tmp == "ADD" ) {
325  item = std::move( attrib.value );
326  } else if ( tmp == "SEL" ) {
327  sel = std::move( attrib.value );
328  } else if ( tmp == "FUN" ) {
329  stmt = std::move( attrib.value );
330  } else if ( tmp == "AUT" ) {
331  aut = std::move( attrib.value );
332  } else if ( tmp == "COL" ) {
333  addr = std::move( attrib.value );
334  }
335  }
336  }
337  // It's now time to parse the criteria for the event selection
338  // The format for the criteria is:
339  // FILE filename1, filename2 ...
340  // JOBID number1-number2, number3, ...
342  string rest = db;
343  files.clear();
344  while ( true ) {
345  int ipos = rest.find_first_not_of( " ," );
346  if ( ipos == -1 ) break;
347  rest = rest.substr( ipos, string::npos ); // remove blanks before
348  int lpos = rest.find_first_of( " ," ); // locate next blank
349  files.push_back( rest.substr( 0, lpos ) ); // insert in list
350  if ( lpos == -1 ) break;
351  rest = rest.substr( lpos, string::npos ); // get the rest
352  }
353  ctxt->setFiles( files );
354  ctxt->setFileIterator( ctxt->files().begin() );
355  return StatusCode::SUCCESS;
356  }
357  return error( "Invalid iteration context." );
358 }
Parse attribute strings allowing iteration over the various attributes.
StatusCode resetCriteria(const std::string &cr, Context &c) const override
Will set a new criteria for the selection of the next list of events and will change the state of the...
constexpr static const auto FAILURE
Definition: StatusCode.h:88
Definition of the MsgStream class used to transmit messages.
Definition: MsgStream.h:24
StatusCode initialize() override
Definition: Service.cpp:63
const Files & files() const
Access to the file container.
Files::const_iterator fileIterator() const
Access to the file iterator.
const std::string & fid() const
Access file id.
RootEvtSelectorContext(const RootEvtSelector *s)
Standard constructor with initialization.
T find_first_not_of(T...args)
StatusCode finalize() override
Definition: Service.cpp:173
std::vector< std::string > Files
Definition of the file container.
bool isSuccess() const
Definition: StatusCode.h:287
std::string m_fid
Connection fid.
sel
Definition: IOTest.py:95
Description:
Definition: RootCnvSvc.h:53
STL namespace.
T end(T...args)
StatusCode last(Context &refContext) const override
Access last item in the iteration.
STL class.
StatusCode releaseContext(Context *&refCtxt) const override
Release existing event iteration context.
T push_back(T...args)
MsgStream & error() const
shortcut for the method msgStream(MSG::ERROR)
Helper class to parse a string of format "type/name".
StatusCode rewind(Context &refCtxt) const override
Rewind the dataset.
ROOT specific event selector context.
T next(T...args)
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:51
StatusCode finalize() override
IService implementation: Service finalization.
T find_first_of(T...args)
void * identifier() const override
Context identifier.
T clear(T...args)
T move(T...args)
const RootEvtSelector * m_sel
Reference to the hosting event selector instance.
constexpr static const auto SUCCESS
Definition: StatusCode.h:87
void setFiles(const Files &f)
Set the file container.
void setBranch(TBranch *b)
Set the top level branch (typically /Event) used to iterate.
StatusCode next(Context &refCtxt) const override
Get next iteration item from the event loop context.
TBranch * getBranch(boost::string_ref section, boost::string_ref branch_name)
Access data branch by name: Get existing branch in read only mode.
Concrete event selector implementation to access ROOT files.
void setEntry(long e)
Set current event entry number.
virtual unsigned long release()=0
Release Interface instance.
StatusCode createAddress(const Context &refCtxt, IOpaqueAddress *&) const override
Create new Opaque address corresponding to the current record.
T begin(T...args)
const std::string & type() const
const std::string & fid() const
Access connection fid.
void setFID(const std::string &fid)
Set connection FID.
Files::const_iterator m_fiter
The iterator to the.
string s
Definition: gaudirun.py:253
T substr(T...args)
Data persistency service interface.
void setFileIterator(Files::const_iterator i)
Set file iterator.
TBranch * branch() const
Access to the top level branch (typically /Event) used to iterate.
StatusCode createContext(Context *&refpCtxt) const override
Create a new event loop context.
const std::string & name() const
Opaque address interface definition.
StatusCode initialize() override
IService implementation: Db event selector override.
Files m_files
The file container managed by this context.
Concrete implementation of the IDataConnection interface to access ROOT files.
TBranch * m_branch
Reference to the top level branch (typically /Event) used to iterate.
StatusCode previous(Context &refCtxt) const override
Get previous iteration item from the event loop context.
Helper functions to set/get the application return code.
Definition: __init__.py:1
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:209
long entry() const
Access to the current event entry number.
T compare(T...args)
long m_entry
Current entry of current file.