ATLAS Offline Software
Loading...
Searching...
No Matches
d3pdReaderMaker.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6// STL include(s):
7#include <iostream>
8#include <vector>
9#include <string>
10
11// Boost include(s):
12#include <boost/program_options.hpp>
13#include <boost/regex.hpp>
14
15// ROOT include(s):
16#include <TFile.h>
17#include <TTree.h>
18#include <TList.h>
19#include <TKey.h>
20#include <TBranch.h>
21#include <TBranchElement.h>
22#include <TLeaf.h>
23#include <TError.h>
24
25// Gaudi/Athena include(s):
26#include "GaudiKernel/StatusCode.h"
27#include "GaudiKernel/Bootstrap.h"
28#include "GaudiKernel/IMessageSvc.h"
29#include "GaudiKernel/ISvcLocator.h"
32
33// Local include(s):
34#include "../CodeGenerator_v2.h"
36
38static const char* const PROGRAM_DESC =
39 "This application can be used to generate D3PDReader classes\n"
40 "based on an existing D3PD file and some additional command\n"
41 "line options";
42
44static const char* const PROGRAM_GREETING =
45 "*******************************************************************\n"
46 "* *\n"
47 "* D3PDReader standalone code generator *\n"
48 "* *\n"
49 "*******************************************************************";
50
52
62template< typename T >
63std::ostream& operator<<( std::ostream& out, const std::vector< T >& vec ) {
64 out << "[";
65 typename std::vector< T >::const_iterator itr = vec.begin();
66 typename std::vector< T >::const_iterator end = vec.end();
67 for( ; itr != end; ++itr ) {
68 out << *itr;
69 if( ++itr != end ) {
70 out << ", ";
71 }
72 --itr;
73 }
74 out << "]";
75 return out;
76}
77
79StatusCode extractVariables( const std::string& file_name, const std::string& tree_name,
80 const std::vector< std::string >& patterns,
81 D3PD::RootObjectMetadata& metadata );
82StatusCode extractVariables( const std::string& file_name, const std::string& tree_name,
83 const std::string& pattern,
84 D3PD::RootObjectMetadata& metadata );
85
87namespace po = boost::program_options;
88//coverity[root_function]
89int main ATLAS_NOT_THREAD_SAFE ( int argc, char* argv[] ) {
90 // Let's disable the ROOT warnings:
91 gErrorIgnoreLevel = kError;
92
93 // Force the Gaudi MessageSvc into existence. This gets rid of the
94 // warnings from getMessageSvc().
95 SmartIF<IMessageSvc> msvc{Gaudi::svcLocator()->service("MessageSvc")};
96 if( !msvc ) {
97 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
98 << "Couldn't set up the Gaudi message service";
99 return 255;
100 }
101
102 //
103 // Construct the object describing all the command line parameters of the
104 // application:
105 //
106 po::options_description desc( PROGRAM_DESC );
107 desc.add_options()
108 ( "help,h", "Give some help with the program usage" )
109 ( "classname,n", po::value< std::string >()->default_value( "D3PDReader" ),
110 "Generated class name" )
111 ( "d3pd-files,f", po::value< std::vector< std::string > >()->multitoken(),
112 "Input D3PD file(s)" )
113 ( "tree,t", po::value< std::string >()->default_value( "" ), "Name of the D3PD tree" )
114 ( "variables,v", po::value< std::vector< std::string > >()->multitoken(),
115 "Variable names to consider (regular expression(s))" )
116 ( "prefixes,p", po::value< std::vector< std::string > >()->multitoken(),
117 "Variable name prefix(es)" )
118 ( "container,c", "Option specifying that a proxy class should be created" )
119 ( "output,o", po::value< std::string >()->default_value( "." ),
120 "Output directory for the generated source files" );
121
122 //
123 // Interpret the command line options provided by the user:
124 //
125 po::variables_map vm;
126 try {
127 po::store( po::parse_command_line( argc, argv, desc ), vm );
128 po::notify( vm );
129 } catch( const std::exception& ex ) {
130 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
131 << "There was a problem with interpreting the command line options.";
132 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
133 << "Message: " << ex.what();
134 return 255;
135 }
136
137 //
138 // If the user asked for some help, let's print it:
139 //
140 if( vm.count( "help" ) ) {
141 std::cout << desc << std::endl;
142 return 0;
143 }
144
145 //
146 // Greet the user, and print all the parameters that the application
147 // received:
148 //
149 std::cout << PROGRAM_GREETING << std::endl;
150 std::cout << "Program will be executed with the options:" << std::endl;
151
152 // Get the name of the class to be generated:
153 std::string class_name = vm[ "classname" ].as< std::string >();
154 std::cout << " Generating class with name: \"D3PDReader::" << class_name << "\""
155 << std::endl;
156
157 // Extract the file name(s) from the command line arguments:
158 if( ! vm.count( "d3pd-files" ) ) {
159 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
160 << "You have to specify at least one D3PD file!";
161 return 255;
162 }
163 const std::vector< std::string > file_names =
164 vm[ "d3pd-files" ].as< std::vector< std::string > >();
165 std::cout << " D3PD file(s): " << file_names << std::endl;
166
167 // Get the name of the D3PD TTree:
168 const std::string d3pd_tree_name = vm[ "tree" ].as< std::string >();
169 if( d3pd_tree_name != "" ) {
170 std::cout << " Using D3PD tree: " << d3pd_tree_name << std::endl;
171 } else {
172 std::cout << " Trying to find D3PD tree automatically" << std::endl;
173 }
174
175 // Get the variable name patterns:
176 std::vector< std::string > var_names;
177 if( vm.count( "variables" ) ) {
178 var_names = vm[ "variables" ].as< std::vector< std::string > >();
179 std::cout << " Only considering variables with names: " << var_names << std::endl;
180 } else {
181 std::cout << " Considering all variables from the input file" << std::endl;
182 }
183
184 // Get the common prefix of the variable names:
185 if( ( ( vm.count( "prefixes" ) != vm.count( "variables" ) ) &&
186 ( vm.count( "prefixes" ) != 1 ) ) ||
187 ( vm.count( "prefixes" ) == 0 ) ) {
188 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
189 << "You have to either specify just one prefix, or the same number as "
190 << "how many variable regular expressions you gave";
191 return 255;
192 }
193 const std::vector< std::string > prefixes =
194 vm[ "prefixes" ].as< std::vector< std::string > >();
195 std::cout << " Common variable prefix(es): " << prefixes << std::endl;
196
197 // Get the parameter describing whether proxy classes should be generated:
198 const bool is_container = vm.count( "container" );
199 std::cout << " Proxy class will " << ( is_container ? "" : "NOT " ) << "be created"
200 << std::endl;
201
202 // Get the output directory:
203 const std::string output = vm[ "output" ].as< std::string >();
204 std::cout << " Output directory: " << output << std::endl;
205
206 //
207 // Create the metadata object(s):
208 //
209 std::map< std::string, D3PD::RootObjectMetadata > metadata;
210 std::vector< std::string >::const_iterator prefix_itr = prefixes.begin();
211 std::vector< std::string >::const_iterator prefix_end = prefixes.end();
212 for( ; prefix_itr != prefix_end; ++prefix_itr ) {
213 metadata[ *prefix_itr ].setName( class_name );
214 metadata[ *prefix_itr ].setPrefix( *prefix_itr );
215 metadata[ *prefix_itr ].setContainer( is_container );
216 }
217
218 //
219 // Iterate over the files:
220 //
221 std::vector< std::string >::const_iterator file_itr = file_names.begin();
222 std::vector< std::string >::const_iterator file_end = file_names.end();
223 for( ; file_itr != file_end; ++file_itr ) {
224
225 if( prefixes.size() > 1 ) {
226
227 for( size_t i = 0; i < prefixes.size(); ++i ) {
228 if( extractVariables( *file_itr, d3pd_tree_name, var_names[ i ],
229 metadata[ prefixes[ i ] ] ).isFailure() ) {
230 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
231 << "Couldn't extract the variable descriptions from file: "
232 << *file_itr;
233 return 255;
234 }
235 }
236
237 } else {
238
239 if( extractVariables( *file_itr, d3pd_tree_name, var_names,
240 metadata[ prefixes[ 0 ] ] ).isFailure() ) {
241 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
242 << "Couldn't extract the variable descriptions from file: "
243 << *file_itr;
244 return 255;
245 }
246
247 }
248
249 }
250
251 //
252 // Merge the information extracted:
253 //
254 D3PD::ObjectMetadata summed_meta;
255 summed_meta.setName( class_name );
256 summed_meta.setPrefix( prefixes[ 0 ] );
257 summed_meta.setContainer( is_container );
258 prefix_itr = prefixes.begin();
259 prefix_end = prefixes.end();
260 for( ; prefix_itr != prefix_end; ++prefix_itr ) {
261 summed_meta += metadata[ *prefix_itr ];
262 }
263
264 //
265 // Generate the D3PDObjectBase class's source:
266 //
267 if( D3PD::Version2::writeD3PDObjectBase( output ).isFailure() ) {
268 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
269 << "Couldn't create the D3PDObjectBase class source";
270 return 255;
271 }
272
273 //
274 // Generate the VarHandle class's source:
275 //
276 if( D3PD::Version2::writeVarHandle( output ).isFailure() ) {
277 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
278 << "Couldn't create the VarHandle class source";
279 return 255;
280 }
281
282 //
283 // Generate the VarProxy class's source if needed:
284 //
285 if( is_container ) {
286 if( D3PD::Version2::writeVarProxy( output ).isFailure() ) {
287 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
288 << "Couldn't create the VarProxy class source";
289 return 255;
290 }
291 }
292
293 //
294 // Generate the VarHandle class's source:
295 //
296 if( D3PD::Version2::writeUserD3PDObject( output ).isFailure() ) {
297 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
298 << "Couldn't create the UserD3PDObject class source";
299 return 255;
300 }
301
302 //
303 // Generate the D3PDReadStats class's source:
304 //
305 if( D3PD::Version2::writeD3PDReadStats( output ).isFailure() ) {
306 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
307 << "Couldn't create the D3PDReadStats class source";
308 return 255;
309 }
310
311 //
312 // Generate the D3PDPerfStats class's source:
313 //
314 if( D3PD::Version2::writeD3PDPerfStats( output ).isFailure() ) {
315 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
316 << "Couldn't create the D3PDPerfStats class source";
317 return 255;
318 }
319
320 //
321 // Generate the Utils source:
322 //
323 if( D3PD::Version2::writeUtils( output ).isFailure() ) {
324 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
325 << "Couldn't create the Utils source";
326 return 255;
327 }
328
329 //
330 // Generate the header of the D3PDReader class:
331 //
332 if( D3PD::Version2::writeHeader( class_name, output, summed_meta ).isFailure() ) {
333 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
334 << "Couldn't generate the D3PDReader class header";
335 return 255;
336 }
337
338 //
339 // Generate the source of the D3PDReader class:
340 //
341 if( D3PD::Version2::writeSource( class_name, output, summed_meta ).isFailure() ) {
342 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker" )
343 << "Couldn't generate the D3PDReader class source";
344 return 255;
345 }
346
347 return 0;
348}
349
361StatusCode extractVariables( const std::string& file_name, const std::string& tree_name,
362 const std::vector< std::string >& patterns,
363 D3PD::RootObjectMetadata& metadata ) {
364
365 //
366 // Try to open the D3PD file:
367 //
368 TFile* file = TFile::Open( file_name.c_str(), "READ" );
369 if( ( ! file ) || file->IsZombie() ) {
370 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
371 << "Failed to open D3PD file: " << file_name;
372 return StatusCode::FAILURE;
373 }
374
375 //
376 // Access the D3PD tree:
377 //
378 TTree* tree = 0;
379 if( tree_name != "" ) {
380 //
381 // If the user specified a tree name, then try to access that tree:
382 //
383 tree = dynamic_cast< TTree* >( file->Get( tree_name.c_str() ) );
384 if( ! tree ) {
385 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
386 << "Couldn't find tree \"" << tree_name << "\" in file: " << file_name;
387 return StatusCode::FAILURE;
388 }
389 } else {
390 //
391 // If the user didn't specify a tree name then search for all the TTree-s
392 // in the file which are not called "CollectionTree", and select the first
393 // one.
394 //
395 std::string name = "";
396 const TList* keys = file->GetListOfKeys();
397 for( Int_t i = 0; i < keys->GetSize(); ++i ) {
398 const TKey* key = dynamic_cast< TKey* >( keys->At( i ) );
399 if (!key) continue;
400 std::string key_name = key->GetName();
401 if( ( key_name.find( "CollectionTree" ) == key_name.npos ) &&
402 ( key->GetClassName() == std::string( "TTree" ) ) ) {
403 name = std::move(key_name);
404 std::cout << "Assuming that the D3PD tree is called: " << name
405 << std::endl;
406 break;
407 }
408 }
409 if( name == "" ) {
410 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
411 << "Couldn't find the name for the D3PD tree!";
412 return StatusCode::FAILURE;
413 }
414 tree = dynamic_cast< TTree* >( file->Get( name.c_str() ) );
415 if( ! tree ) {
416 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
417 << "Couldn't find tree \"" << tree_name << "\" in file: "
418 << file_name;
419 return StatusCode::FAILURE;
420 }
421 }
422
423 //
424 // Look at all the branches in the tree:
425 //
426 TObjArray* branches = tree->GetListOfBranches();
427 for( Int_t i = 0; i < branches->GetSize(); ++i ) {
428
429 //
430 // Only consider the branches that match one of the specified patterns:
431 //
432 bool pass = false;
433 std::vector< std::string >::const_iterator itr = patterns.begin();
434 std::vector< std::string >::const_iterator end = patterns.end();
435 for( ; itr != end; ++itr ) {
436 if( boost::regex_search( branches->At( i )->GetName(),
437 boost::basic_regex< char >( itr->c_str() ) ) ) {
438 pass = true;
439 break;
440 }
441 }
442 if( ! pass ) continue;
443
444 //
445 // Extract the type of the variable. This has to be done a little differently for
446 // primitive types and objects (STL containers).
447 //
448 const std::string branch_type( branches->At( i )->IsA()->GetName() );
449 if( branch_type == "TBranch" ) {
450
451 //
452 // If it's a branch describing a primitive type, then the name of the type
453 // has to be extracted from the single leaf of the branch:
454 //
455 TBranch* branch = dynamic_cast< TBranch* >( branches->At( i ) );
456 if( ! branch ) {
457 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
458 << "Could not cast object into a TBranch!";
459 return StatusCode::FAILURE;
460 }
461 TObjArray* leaves = branch->GetListOfLeaves();
462 bool leaf_found = false;
463 for( Int_t j = 0; j < leaves->GetSize(); ++j ) {
464 TLeaf* leaf = dynamic_cast< TLeaf* >( leaves->At( 0 ) );
465 if( leaf ) {
466 if( metadata.addVariable( branch->GetName(),
467 leaf->GetTypeName(),
468 branch->GetTitle() ).isFailure() ) {
469 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
470 << "Failed adding variable: " << branch->GetName();
471 return StatusCode::FAILURE;
472 }
473 leaf_found = true;
474 break;
475 }
476 }
477 if( ! leaf_found ) {
478 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
479 << "Couldn't find the TLeaf describing a primitive type!";
480 return StatusCode::FAILURE;
481 }
482
483 } else if( branch_type == "TBranchElement" ) {
484
485 //
486 // If it's a branch describing an object, then the situation is a bit simpler.
487 // In this case we can ask the TBranchElement directly for the class name.
488 // Remember however that this is a "simplified" name of the class, like
489 // "vector<float>".
490 //
491 TBranchElement* branch = dynamic_cast< TBranchElement* >( branches->At( i ) );
492 if( ! branch ) {
493 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
494 << "Could not cast object into a TBranchElement!";
495 return StatusCode::FAILURE;
496 }
497 if( metadata.addVariable( branch->GetName(),
498 branch->GetClassName(),
499 branch->GetTitle() ).isFailure() ) {
500 REPORT_MESSAGE_WITH_CONTEXT( MSG::ERROR, "d3pdReaderMaker::extractVariables" )
501 << "Failed adding variable: " << branch->GetName();
502 return StatusCode::FAILURE;
503 }
504 }
505 }
506
507 // Clean up after ourselves:
508 delete file;
509 return StatusCode::SUCCESS;
510}
511
512StatusCode extractVariables( const std::string& file_name, const std::string& tree_name,
513 const std::string& pattern,
514 D3PD::RootObjectMetadata& metadata ) {
515
516 std::vector< std::string > patterns; patterns.push_back( pattern );
517 return extractVariables( file_name, tree_name, patterns, metadata );
518}
std::vector< size_t > vec
Helpers for checking error return status codes and reporting errors.
#define REPORT_MESSAGE_WITH_CONTEXT(LVL, CONTEXT_NAME)
Report a message, with an explicitly specified context name.
int main(int, char **)
Main class for all the CppUnit test classes.
Define macros for attributes used to control the static checker.
#define ATLAS_NOT_THREAD_SAFE
getNoisyStrip() Find noisy strips from hitmaps and write out into xml/db formats
D3PD variable metadata handling class.
void setPrefix(const std::string &prefix)
Set the prefix given to variables in this D3PDObject.
void setName(const std::string &name)
Set the name of the D3PDObject that this object describes.
void setContainer(bool container)
Set whether the D3PDObject describes a container or not.
Extension of the ObjectMetadata class for reading D3PD files.
std::ostream & operator<<(std::ostream &out, const std::vector< T > &vec)
Formatted printing for vector objects.
static const char *const PROGRAM_GREETING
A greeting text that's printed when the application starts up.
StatusCode extractVariables(const std::string &file_name, const std::string &tree_name, const std::vector< std::string > &patterns, D3PD::RootObjectMetadata &metadata)
This function is used to extract the variable descriptions from the D3PD file.
static const char *const PROGRAM_DESC
The program description that's printed with the available parameters.
std::vector< std::string > patterns
Definition listroot.cxx:187
StatusCode writeUtils(const std::string &dir)
This function can be used to create source files containing some utility functions.
StatusCode writeHeader(const std::string &classname, const std::string &dir, const ObjectMetadata &metadata)
This function is used to create the header of the class describing a set of D3PD variables.
StatusCode writeD3PDReadStats(const std::string &dir)
This function can be used to create the D3PDReader::D3PDReadStats class's source files.
StatusCode writeSource(const std::string &classname, const std::string &dir, const ObjectMetadata &metadata)
This function is used to generate the source file of a D3PDReader class.
StatusCode writeVarHandle(const std::string &dir)
This function can be used to create the D3PDReader::VarHandle class's source files.
StatusCode writeVarProxy(const std::string &dir)
This function can be used to create the D3PDReader::VarProxy class's source files.
StatusCode writeD3PDObjectBase(const std::string &dir)
This function can be used to create the D3PDReader::D3PDObjectBase class's source files.
StatusCode writeUserD3PDObject(const std::string &dir)
This function can be used to create the D3PDReader::UserD3PDObject class's source files.
StatusCode writeD3PDPerfStats(const std::string &dir)
This function can be used to create the D3PDReader::D3PDPerfStats class's source files.
TChain * tree
TFile * file