DD4hep  1.30.0
Detector Description Toolkit for High Energy Physics
Geant4EventReaderHepMC.cpp
Go to the documentation of this file.
1 //==========================================================================
2 // AIDA Detector description implementation
3 //--------------------------------------------------------------------------
4 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
5 // All rights reserved.
6 //
7 // For the licensing terms see $DD4hepINSTALL/LICENSE.
8 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
9 //
10 // Author : M.Frank
11 //
12 //==========================================================================
13 
24 // Framework include files
25 #include <DDG4/IoStreams.h>
26 #include <DDG4/Geant4InputAction.h>
27 
28 // C/C++ include files
29 
31 namespace dd4hep {
32 
34  namespace sim {
35 
37  namespace HepMC {
39  class EventStream;
40  }
41 
43 
56  typedef boost::iostreams::stream<dd4hep_file_source<int> > in_stream;
57  //typedef boost::iostreams::stream<dd4hep_file_source<TFile*> > in_stream;
59  protected:
62  public:
64  explicit Geant4EventReaderHepMC(const std::string& nam);
68  virtual EventReaderStatus readParticles(int event_number,
69  Vertices& vertices,
70  std::vector<Particle*>& particles) override;
71  virtual EventReaderStatus moveToEvent(int event_number) override;
72  virtual EventReaderStatus skipEvent() override { return EVENT_READER_OK; }
73 
74  };
75  } /* End namespace sim */
76 } /* End namespace dd4hep */
77 
78 //====================================================================
79 // AIDA Detector description implementation
80 //--------------------------------------------------------------------
81 //
82 //====================================================================
83 // #include <DDG4/Geant4EventReaderHepMC"
84 
85 // Framework include files
86 #include <DDG4/Factories.h>
87 #include <DD4hep/Printout.h>
88 #include <DDG4/Geant4Primary.h>
89 #include <CLHEP/Units/SystemOfUnits.h>
90 #include <CLHEP/Units/PhysicalConstants.h>
91 
92 // C/C++ include files
93 #include <cerrno>
94 #include <climits>
95 #include <algorithm>
96 
97 using namespace dd4hep::sim;
98 using PropertyMask = dd4hep::detail::ReferenceBitMask<int>;
99 
100 // Factory entry
102 
103 namespace dd4hep {
105 
107  namespace sim {
108 
110  namespace HepMC {
111 
113  /*
114  * \author P.Kostka (main author)
115  * \author M.Frank (code reshuffeling into new DDG4 scheme)
116  * \version 1.0
117  * \ingroup DD4HEP_SIMULATION
118  */
119  class EventHeader {
120  public:
121  int id;
123  int bp1, bp2;
126  float scale;
127  float alpha_qcd;
128  float alpha_qed;
129  std::vector<float> weights;
130  std::vector<long> random;
132  EventHeader() : id(0), num_vertices(0), bp1(0), bp2(0),
133  signal_process_id(0), signal_process_vertex(0),
134  scale(0.0), alpha_qcd(0.0), alpha_qed(0.0), weights(), random() {}
135  };
136 
139 
141  /*
142  * \author P.Kostka (main author)
143  * \author M.Frank (code reshuffeling into new DDG4 scheme)
144  * \version 1.0
145  * \ingroup DD4HEP_SIMULATION
146  */
147  class EventStream {
148  public:
149  typedef std::map<int,Geant4Vertex*> Vertices;
150  typedef std::map<int,Geant4Particle*> Particles;
151 
152  std::istream& instream;
153 
154  // io information
155  std::string key;
156  double mom_unit, pos_unit;
157  int io_type;
158 
159  float xsection, xsection_err;
163 
165  EventStream(std::istream& in) : instream(in), mom_unit(0.0), pos_unit(0.0),
166  io_type(0), xsection(0.0), xsection_err(0.0)
167  { use_default_units(); }
169  bool ok() const;
171  Particles& particles() { return m_particles; }
172  Vertices& vertices() { return m_vertices; }
173  void set_io(int typ, const std::string& k)
174  { io_type = typ; key = k; }
176  { mom_unit = CLHEP::MeV; pos_unit = CLHEP::mm; }
177  bool read();
178  void clear();
179  };
180 
181  char get_input(std::istream& is, std::istringstream& iline);
182  int read_until_event_end(std::istream & is);
183  int read_weight_names(EventStream &, std::istringstream& iline);
184  int read_particle(EventStream &info, std::istringstream& iline, Geant4Particle * p);
185  int read_vertex(EventStream &info, std::istream& is, std::istringstream & iline);
186  int read_event_header(EventStream &info, std::istringstream & input, EventHeader& header);
187  int read_cross_section(EventStream &info, std::istringstream & input);
188  int read_units(EventStream &info, std::istringstream & input);
189  int read_heavy_ion(EventStream &, std::istringstream & input);
190  int read_pdf(EventStream &, std::istringstream & input);
191  Geant4Vertex* vertex(EventStream& info, int i);
192  void fix_particles(EventStream &info);
193  }
194  }
195 }
196 
197 // For debugging specify here vertex numbers
198 //#define DD4HEP_DEBUG_HEP_MC_VERTEX -390
199 //#define DD4HEP_DEBUG_HEP_MC_PARTICLE 418
200 
203  : Geant4EventReader(nam), m_input(), m_events(0)
204 {
205  // Now open the input file:
206  m_input.open(nam.c_str(), BOOST_IOS::in|BOOST_IOS::binary);
207  if ( not m_input.is_open() ) {
208  except("+++ Failed to open input stream: %s Error:%s.", nam.c_str(), ::strerror(errno));
209  }
210  m_events = new HepMC::EventStream(m_input);
211 }
212 
215  delete m_events;
216  m_events = 0;
217  m_input.close();
218 }
219 
222 Geant4EventReaderHepMC::moveToEvent(int event_number) {
223  if( m_currEvent < event_number && event_number != 0 ) {
224  printout(INFO,"EventReaderHepMC::moveToEvent","Current event:%d Skipping the next %d events",
225  m_currEvent, event_number);
226  while ( m_currEvent < event_number ) {
227  if ( not m_events->read() ) return EVENT_READER_ERROR;
228  ++m_currEvent;
229  }
230  }
231  printout(DEBUG,"EventReaderHepMC::moveToEvent","Current event number: %d",m_currEvent);
232  return EVENT_READER_OK;
233 }
234 
238  Vertices& vertices,
239  Particles& output) {
240 
241  //fg: for now we create exactly one event vertex here ( as before )
242  // this needs revisiting as HepMC allows to have more than one vertex ...
243  Geant4Vertex* primary_vertex = new Geant4Vertex ;
244  vertices.emplace_back( primary_vertex );
245  primary_vertex->x = 0;
246  primary_vertex->y = 0;
247  primary_vertex->z = 0;
248 
249  if ( !m_events->ok() ) {
250  vertices.clear();
251  output.clear();
252  return EVENT_READER_EOF;
253  }
254  else if ( m_events->read() ) {
255  EventStream::Particles& parts = m_events->particles();
256 
257  Position pos(primary_vertex->x,primary_vertex->y,primary_vertex->z);
258 
259  output.reserve(parts.size());
260  transform(parts.begin(),parts.end(),back_inserter(output),detail::reference2nd(parts));
261  m_events->clear();
262  if (pos.mag2() > std::numeric_limits<double>::epsilon() ) {
263  for(Particles::iterator k=output.begin(); k != output.end(); ++k) {
264  Geant4ParticleHandle p(*k);
265  p->vsx += pos.x();
266  p->vsy += pos.y();
267  p->vsz += pos.z();
268  p->vex += pos.x();
269  p->vey += pos.y();
270  p->vez += pos.z();
271  }
272  }
273  for(Particles::const_iterator k=output.begin(); k != output.end(); ++k) {
274  Geant4ParticleHandle p(*k);
275  printout(VERBOSE,m_name,
276  "+++ %s ID:%3d status:%08X typ:%9d Mom:(%+.2e,%+.2e,%+.2e)[MeV] "
277  "time: %+.2e [ns] #Dau:%3d #Par:%1d",
278  "",p->id,p->status,p->pdgID,
279  p->psx/CLHEP::MeV,p->psy/CLHEP::MeV,p->psz/CLHEP::MeV,p->time/CLHEP::ns,
280  p->daughters.size(),
281  p->parents.size());
282  //output.emplace_back(p);
283 
284  //add particles to the 'primary vertex'
285  if ( p->parents.size() == 0 ) {
286  PropertyMask status(p->status);
287  if ( status.isSet(G4PARTICLE_GEN_EMPTY) || status.isSet(G4PARTICLE_GEN_DOCUMENTATION) )
288  primary_vertex->in.insert(p->id); // Beam particles and primary quarks etc.
289  else
290  primary_vertex->out.insert(p->id); // Stuff, to be given to Geant4 together with daughters
291  }
292  }
293  ++m_currEvent;
294  return EVENT_READER_OK;
295  }
296  vertices.clear();
297  output.clear();
298  return EVENT_READER_EOF;
299 }
300 
302  EventStream::Particles& parts = info.particles();
303  EventStream::Vertices& verts = info.vertices();
304  EventStream::Particles::iterator i;
305  std::set<int>::const_iterator id, ip;
306  for(i=parts.begin(); i != parts.end(); ++i) {
307  Geant4ParticleHandle p((*i).second);
308  int end_vtx_id = p->secondaries;
309  p->secondaries = 0;
310  Geant4Vertex* v = vertex(info,end_vtx_id);
311 #if defined(DD4HEP_DEBUG_HEP_MC_VERTEX)
312  if ( end_vtx_id == DD4HEP_DEBUG_HEP_MC_VERTEX ) {
313  std::cout << "End-vertex:" << end_vtx_id << std::endl;
314  }
315 #endif
316  if ( v ) {
317  p->vex = v->x;
318  p->vey = v->y;
319  p->vez = v->z;
320  v->in.insert(p->id);
321  for(id=v->out.begin(); id!=v->out.end();++id) {
322  EventStream::Particles::iterator ipp = parts.find(*id);
323  Geant4Particle* dau = ipp != parts.end() ? (*ipp).second : 0;
324  if ( !dau )
325  std::cout << "ERROR: Invalid daughter particle: " << *id << std::endl;
326  else
327  dau->parents.insert(p->id);
328  p->daughters.insert(*id);
329  }
330  }
331  }
332  for(const auto& iv : verts) {
333  Geant4Vertex* v = iv.second;
334  for (int pout : v->out) {
335  EventStream::Particles::iterator ipp = parts.find(pout);
336  if ( ipp != parts.end() ) {
337  Geant4Particle* p = (*ipp).second;
338  for (int d : v->in) {
339  p->parents.insert(d);
340  }
341  }
342  }
343  }
346  std::vector<Geant4Particle*> beam;
347  for(const auto& ipart : parts) {
348  Geant4ParticleHandle p(ipart.second);
349  if ( p->parents.size() == 0 ) {
350  for(int d : p->daughters) {
351  Geant4Particle *pp = parts[d];
352  beam.emplace_back(pp);
353  }
354  }
355  }
356  for(auto* ipp : beam) {
357  //std::cout << "Clear parents of " << (*ipp)->id << std::endl;
358  ipp->parents.clear();
359  ipp->status = G4PARTICLE_GEN_DECAYED;
360  }
361 }
362 
364  EventStream::Vertices::iterator it=info.vertices().find(i);
365  return (it==info.vertices().end()) ? 0 : (*it).second;
366 }
367 
368 char HepMC::get_input(std::istream& is, std::istringstream& iline) {
369  char value = is.peek();
370  if ( !is ) { // make sure the stream is valid
371  std::cerr << "StreamHelpers: setting badbit." << std::endl;
372  is.clear(std::ios::badbit);
373  return -1;
374  }
375  std::string line, firstc;
376  getline(is,line);
377  if ( !is ) { // make sure the stream is valid
378  std::cerr << "StreamHelpers: setting badbit." << std::endl;
379  is.clear(std::ios::badbit);
380  return -1;
381  }
382  iline.clear();
383  iline.str(line);
384  if ( line.length()>0 && line[1] == ' ' ) iline >> firstc;
385  return iline ? value : -1;
386 }
387 
388 int HepMC::read_until_event_end(std::istream & is) {
389  std::string line;
390  while ( is ) {
391  char val = is.peek();
392  if( val == 'E' ) { // next event
393  return 1;
394  }
395  getline(is,line);
396  if ( !is.good() ) return 0;
397  }
398  return 0;
399 }
400 
401 int HepMC::read_weight_names(EventStream&, std::istringstream&) {
402 #if 0
403  int HepMC::read_weight_names(EventStream& info, std::istringstream& iline)
404  size_t name_size = 0;
405  iline >> name_size;
406  info.weights.names.clear();
407  info.weights.weights.clear();
408 
409  std::string name;
410  WeightContainer namedWeight;
411  std::string::size_type i1 = line.find("\""), i2, len = line.size();
412  for(size_t ii = 0; ii < name_size; ++ii) {
413  // weight names may contain blanks
414  if(i1 >= len) {
415  std::cout << "debug: attempting to read past the end of the named weight line " << std::endl;
416  std::cout << "debug: We should never get here" << std::endl;
417  std::cout << "debug: Looking for the end of this event" << std::endl;
419  }
420  i2 = line.find("\"",i1+1);
421  name = line.substr(i1+1,i2-i1-1);
422  if ( namedWeight.names.find(name) == namedWeight.names.end() )
423  namedWeight.names[name] = namedWeight.names.size();
424  namedWeight.weights[ii] = info.weights.weights[ii];
425  i1 = line.find("\"",i2+1);
426  }
427 
428  info.weights.weights = namedWeight.weights;
429  info.weights.names = namedWeight.names;
430 #endif
431  return 1;
432 }
433 
434 int HepMC::read_particle(EventStream &info, std::istringstream& input, Geant4Particle * p) {
435  float ene = 0., theta = 0., phi = 0;
436  int size = 0, stat=0;
437  PropertyMask status(p->status);
438 
439  // check that the input stream is still OK after reading item
440  input >> p->id >> p->pdgID >> p->psx >> p->psy >> p->psz >> ene;
441  p->id = info.particles().size();
442 #if defined(DD4HEP_DEBUG_HEP_MC_PARTICLE)
443  if ( p->id == DD4HEP_DEBUG_HEP_MC_PARTICLE ) {
444  std::cout << "Particle id: " << p->id << std::endl;
445  }
446 #endif
447  p->charge = 0;
448  p->psx *= info.mom_unit;
449  p->psy *= info.mom_unit;
450  p->psz *= info.mom_unit;
451  ene *= info.mom_unit;
452  if ( !input )
453  return 0;
454  else if ( info.io_type != ascii ) {
455  input >> p->mass;
456  p->mass *= info.mom_unit;
457  }
458  else {
459  p->mass = std::sqrt(fabs(ene*ene - (p->psx*p->psx + p->psy*p->psy + p->psz*p->psz)));
460  }
461  // Reuse here the secondaries to store the end-vertex ID
462  input >> stat >> theta >> phi >> p->secondaries >> size;
463  if( !input ) {
464  return 0;
465  }
466  //
467  // Generator status
468  // Simulator status 0 until simulator acts on it
469  status.clear();
470  if ( stat == 0 ) status.set(G4PARTICLE_GEN_EMPTY);
471  else if ( stat == 0x1 ) status.set(G4PARTICLE_GEN_STABLE);
472  else if ( stat == 0x2 ) status.set(G4PARTICLE_GEN_DECAYED);
473  else if ( stat == 0x3 ) status.set(G4PARTICLE_GEN_DOCUMENTATION);
474  else if ( stat == 0x4 ) status.set(G4PARTICLE_GEN_DOCUMENTATION);
475  else if ( stat == 0xB ) status.set(G4PARTICLE_GEN_DOCUMENTATION);
476  else status.set(G4PARTICLE_GEN_OTHER);
478  if ( p->secondaries != 0 ) {
479  status.set(G4PARTICLE_GEN_DECAYED);
480  }
483 
484  // read flow patterns if any exist. Protect against tainted readings.
485  size = std::min(size,100);
486  for (int i = 0; i < size; ++i ) {
487  input >> p->colorFlow[0] >> p->colorFlow[1];
488  if(!input) return 0;
489  }
490  return 1;
491 }
492 
493 int HepMC::read_vertex(EventStream &info, std::istream& is, std::istringstream & input) {
494  int id=0, dummy = 0, num_orphans_in=0, num_particles_out=0, weights_size=0;
495  std::vector<float> weights;
496  Geant4Vertex* v = new Geant4Vertex();
497  Geant4Particle* p;
498 
499  if( !input ) {
500  delete v;
501  return 0;
502  }
503  input >> id >> dummy >> v->x >> v->y >> v->z >> v->time
504  >> num_orphans_in >> num_particles_out >> weights_size;
505  if( !input ) {
506  delete v;
507  return 0;
508  }
509  if ( weights_size < 0 || weights_size > USHRT_MAX ) {
510  delete v;
511  return 0;
512  }
513 #if defined(DD4HEP_DEBUG_HEP_MC_VERTEX)
514  if ( id == DD4HEP_DEBUG_HEP_MC_VERTEX ) {
515  printout(ALWAYS,"HepMC","++ Created Vertex ID=%d",id);
516  }
517 #endif
518  v->x *= info.pos_unit;
519  v->y *= info.pos_unit;
520  v->z *= info.pos_unit;
521  for (int i1 = 0; i1 < weights_size; ++i1) {
522  float value = 0e0;
523  input >> value;
524  weights.emplace_back(value);
525  if( !input ) {
526  delete v;
527  return 0;
528  }
529  }
530  info.vertices().emplace(id,v);
531  for(char value = is.peek(); value=='P'; value=is.peek()) {
532  value = get_input(is,input);
533  if( !input || value < 0 )
534  return 0;
535 
536  read_particle(info, input,p = new Geant4Particle());
537  if( !input ) {
538  printout(ERROR,"HepMC","++ Vertex %d Failed to daughter read particle!",id);
539  delete p;
540  return 0;
541  }
542  info.particles().emplace(p->id,p);
543  p->pex = p->psx;
544  p->pey = p->psy;
545  p->pez = p->psz;
546  if ( --num_orphans_in >= 0 ) {
547  v->in.insert(p->id);
548  p->vex = v->x;
549  p->vey = v->y;
550  p->vez = v->z;
551 #if defined(DD4HEP_DEBUG_HEP_MC_VERTEX)
552  if ( id == DD4HEP_DEBUG_HEP_MC_VERTEX ) {
553  printout(ALWAYS,"HepMC","++ Vertex %d Add INGOING Particle: %d",id,p->id);
554  }
555 #endif
556  }
557  else if ( num_particles_out >= 0 ) {
558  v->out.insert(p->id);
559  p->vsx = v->x;
560  p->vsy = v->y;
561  p->vsz = v->z;
562 #if defined(DD4HEP_DEBUG_HEP_MC_VERTEX)
563  if ( id == DD4HEP_DEBUG_HEP_MC_VERTEX ) {
564  printout(ALWAYS,"HepMC","++ Vertex %d Add OUTGOING Particle: %d",id,p->id);
565  }
566 #endif
567  }
568  else {
569  delete p;
570  except("HepMC", "Invalid number of particles....");
571  }
572  }
573  return 1;
574 }
575 
576 int HepMC::read_event_header(EventStream &info, std::istringstream & input, EventHeader& header) {
577  // read values into temp variables, then fill GenEvent
578  int size = 0;
579  input >> header.id;
580  if( info.io_type == gen || info.io_type == extascii ) {
581  int nmpi = -1;
582  input >> nmpi;
583  if( input.fail() ) return 0;
584  //MSF set_mpi( nmpi );
585  }
586  input >> header.scale;
587  input >> header.alpha_qcd;
588  input >> header.alpha_qed;
589  input >> header.signal_process_id;
590  input >> header.signal_process_vertex;
591  input >> header.num_vertices;
592  if( info.io_type == gen || info.io_type == extascii )
593  input >> header.bp1 >> header.bp2;
594 
595  printout(DEBUG,"HepMC","++ Event header: %s",input.str().c_str());
596  input >> size;
597  input.clear();
598  if( input.fail() )
599  return 0;
600  if( size < 0 || size > USHRT_MAX )
601  return 0;
602 
603  for(int i = 0; i < size; ++i ) {
604  long val = 0e0;
605  input >> val;
606  header.random.emplace_back(val);
607  if( input.fail() ) return 0;
608  }
609 
610  input >> size;
611  if( input.fail() )
612  return 0;
613  if( size < 0 || size > USHRT_MAX )
614  return 0;
615 
616  std::vector<float> wgt;
617  for( int ii = 0; ii < size; ++ii ) {
618  float val = 0e0;
619  input >> val;
620  wgt.emplace_back(val);
621  if( input.fail() ) return 0;
622  }
623  // weight names will be added later if they exist
624  if( !wgt.empty() ) header.weights = std::move(wgt);
625  return 1;
626 }
627 
628 int HepMC::read_cross_section(EventStream &info, std::istringstream & input) {
629  input >> info.xsection >> info.xsection_err;
630  return input.fail() ? 0 : 1;
631 }
632 
633 int HepMC::read_units(EventStream &info, std::istringstream & input) {
634  if( info.io_type == gen ) {
635  std::string mom, pos;
636  input >> mom >> pos;
637  if ( !input.fail() ) {
638  if ( mom == "KEV" ) info.mom_unit = CLHEP::keV;
639  else if ( mom == "MEV" ) info.mom_unit = CLHEP::MeV;
640  else if ( mom == "GEV" ) info.mom_unit = CLHEP::GeV;
641  else if ( mom == "TEV" ) info.mom_unit = CLHEP::TeV;
642 
643  if ( pos == "MM" ) info.pos_unit = CLHEP::mm;
644  else if ( pos == "CM" ) info.pos_unit = CLHEP::cm;
645  else if ( pos == "M" ) info.pos_unit = CLHEP::m;
646  }
647  }
648  return input.fail() ? 0 : 1;
649 }
650 
651 int HepMC::read_heavy_ion(EventStream &, std::istringstream & input) {
652  // read values into temp variables, then create a new HeavyIon object
653  int nh =0, np =0, nt =0, nc =0,
654  neut = 0, prot = 0, nw =0, nwn =0, nwnw =0;
655  float impact = 0., plane = 0., xcen = 0., inel = 0.;
656  input >> nh >> np >> nt >> nc >> neut >> prot >> nw >> nwn >> nwnw;
657  input >> impact >> plane >> xcen >> inel;
658  /*
659  std::cerr << "Reading heavy ion, but igoring data!" << std::endl;
660  ion->set_Ncoll_hard(nh);
661  ion->set_Npart_proj(np);
662  ion->set_Npart_targ(nt);
663  ion->set_Ncoll(nc);
664  ion->set_spectator_neutrons(neut);
665  ion->set_spectator_protons(prot);
666  ion->set_N_Nwounded_collisions(nw);
667  ion->set_Nwounded_N_collisions(nwn);
668  ion->set_Nwounded_Nwounded_collisions(nwnw);
669  ion->set_impact_parameter(impact);
670  ion->set_event_plane_angle(plane);
671  ion->set_eccentricity(xcen);
672  ion->set_sigma_inel_NN(inel);
673  */
674  return input.fail() ? 0 : 1;
675 }
676 
677 int HepMC::read_pdf(EventStream &, std::istringstream & input) {
678  // read values into temp variables, then create a new PdfInfo object
679  int id1 =0, id2 =0;
680  double x1 = 0., x2 = 0., scale = 0., pdf1 = 0., pdf2 = 0.;
681  input >> id1 ;
682  if ( input.fail() )
683  return 0;
684  // check now for empty PdfInfo line
685  if( id1 == 0 )
686  return 0;
687  // continue reading
688  input >> id2 >> x1 >> x2 >> scale >> pdf1 >> pdf2;
689  if ( input.fail() )
690  return 0;
691  /*
692  std::cerr << "Reading pdf, but igoring data!" << std::endl;
693  pdf->set_id1( id1 );
694  pdf->set_id2( id2 );
695  pdf->set_x1( x1 );
696  pdf->set_x2( x2 );
697  pdf->set_scalePDF( scale );
698  pdf->set_pdf1( pdf1 );
699  pdf->set_pdf2( pdf2 );
700  */
701  // check to see if we are at the end of the line
702  if( !input.eof() ) {
703  int pdf_id1=0, pdf_id2=0;
704  input >> pdf_id1 >> pdf_id2;
705  /*
706  pdf->set_pdf_id1( pdf_id1 );
707  pdf->set_pdf_id2( pdf_id2 );
708  */
709  }
710  return input.fail() ? 0 : 1;
711 }
712 
715  // make sure the stream is good
716  if ( instream.eof() || instream.fail() ) {
717  instream.clear(std::ios::badbit);
718  return false;
719  }
720  return true;
721 }
722 
724  detail::releaseObjects(m_vertices);
725  detail::releaseObjects(m_particles);
726 }
727 
729  EventStream& info = *this;
730  bool event_read = false;
731 
732  detail::releaseObjects(vertices());
733  detail::releaseObjects(particles());
734 
735  while( instream.good() ) {
736  char value = instream.peek();
737  std::istringstream input_line;
738  if ( value == 'E' && event_read )
739  break;
740  else if ( instream.eof() && event_read )
741  goto Done;
742  else if ( instream.eof() )
743  return false;
744  else if ( value=='#' || ::isspace(value) ) {
745  get_input(instream,input_line);
746  continue;
747  }
748  value = get_input(instream,input_line);
749 
750  // On failure switch to end
751  if( !input_line || value < 0 )
752  goto Skip;
753 
754  switch( value ) {
755  case 'H': {
756  int iotype = 0;
757  std::string key_value;
758  input_line >> key_value;
759  // search for event listing key before first event only.
760  key_value = key_value.substr(0,key_value.find('\r'));
761  if ( key_value == "H" && (this->io_type == gen || this->io_type == extascii) ) {
762  read_heavy_ion(info, input_line);
763  break;
764  }
765  else if( key_value == "HepMC::IO_GenEvent-START_EVENT_LISTING" )
766  this->set_io(gen,key_value);
767  else if( key_value == "HepMC::IO_Ascii-START_EVENT_LISTING" )
768  this->set_io(ascii,key_value);
769  else if( key_value == "HepMC::IO_ExtendedAscii-START_EVENT_LISTING" )
770  this->set_io(extascii,key_value);
771  else if( key_value == "HepMC::IO_Ascii-START_PARTICLE_DATA" )
772  this->set_io(ascii_pdt,key_value);
773  else if( key_value == "HepMC::IO_ExtendedAscii-START_PARTICLE_DATA" )
774  this->set_io(extascii_pdt,key_value);
775  else if( key_value == "HepMC::IO_GenEvent-END_EVENT_LISTING" )
776  iotype = gen;
777  else if( key_value == "HepMC::IO_Ascii-END_EVENT_LISTING" )
778  iotype = ascii;
779  else if( key_value == "HepMC::IO_ExtendedAscii-END_EVENT_LISTING" )
780  iotype = extascii;
781  else if( key_value == "HepMC::IO_Ascii-END_PARTICLE_DATA" )
782  iotype = ascii_pdt;
783  else if( key_value == "HepMC::IO_ExtendedAscii-END_PARTICLE_DATA" )
784  iotype = extascii_pdt;
785 
786  if( iotype != 0 && this->io_type != iotype ) {
787  std::cerr << "GenEvent::find_end_key: iotype keys have changed. "
788  << "MALFORMED INPUT" << std::endl;
789  instream.clear(std::ios::badbit);
790  return false;
791  }
792  else if ( iotype != 0 ) {
793  break;
794  }
795  continue;
796  }
797  case 'E': // deal with the event line
798  if ( !read_event_header(info, input_line, this->header) )
799  goto Skip;
800  event_read = true;
801  continue;
802 
803  case 'N': // get weight names
804  if ( !read_weight_names(info, input_line) )
805  goto Skip;
806  continue;
807 
808  case 'U': // get unit information if it exists
809  if ( !read_units(info, input_line) )
810  goto Skip;
811  continue;
812 
813  case 'C': // we have a GenCrossSection line
814  if ( !read_cross_section(info, input_line) )
815  goto Skip;
816  continue;
817 
818  case 'V': // Read vertex with particles
819  if ( !read_vertex(info, instream, input_line) )
820  goto Skip;
821  continue;
822 
823  case 'F': // Read PDF
824  if ( !read_pdf(info, input_line) )
825  goto Skip;
826  continue;
827 
828  case 'P': // we should not find this line
829  std::cerr << "streaming input: found unexpected Particle line." << std::endl;
830  continue;
831 
832  default: // ignore everything else
833  continue;
834  }
835  continue;
836  Skip:
837  printout(WARNING,"HepMC::EventStream","+++ Skip event with ID: %d",this->header.id);
838  detail::releaseObjects(vertices());
839  detail::releaseObjects(particles());
840  read_until_event_end(instream);
841  event_read = false;
842  if ( instream.eof() ) return false;
843  }
844 
845  if( not instream.good() ) return false;
846  Done:
848  detail::releaseObjects(vertices());
849  return true;
850 }
851 
dd4hep::sim::HepMC::read_weight_names
int read_weight_names(EventStream &, std::istringstream &iline)
Definition: Geant4EventReaderHepMC.cpp:401
dd4hep::sim::HepMC::read_pdf
int read_pdf(EventStream &, std::istringstream &input)
Definition: Geant4EventReaderHepMC.cpp:677
dd4hep::sim::Geant4EventReaderHepMC::~Geant4EventReaderHepMC
virtual ~Geant4EventReaderHepMC()
Default destructor.
dd4hep::sim::HepMC::EventStream::Vertices
std::map< int, Geant4Vertex * > Vertices
Definition: Geant4EventReaderHepMC.cpp:149
dd4hep::sim::HepMC::EventStream
HepMC EventStream class used internally by the Geant4EventReaderHepMC plugin.
Definition: Geant4EventReaderHepMC.cpp:147
dd4hep::sim::Geant4Particle::vsz
double vsz
Definition: Geant4Particle.h:124
dd4hep::sim::G4PARTICLE_GEN_STABLE
@ G4PARTICLE_GEN_STABLE
Definition: Geant4Particle.h:71
dd4hep::sim::Geant4Vertex::y
double y
Definition: Geant4Vertex.h:53
v
View * v
Definition: MultiView.cpp:28
dd4hep::sim::Geant4Particle::mass
double mass
Particle mass.
Definition: Geant4Particle.h:132
dd4hep::sim::HepMC::EventStream::clear
void clear()
Definition: Geant4EventReaderHepMC.cpp:723
dd4hep::sim::Geant4Particle::vey
double vey
Definition: Geant4Particle.h:126
dd4hep::sim::Geant4EventReader::EventReaderStatus
EventReaderStatus
Status codes of the event reader object. Anything with NOT low-bit set is an error.
Definition: Geant4InputAction.h:68
dd4hep::info
std::size_t info(const std::string &src, const std::string &msg)
Definition: RootDictionary.h:65
dd4hep::sim::Geant4Particle::pez
double pez
Definition: Geant4Particle.h:130
dd4hep::sim::HepMC::get_input
char get_input(std::istream &is, std::istringstream &iline)
Definition: Geant4EventReaderHepMC.cpp:368
dd4hep::sim::HepMC::EventStream::vertex
Geant4Vertex * vertex(int i)
dd4hep::sim::Geant4Particle::id
int id
not persistent
Definition: Geant4Particle.h:108
dd4hep::sim::HepMC::read_event_header
int read_event_header(EventStream &info, std::istringstream &input, EventHeader &header)
Definition: Geant4EventReaderHepMC.cpp:576
DECLARE_GEANT4_EVENT_READER
#define DECLARE_GEANT4_EVENT_READER(name)
Plugin defintion to create event reader objects.
Definition: Factories.h:237
dd4hep::sim::HepMC::EventStream::use_default_units
void use_default_units()
Definition: Geant4EventReaderHepMC.cpp:175
dd4hep::sim::HepMC::EventHeader::bp2
int bp2
Definition: Geant4EventReaderHepMC.cpp:123
dd4hep::sim::Geant4Particle::pdgID
int pdgID
Definition: Geant4Particle.h:115
dd4hep::sim::HepMC::read_heavy_ion
int read_heavy_ion(EventStream &, std::istringstream &input)
Definition: Geant4EventReaderHepMC.cpp:651
dd4hep::sim::Geant4Particle::psz
double psz
Definition: Geant4Particle.h:128
dd4hep::sim::HepMC::EventHeader::num_vertices
int num_vertices
Definition: Geant4EventReaderHepMC.cpp:122
dd4hep::sim::Geant4EventReaderHepMC::Geant4EventReaderHepMC
Geant4EventReaderHepMC(const std::string &nam)
Initializing constructor.
dd4hep::sim::G4PARTICLE_GEN_STATUS_MASK
@ G4PARTICLE_GEN_STATUS_MASK
Definition: Geant4Particle.h:83
dd4hep::sim::HepMC::EventHeader::signal_process_id
int signal_process_id
Definition: Geant4EventReaderHepMC.cpp:124
dd4hep::sim::Geant4EventReaderHepMC::m_input
in_stream m_input
Definition: Geant4EventReaderHepMC.cpp:60
dd4hep::sim::HepMC::EventHeader::EventHeader
EventHeader()
Default constructor.
Definition: Geant4EventReaderHepMC.cpp:132
dd4hep::sim::Geant4Particle::secondaries
int secondaries
Definition: Geant4Particle.h:114
dd4hep::sim::G4PARTICLE_GEN_DOCUMENTATION
@ G4PARTICLE_GEN_DOCUMENTATION
Definition: Geant4Particle.h:73
dd4hep::sim::HepMC::ascii_pdt
@ ascii_pdt
Definition: Geant4EventReaderHepMC.cpp:138
dd4hep::sim::HepMC::EventStream::io_type
int io_type
Definition: Geant4EventReaderHepMC.cpp:157
IoStreams.h
epsilon
const double epsilon
Definition: test_cellid_position_converter.cpp:41
dd4hep::sim::HepMC::EventStream::key
std::string key
Definition: Geant4EventReaderHepMC.cpp:155
Geant4EventReaderHepMC
Plugin to read HepMC2 ASCII files.
dd4hep::sim::HepMC::fix_particles
void fix_particles(EventStream &info)
Definition: Geant4EventReaderHepMC.cpp:301
dd4hep::sim::Geant4EventReaderHepMC::in_stream
boost::iostreams::stream< dd4hep_file_source< int > > in_stream
Definition: Geant4EventReaderHepMC.cpp:56
dd4hep::sim::HepMC::EventStream::header
EventHeader header
Definition: Geant4EventReaderHepMC.cpp:160
dd4hep::sim::HepMC::extascii_pdt
@ extascii_pdt
Definition: Geant4EventReaderHepMC.cpp:138
dd4hep::sim::G4PARTICLE_GEN_EMPTY
@ G4PARTICLE_GEN_EMPTY
Definition: Geant4Particle.h:70
dd4hep::sim::HepMC::EventStream::pos_unit
double pos_unit
Definition: Geant4EventReaderHepMC.cpp:156
dd4hep::sim::HepMC::EventHeader::signal_process_vertex
int signal_process_vertex
Definition: Geant4EventReaderHepMC.cpp:125
dd4hep::sim::Geant4EventReaderHepMC::readParticles
virtual EventReaderStatus readParticles(int event_number, Vertices &vertices, std::vector< Particle * > &particles) override
Read an event and fill a vector of MCParticles.
dd4hep::sim::Geant4EventReaderHepMC::EventStream
HepMC::EventStream EventStream
Definition: Geant4EventReaderHepMC.cpp:58
dd4hep::sim::HepMC::EventStream::particles
Particles & particles()
Definition: Geant4EventReaderHepMC.cpp:171
dd4hep::sim::Geant4Particle::vez
double vez
Definition: Geant4Particle.h:126
dd4hep::sim::HepMC::EventHeader::alpha_qed
float alpha_qed
Definition: Geant4EventReaderHepMC.cpp:128
dd4hep::sim::HepMC::EventStream::xsection_err
float xsection_err
Definition: Geant4EventReaderHepMC.cpp:159
dd4hep::sim::Geant4EventReaderHepMC::m_events
EventStream * m_events
Definition: Geant4EventReaderHepMC.cpp:61
dd4hep::sim::Geant4Vertex::x
double x
Vertex data.
Definition: Geant4Vertex.h:53
dd4hep::sim::HepMC::ascii
@ ascii
Definition: Geant4EventReaderHepMC.cpp:138
dd4hep::sim::HepMC::known_io
known_io
The known_io enum is used to track which type of input is being read.
Definition: Geant4EventReaderHepMC.cpp:138
dd4hep::sim::Geant4Particle::colorFlow
int colorFlow[2]
Definition: Geant4Particle.h:117
dd4hep::sim::HepMC::extascii
@ extascii
Definition: Geant4EventReaderHepMC.cpp:138
dd4hep::sim::HepMC::EventHeader::random
std::vector< long > random
Definition: Geant4EventReaderHepMC.cpp:130
dd4hep::sim::HepMC::EventStream::Particles
std::map< int, Geant4Particle * > Particles
Definition: Geant4EventReaderHepMC.cpp:150
dd4hep::sim::Geant4Particle::genStatus
unsigned short genStatus
Definition: Geant4Particle.h:118
dd4hep::sim::HepMC::gen
@ gen
Definition: Geant4EventReaderHepMC.cpp:138
dd4hep::sim::Geant4Particle::status
int status
Definition: Geant4Particle.h:116
dd4hep::sim::HepMC::EventHeader::weights
std::vector< float > weights
Definition: Geant4EventReaderHepMC.cpp:129
dd4hep::sim::G4PARTICLE_GEN_DECAYED
@ G4PARTICLE_GEN_DECAYED
Definition: Geant4Particle.h:72
dd4hep::sim::HepMC::EventHeader::scale
float scale
Definition: Geant4EventReaderHepMC.cpp:126
dd4hep::sim::Geant4Particle::daughters
Particles daughters
The list of daughters of this MC particle.
Definition: Geant4Particle.h:140
dd4hep::sim::HepMC::EventHeader::id
int id
Definition: Geant4EventReaderHepMC.cpp:121
dd4hep::sim::HepMC::EventStream::m_vertices
Vertices m_vertices
Definition: Geant4EventReaderHepMC.cpp:161
dd4hep::sim::Geant4Particle::pex
double pex
The track momentum at the end vertex.
Definition: Geant4Particle.h:130
dd4hep::sim::Geant4EventReader
Basic geant4 event reader class. This interface/base-class must be implemented by concrete readers.
Definition: Geant4InputAction.h:60
dd4hep::sim::Geant4EventReaderHepMC::moveToEvent
virtual EventReaderStatus moveToEvent(int event_number) override
Move to the indicated event number.
dd4hep::sim::Geant4Vertex::out
Particles out
The list of outgoing particles.
Definition: Geant4Vertex.h:55
dd4hep::sim::HepMC::read_until_event_end
int read_until_event_end(std::istream &is)
Definition: Geant4EventReaderHepMC.cpp:388
dd4hep::sim::HepMC::EventStream::m_particles
Particles m_particles
Definition: Geant4EventReaderHepMC.cpp:162
dd4hep::sim::HepMC::EventStream::ok
bool ok() const
Check if data stream is in proper state and has data.
Definition: Geant4EventReaderHepMC.cpp:714
dd4hep::sim::Geant4Particle::vex
double vex
The end vertex.
Definition: Geant4Particle.h:126
dd4hep::sim::Geant4Particle::vsx
double vsx
The starting vertex.
Definition: Geant4Particle.h:124
dd4hep::sim::Geant4Particle::psx
double psx
The track momentum at the start vertex.
Definition: Geant4Particle.h:128
dd4hep::sim::Geant4Particle::pey
double pey
Definition: Geant4Particle.h:130
dd4hep::Position
ROOT::Math::XYZVector Position
Definition: Objects.h:81
dd4hep::sim::HepMC::EventHeader::bp1
int bp1
Definition: Geant4EventReaderHepMC.cpp:123
Factories.h
dd4hep::sim::HepMC::EventHeader::alpha_qcd
float alpha_qcd
Definition: Geant4EventReaderHepMC.cpp:127
key
unsigned char key
Definition: AlignmentsCalculator.cpp:69
dd4hep::sim
Namespace for the Geant4 based simulation part of the AIDA detector description toolkit.
Definition: Geant4Output2EDM4hep.cpp:49
PropertyMask
dd4hep::detail::ReferenceBitMask< int > PropertyMask
Definition: HepMC3EventReader.cpp:37
Geant4InputAction.h
Vertices
Geant4InputAction::Vertices Vertices
Definition: Geant4InputAction.cpp:27
dd4hep::sim::HepMC::EventStream::vertices
Vertices & vertices()
Definition: Geant4EventReaderHepMC.cpp:172
dd4hep::sim::Geant4EventReader::EVENT_READER_OK
@ EVENT_READER_OK
Definition: Geant4InputAction.h:70
Geant4Primary.h
dd4hep
Namespace for the AIDA detector description toolkit.
Definition: AlignmentsCalib.h:28
dd4hep::sim::Geant4ParticleHandle
Data structure to access derived MC particle information.
Definition: Geant4Particle.h:181
dd4hep::sim::Geant4Particle::parents
Particles parents
The list of parents of this MC particle.
Definition: Geant4Particle.h:138
dd4hep::sim::HepMC::read_cross_section
int read_cross_section(EventStream &info, std::istringstream &input)
Definition: Geant4EventReaderHepMC.cpp:628
dd4hep::sim::Geant4Vertex::in
Particles in
The list of incoming particles.
Definition: Geant4Vertex.h:57
dd4hep::sim::Geant4Particle::vsy
double vsy
Definition: Geant4Particle.h:124
dd4hep::sim::Geant4Vertex::z
double z
Definition: Geant4Vertex.h:53
dd4hep::sim::HepMC::vertex
Geant4Vertex * vertex(EventStream &info, int i)
Definition: Geant4EventReaderHepMC.cpp:363
dd4hep::sim::Geant4Particle
Data structure to store the MC particle information.
Definition: Geant4Particle.h:103
dd4hep::sim::Geant4Vertex
Data structure to store the MC vertex information.
Definition: Geant4Vertex.h:45
dd4hep::sim::Geant4Particle::charge
char charge
Definition: Geant4Particle.h:119
dd4hep::sim::HepMC::EventStream::EventStream
EventStream(std::istream &in)
Default constructor.
Definition: Geant4EventReaderHepMC.cpp:165
dd4hep::sim::HepMC::EventStream::instream
std::istream & instream
Definition: Geant4EventReaderHepMC.cpp:152
dd4hep::sim::HepMC::EventStream::set_io
void set_io(int typ, const std::string &k)
Definition: Geant4EventReaderHepMC.cpp:173
dd4hep::sim::Geant4EventReader::Vertices
std::vector< Vertex * > Vertices
Definition: Geant4InputAction.h:66
dd4hep::sim::Geant4Particle::psy
double psy
Definition: Geant4Particle.h:128
dd4hep::sim::HepMC::read_particle
int read_particle(EventStream &info, std::istringstream &iline, Geant4Particle *p)
Definition: Geant4EventReaderHepMC.cpp:434
dd4hep::sim::HepMC::EventHeader
HepMC EventHeader class used internally by the Geant4EventReaderHepMC plugin.
Definition: Geant4EventReaderHepMC.cpp:119
dd4hep::sim::Geant4EventReaderHepMC::skipEvent
virtual EventReaderStatus skipEvent() override
Skip event. To be implemented for sequential sources.
Definition: Geant4EventReaderHepMC.cpp:72
dd4hep::sim::G4PARTICLE_GEN_OTHER
@ G4PARTICLE_GEN_OTHER
Definition: Geant4Particle.h:76
Printout.h
dd4hep::sim::HepMC::read_vertex
int read_vertex(EventStream &info, std::istream &is, std::istringstream &iline)
Definition: Geant4EventReaderHepMC.cpp:493
dd4hep::sim::HepMC::read_units
int read_units(EventStream &info, std::istringstream &input)
Definition: Geant4EventReaderHepMC.cpp:633
dd4hep::sim::HepMC::EventStream::read
bool read()
Definition: Geant4EventReaderHepMC.cpp:728