Project
Loading...
Searching...
No Matches
GeneratorHepMC.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13
18#include "SimConfig/SimConfig.h"
19#include "HepMC3/ReaderFactory.h"
20#include "HepMC3/GenEvent.h"
21#include "HepMC3/GenParticle.h"
22#include "HepMC3/GenVertex.h"
23#include "HepMC3/FourVector.h"
24#include "HepMC3/Version.h"
25#include "TParticle.h"
26
27#include <fairlogger/Logger.h>
28#include "FairPrimaryGenerator.h"
29#include <cmath>
30#include <sstream>
31
32namespace o2
33{
34namespace eventgen
35{
36
37/*****************************************************************/
38/*****************************************************************/
39
41 : GeneratorHepMC("ALICEo2", "ALICEo2 HepMC Generator")
42{
43}
44
45/*****************************************************************/
46
47GeneratorHepMC::GeneratorHepMC(const Char_t* name, const Char_t* title)
48 : Generator(name, title)
49{
52 mEvent = new HepMC3::GenEvent();
53 mInterface = reinterpret_cast<void*>(mEvent);
54 mInterfaceName = "hepmc";
55}
56
57/*****************************************************************/
58
60{
62 LOG(info) << "Destructing GeneratorHepMC";
63 if (mReader) {
64 mReader->close();
65 }
66 if (mEvent) {
67 delete mEvent;
68 }
69 if (not mCmd.empty()) {
70 // Must be executed before removing the temporary file
71 // otherwise the current child process might still be writing on it
72 // causing unwanted stdout messages which could slow down the system
74 }
75 removeTemp();
76}
77
78/*****************************************************************/
81 const conf::SimConfig& config)
82{
83 if (not param.fileName.empty()) {
84 LOG(warn) << "The use of the key \"HepMC.fileName\" is "
85 << "deprecated, use \"GeneratorFileOrCmd.fileNames\" instead";
86 }
87
88 GeneratorFileOrCmd::setup(param0, config);
89 if (not param.fileName.empty()) {
90 setFileNames(param.fileName);
91 }
92
93 mVersion = param.version;
94 mPrune = param.prune;
95 setEventsToSkip(param.eventsToSkip);
96
97 // we are skipping ahead in the HepMC stream now
98 for (int i = 0; i < mEventsToSkip; ++i) {
100 }
101
102 if (param.version != 0 and mCmd.empty()) {
103 LOG(warn) << "The key \"HepMC.version\" is no longer needed when "
104 << "reading from files. The format version of the input files "
105 << "are automatically deduced. However, it is mandatory when reading "
106 << "from a pipe containing HepMC2 data.";
107 }
108}
109
110/*****************************************************************/
112 const HepMCGenConfig& param,
113 const conf::SimConfig& config)
114{
115 if (not param.fileName.empty()) {
116 LOG(warn) << "The use of the key \"HepMC.fileName\" is "
117 << "deprecated, use \"GeneratorFileOrCmd.fileNames\" instead";
118 }
119
120 GeneratorFileOrCmd::setup(param0, config);
121 if (not param.fileName.empty()) {
122 setFileNames(param.fileName);
123 }
124
125 mVersion = param.version;
126 mPrune = param.prune;
127 setEventsToSkip(param.eventsToSkip);
128
129 // we are skipping ahead in the HepMC stream now
130 for (int i = 0; i < mEventsToSkip; ++i) {
132 }
133
134 if (param.version != 0 and mCmd.empty()) {
135 LOG(warn) << "The key \"HepMC.version\" is no longer needed when "
136 << "reading from files. The format version of the input files "
137 << "are automatically deduced. However, it is mandatory when reading "
138 << "from a pipe containing HepMC2 data.";
139 }
140}
141
142/*****************************************************************/
144{
145 LOG(debug) << "Generating an event";
147 int tries = 0;
148 constexpr int max_tries = 3;
149 do {
150 LOG(debug) << " try # " << tries;
151 if (not mReader and not makeReader()) {
152 return false;
153 }
154
156 mEvent->clear();
157 mReader->read_event(*mEvent);
158 if (not mReader->failed()) {
160 mEvent->set_units(HepMC3::Units::GEV, HepMC3::Units::MM);
161 LOG(debug) << "Read one event " << mEvent->event_number();
162 return true;
163 } else {
164 LOG(error) << "Event reading from HepMC failed ...";
165 }
166 tries++;
167 } while (tries < max_tries);
168
169 LOG(error) << "HepMC event gen failed (Does the file/stream have enough events)?";
170
172 return false;
173}
174
175/*****************************************************************/
177{
178 HepMC3::GenEvent& event = *mEvent;
179
180 auto particles = event.particles();
181 auto vertices = event.vertices();
182 std::list<HepMC3::GenParticlePtr> toRemove;
183
184 LOG(debug) << "HepMC events has " << particles.size()
185 << " particles and " << vertices.size()
186 << " vertices" << std::endl;
187
188 size_t nSelect = 0;
189 for (size_t i = 0; i < particles.size(); ++i) {
190 auto particle = particles[i];
191 if (select(particle)) {
192 nSelect++;
193 continue;
194 }
195
196 // Remove particle from the event
197 toRemove.push_back(particle);
198 LOG(debug) << " Remove " << std::setw(3) << particle->id();
199
200 auto endVtx = particle->end_vertex();
201 auto prdVtx = particle->production_vertex();
202 if (endVtx) {
203 // Disconnect this particle from its out going vertex
204 endVtx->remove_particle_in(particle);
205 LOG(debug) << " end " << std::setw(3) << endVtx->id();
206
207 if (prdVtx and prdVtx->id() != endVtx->id()) {
208 auto outbound = endVtx->particles_out();
209 auto inbound = endVtx->particles_in();
210 LOG(debug) << " prd " << std::setw(3) << prdVtx->id() << " "
211 << std::setw(3) << outbound.size() << " out "
212 << " "
213 << std::setw(3) << inbound.size() << " in ";
214
215 // Other out-bound particles of the end vertex are attached as
216 // out-going to the production vertex of this particle.
217 for (auto outgoing : outbound) {
218 // This should also detach the particle from its old
219 // end-vertex.
220 if (outgoing) {
221 auto ee = outgoing->end_vertex();
222 if (not ee or ee->id() != prdVtx->id()) {
223 prdVtx->add_particle_out(outgoing);
224 }
225 LOG(debug) << " " << std::setw(3) << outgoing->id();
226 }
227 }
228
229 // Other incoming particles to the end vertex of this
230 // particles are attached incoming particles to the production
231 // vertex of this particle.
232 for (auto incoming : inbound) {
233 if (incoming) {
234 auto pp = incoming->production_vertex();
235 if (not pp or pp->id() != prdVtx->id()) {
236 prdVtx->add_particle_in(incoming);
237 }
238
239 LOG(debug) << " " << std::setw(3) << incoming->id();
240 }
241 }
242 }
243 }
244 if (prdVtx) {
245 prdVtx->remove_particle_out(particle);
246 }
247 }
248
249 LOG(debug) << "Selected " << nSelect << " particles\n"
250 << "Removing " << toRemove.size() << " particles";
251 size_t oldSize = particles.size();
252 for (auto particle : toRemove) {
253 event.remove_particle(particle);
254 }
255
256 std::list<HepMC3::GenVertexPtr> remVtx;
257 for (auto vtx : event.vertices()) {
258 if (not vtx or
259 (vtx->particles_out().empty() and
260 vtx->particles_in().empty())) {
261 remVtx.push_back(vtx);
262 }
263 }
264 LOG(debug) << "Removing " << remVtx.size() << " vertexes";
265 for (auto vtx : remVtx) {
266 event.remove_vertex(vtx);
267 }
268
269 LOG(debug) << "HepMC events was pruned from " << oldSize
270 << " particles to " << event.particles().size()
271 << " particles and " << event.vertices().size()
272 << " vertices";
273}
274
275/*****************************************************************/
276
278{
280 if (mPrune) {
281 auto select = [](HepMC3::ConstGenParticlePtr particle) {
282 switch (particle->status()) {
283 case 1: // Final st
284 case 2: // Decayed
285 case 4: // Beam
286 return true;
287 }
288 // To also keep diffractive particles
289 // if (particle->pid() == 9902210) return true;
290 return false;
291 };
293 }
294
296 mParticles.clear();
297 auto particles = mEvent->particles();
298 for (int i = 0; i < particles.size(); ++i) {
299
301 auto particle = particles.at(i);
302 auto momentum = particle->momentum();
303 auto vertex = particle->production_vertex()->position();
304 auto parents = particle->parents();
305 auto children = particle->children();
306
308 auto m1 = parents.empty() ? -1 : parents.front()->id() - 1;
309 auto m2 = parents.empty() ? -1 : parents.back()->id() - 1;
310
312 auto d1 = children.empty() ? -1 : children.front()->id() - 1;
313 auto d2 = children.empty() ? -1 : children.back()->id() - 1;
314
316 mParticles.push_back(TParticle(particle->pid(), // Particle type
317 particle->status(), // Status code
318 m1, // First mother
319 m2, // Second mother
320 d1, // First daughter
321 d2, // Last daughter
322 momentum.x(), // X-momentum
323 momentum.y(), // Y-momentum
324 momentum.z(), // Z-momentum
325 momentum.t(), // Energy
326 vertex.x(), // Production X
327 vertex.y(), // Production Y
328 vertex.z(), // Production Z
329 vertex.t())); // Production time
331 mParticles.back(), // Add to back
332 particle->status() == 1); // only final state are to be propagated
333
334 }
337 return kTRUE;
338}
339
340namespace
341{
342template <typename AttributeType, typename TargetType>
343bool putAttributeInfoImpl(o2::dataformats::MCEventHeader* eventHeader,
344 const std::string& name,
345 const std::shared_ptr<HepMC3::Attribute>& a)
346{
347 if (auto* p = dynamic_cast<AttributeType*>(a.get())) {
348 eventHeader->putInfo<TargetType>(name, p->value());
349 return true;
350 }
351 return false;
352}
353
354void putAttributeInfo(o2::dataformats::MCEventHeader* eventHeader,
355 const std::string& name,
356 const std::shared_ptr<HepMC3::Attribute>& a)
357{
358 using IntAttribute = HepMC3::IntAttribute;
359 using LongAttribute = HepMC3::LongAttribute;
360 using FloatAttribute = HepMC3::FloatAttribute;
361 using DoubleAttribute = HepMC3::DoubleAttribute;
362 using StringAttribute = HepMC3::StringAttribute;
363 using CharAttribute = HepMC3::CharAttribute;
364 using LongLongAttribute = HepMC3::LongLongAttribute;
365 using LongDoubleAttribute = HepMC3::LongDoubleAttribute;
366 using UIntAttribute = HepMC3::UIntAttribute;
367 using ULongAttribute = HepMC3::ULongAttribute;
368 using ULongLongAttribute = HepMC3::ULongLongAttribute;
369 using BoolAttribute = HepMC3::BoolAttribute;
370
371 if (putAttributeInfoImpl<IntAttribute, int>(eventHeader, name, a)) {
372 return;
373 }
374 if (putAttributeInfoImpl<LongAttribute, int>(eventHeader, name, a)) {
375 return;
376 }
377 if (putAttributeInfoImpl<FloatAttribute, float>(eventHeader, name, a)) {
378 return;
379 }
380 if (putAttributeInfoImpl<DoubleAttribute, float>(eventHeader, name, a)) {
381 return;
382 }
383 if (putAttributeInfoImpl<StringAttribute, std::string>(eventHeader, name, a)) {
384 return;
385 }
386 if (putAttributeInfoImpl<CharAttribute, char>(eventHeader, name, a)) {
387 return;
388 }
389 if (putAttributeInfoImpl<LongLongAttribute, int>(eventHeader, name, a)) {
390 return;
391 }
392 if (putAttributeInfoImpl<LongDoubleAttribute, float>(eventHeader, name, a)) {
393 return;
394 }
395 if (putAttributeInfoImpl<UIntAttribute, int>(eventHeader, name, a)) {
396 return;
397 }
398 if (putAttributeInfoImpl<ULongAttribute, int>(eventHeader, name, a)) {
399 return;
400 }
401 if (putAttributeInfoImpl<ULongLongAttribute, int>(eventHeader, name, a)) {
402 return;
403 }
404 if (putAttributeInfoImpl<BoolAttribute, bool>(eventHeader, name, a)) {
405 return;
406 }
407}
408} // namespace
409
410/*****************************************************************/
411
413{
416
417 eventHeader->putInfo<std::string>(Key::generator, "hepmc");
418 eventHeader->putInfo<int>(Key::generatorVersion, HEPMC3_VERSION_CODE);
419
420 auto xSection = mEvent->cross_section();
421 auto pdfInfo = mEvent->pdf_info();
422 auto hiInfo = mEvent->heavy_ion();
423
424 // Workaround for a bug in HepMC3 (3.3.1 on 23/02/2026): GenHeavyIon::from_string() for the "v0"
425 // format skips reading user_cent_estimate, but to_string() always writes it.
426 // This shifts all subsequent fields by one, causing a istringstream failure and and heavy_ion()
427 // to return null even when the attribute is present and well-formed.
428 // For now we use this manual parser in case the infos are available
429 if (!hiInfo) {
430 auto attStr = mEvent->attribute_as_string("GenHeavyIon");
431 if (!attStr.empty() && attStr[0] == 'v') {
432 std::istringstream is(attStr);
433 std::string version;
434 is >> version;
435 if (version == "v0") {
436 auto hi = std::make_shared<HepMC3::GenHeavyIon>();
437 double spectNeutrons, spectProtons, eccentricity, userCentEst;
438 is >> hi->Ncoll_hard >> hi->Npart_proj >> hi->Npart_targ >> hi->Ncoll >> spectNeutrons >> spectProtons // deprecated v0 fields
439 >> hi->N_Nwounded_collisions >> hi->Nwounded_N_collisions >> hi->Nwounded_Nwounded_collisions >> hi->impact_parameter >> hi->event_plane_angle >> eccentricity // deprecated v0 field
440 >> hi->sigma_inel_NN >> hi->centrality >> userCentEst // GenHeavyIon::to_string always writes this, but GenHeavyIon::from_string skips it for v0 (HepMC3 bug to fix)
441 >> hi->Nspec_proj_n >> hi->Nspec_targ_n >> hi->Nspec_proj_p >> hi->Nspec_targ_p;
442 if (!is.fail()) {
443 LOG(debug) << "GenHeavyIon: using manual v0 parser (workaround for HepMC3 from_string bug)";
444 hiInfo = hi;
445 } else {
446 LOG(warn) << "GenHeavyIon: manual v0 parser also failed on: [" << attStr << "]";
447 }
448 }
449 }
450 }
451
452 // Set default cross-section
453 if (xSection) {
454 eventHeader->putInfo<float>(Key::xSection, xSection->xsec());
455 eventHeader->putInfo<float>(Key::xSectionError, xSection->xsec_err());
456 eventHeader->putInfo<int>(Key::acceptedEvents,
457 xSection->get_accepted_events());
458 eventHeader->putInfo<int>(Key::attemptedEvents,
459 xSection->get_attempted_events());
460 }
461
462 // Set weights and cross sections
463 size_t iw = 0;
464 for (auto w : mEvent->weights()) {
465 std::string post = (iw > 0 ? "_" + std::to_string(iw) : "");
466 eventHeader->putInfo<float>(Key::weight + post, w);
467 if (xSection) {
468 eventHeader->putInfo<float>(Key::xSection, xSection->xsec(iw));
469 eventHeader->putInfo<float>(Key::xSectionError, xSection->xsec_err(iw));
470 }
471 iw++;
472 }
473
474 // Set the PDF information
475 if (pdfInfo) {
476 eventHeader->putInfo<int>(Key::pdfParton1Id, pdfInfo->parton_id[0]);
477 eventHeader->putInfo<int>(Key::pdfParton2Id, pdfInfo->parton_id[1]);
478 eventHeader->putInfo<float>(Key::pdfX1, pdfInfo->x[0]);
479 eventHeader->putInfo<float>(Key::pdfX2, pdfInfo->x[1]);
480 eventHeader->putInfo<float>(Key::pdfScale, pdfInfo->scale);
481 eventHeader->putInfo<float>(Key::pdfXF1, pdfInfo->xf[0]);
482 eventHeader->putInfo<float>(Key::pdfXF2, pdfInfo->xf[1]);
483 eventHeader->putInfo<int>(Key::pdfCode1, pdfInfo->pdf_id[0]);
484 eventHeader->putInfo<int>(Key::pdfCode2, pdfInfo->pdf_id[1]);
485 }
486
487 // Set heavy-ion information
488 if (hiInfo) {
489 eventHeader->SetB(hiInfo->impact_parameter); // sets the impact parameter to the FairMCEventHeader field for quick access in the AO2D
490 eventHeader->putInfo<float>(Key::impactParameter,
491 hiInfo->impact_parameter);
492 eventHeader->putInfo<int>(Key::nPart,
493 hiInfo->Npart_proj + hiInfo->Npart_targ);
494 eventHeader->putInfo<int>(Key::nPartProjectile, hiInfo->Npart_proj);
495 eventHeader->putInfo<int>(Key::nPartTarget, hiInfo->Npart_targ);
496 eventHeader->putInfo<int>(Key::nColl, hiInfo->Ncoll);
497 eventHeader->putInfo<int>(Key::nCollHard, hiInfo->Ncoll_hard);
498 eventHeader->putInfo<int>(Key::nCollNNWounded,
499 hiInfo->N_Nwounded_collisions);
500 eventHeader->putInfo<int>(Key::nCollNWoundedN,
501 hiInfo->Nwounded_N_collisions);
502 eventHeader->putInfo<int>(Key::nCollNWoundedNwounded,
503 hiInfo->Nwounded_Nwounded_collisions);
504 eventHeader->putInfo<double>(Key::planeAngle, hiInfo->event_plane_angle);
505 eventHeader->putInfo<float>(Key::sigmaInelNN, hiInfo->sigma_inel_NN);
506 eventHeader->putInfo<float>(Key::centrality, hiInfo->centrality);
507 eventHeader->putInfo<int>(Key::nSpecProjectileProton, hiInfo->Nspec_proj_p);
508 eventHeader->putInfo<int>(Key::nSpecProjectileNeutron, hiInfo->Nspec_proj_n);
509 eventHeader->putInfo<int>(Key::nSpecTargetProton, hiInfo->Nspec_targ_p);
510 eventHeader->putInfo<int>(Key::nSpecTargetNeutron, hiInfo->Nspec_targ_n);
511 }
512
513 for (auto na : mEvent->attributes()) {
514 std::string name = na.first;
515 if (name == "GenPdfInfo" ||
516 name == "GenCrossSection" ||
517 name == "GenHeavyIon") {
518 continue;
519 }
520
521 for (auto ia : na.second) {
522 int no = ia.first;
523 auto at = ia.second;
524 std::string post = (no == 0 ? "" : std::to_string(no));
525
526 putAttributeInfo(eventHeader, name + post, at);
527 }
528 }
529}
530
531/*****************************************************************/
532
534{
535 // Reset the reader smart pointer
536 LOG(debug) << "Reseting the reader";
537 mReader.reset();
538
539 // Check that we have any file names left
540 if (mFileNames.size() < 1) {
541 LOG(debug) << "No more files to read, return false";
542 return false;
543 }
544
545 // If we have file names left, pop the top of the list (LIFO)
546 auto filename = mFileNames.front();
547 mFileNames.pop_front();
548
549 LOG(debug) << "Next file to read: \"" << filename << "\" "
550 << mFileNames.size() << " left";
551
552 if (not mCmd.empty()) {
553 // For FIFO reading, we assume straight ASCII output always.
554 // Unfortunately, the HepMC3::deduce_reader `stat`s the filename
555 // which isn't supported on a FIFO, so we have to use the reader
556 // directly. Here, we allow for version 2 formats if the user
557 // specifies that
558 LOG(info) << "Creating ASCII reader of " << filename;
559 if (mVersion == 2) {
560 mReader = std::make_shared<HepMC3::ReaderAsciiHepMC2>(filename);
561 } else {
562 mReader = std::make_shared<HepMC3::ReaderAscii>(filename);
563 }
564 } else {
565 LOG(info) << "Deduce a reader of " << filename;
566 mReader = HepMC3::deduce_reader(filename);
567 }
568
569 bool ret = bool(mReader) and not mReader->failed();
570 LOG(info) << "Reader is " << mReader.get() << " " << ret;
571 return ret;
572}
573
574/*****************************************************************/
575
577{
582
583 // If a EG command line is given, then we make a fifo on a temporary
584 // file, and directs the EG to write to that fifo. We will then set
585 // up the HepMC3 reader to read from that fifo.
586 //
587 // o2-sim -g hepmc --configKeyValues "HepMC.progCmd=<cmd>" ...
588 //
589 // where <cmd> is the command line to run an event generator. The
590 // event generator should output HepMC event records to standard
591 // output. Nothing else, but the HepMC event record may be output
592 // to standard output. If the EG has other output to standard
593 // output, then a filter can be set-up. For example
594 //
595 // crmc -n 3 -o hepmc3 -c /optsw/inst/etc/crmc.param -f /dev/stdout \
596 // | sed -n 's/^\‍(HepMC::\|[EAUWVP] \‍)/\1/p'
597 //
598 // What's more, the event generator program _must_ accept the
599 // following command line argument
600 //
601 // `-n NEVENTS` to set the number of events to produce.
602 //
603 // Optionally, the command line should also accept
604 //
605 // `-s SEED` to set the random number seed
606 // `-b FM` to set the maximum impact parameter to sample
607 // `-o OUTPUT` to set the output file name
608 //
609 // All of this can conviniently be achieved via a wrapper script
610 // around the actual EG program.
611 if (not mCmd.empty()) {
612 if (mFileNames.empty()) {
613 // Set filename to be a temporary name
614 if (not makeTemp(false)) {
615 return false;
616 }
617 } else {
618 // Use the first filename as output for cmd line
619 if (not makeTemp(true)) {
620 return false;
621 }
622 }
623
624 // Make a fifo
625 if (not makeFifo()) {
626 return false;
627 }
628
629 // Build command line, rediret stdout to our fifo and put
630 std::string cmd = makeCmdLine();
631 LOG(debug) << "EG command line is \"" << cmd << "\"";
632
633 // Execute the command line
634 if (not executeCmdLine(cmd)) {
635 LOG(fatal) << "Failed to spawn \"" << cmd << "\"";
636 return false;
637 }
638 } else {
639 // If no command line was given, ensure that all files are present
640 // on the system. Note, in principle, HepMC3 can read from remote
641 // files
642 //
643 // root:// XRootD served
644 // http[s]:// Web served
645 // gsidcap:// DCap served
646 //
647 // These will all be handled in HepMC3 via ROOT's TFile protocol
648 // and the files are assumed to contain a TTree named
649 // `hepmc3_tree` and that tree has the branches
650 //
651 // `hepmc3_event` with object of type `HepMC3::GenEventData`
652 // `GenRunInfo` with object of type `HepMC3::GenRunInfoData`
653 //
654 // where the last branch is optional.
655 //
656 // However, here we will assume system local files. If _any_ of
657 // the listed files do not exist, then we fail.
658 if (not ensureFiles()) {
659 return false;
660 }
661 }
662
663 // Create reader for current (first) file
664 return true;
665}
666
667/*****************************************************************/
668/*****************************************************************/
669
670} /* namespace eventgen */
671} /* namespace o2 */
o2::monitoring::tags::Key Key
std::ostringstream debug
uint64_t vertex
Definition RawEventData.h:9
int32_t i
Utility functions for MC particles.
uint32_t version
Definition RawData.h:8
void putInfo(std::string const &key, T const &value)
void updateHeader(o2::dataformats::MCEventHeader *eventHeader) override
void setEventsToSkip(uint64_t val)
std::shared_ptr< HepMC3::Reader > mReader
void setup(const GeneratorFileOrCmdParam &param0, const GeneratorHepMCParam &param, const conf::SimConfig &config)
std::string mInterfaceName
Definition Generator.h:127
std::vector< TParticle > mParticles
Definition Generator.h:147
Bool_t Init() override
static void encodeParticleStatusAndTracking(TParticle &particle, bool wanttracking=true)
Definition MCUtils.cxx:211
GLuint const GLchar * name
Definition glcorearb.h:781
GLenum GLfloat param
Definition glcorearb.h:271
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLubyte GLubyte GLubyte GLubyte w
Definition glcorearb.h:852
std::vector< InputSpec > select(char const *matcher="")
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string to_string(gsl::span< T, Size > span)
Definition common.h:52
std::string filename()
void setFileNames(const std::string &filenames)
void setup(const GeneratorFileOrCmdParam &param, const conf::SimConfig &config)
std::list< std::string > mFileNames
virtual std::string makeCmdLine() const
virtual bool executeCmdLine(const std::string &cmd)
virtual bool makeTemp(const bool &)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"