DD4hep  1.30.0
Detector Description Toolkit for High Energy Physics
Geant4EventReaderHepEvt.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 
25 // Framework include files
26 #include <DDG4/Geant4InputAction.h>
27 
28 // C/C++ include files
29 #include <fstream>
30 
32 namespace dd4hep {
33 
35  namespace sim {
36 
38 
48 
49  protected:
50  std::ifstream m_input;
51  int m_format;
52 
53  public:
55  explicit Geant4EventReaderHepEvt(const std::string& nam, int format);
59  virtual EventReaderStatus readParticles(int event_number,
60  Vertices& vertices,
61  std::vector<Particle*>& particles);
62  virtual EventReaderStatus moveToEvent(int event_number);
64  };
65  } /* End namespace sim */
66 } /* End namespace dd4hep */
67 
68 //====================================================================
69 // AIDA Detector description implementation
70 //--------------------------------------------------------------------
71 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
72 // All rights reserved.
73 //
74 // For the licensing terms see $DD4hepINSTALL/LICENSE.
75 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
76 //
77 // Author : M.Frank
78 //
79 //====================================================================
80 // #include "DDG4/Geant4EventReaderHepEvt"
81 
82 // Framework include files
83 #include <DDG4/Factories.h>
84 #include <DD4hep/Printout.h>
85 #include <CLHEP/Units/SystemOfUnits.h>
86 #include <CLHEP/Units/PhysicalConstants.h>
87 
88 // C/C++ include files
89 #include <cerrno>
90 
91 using namespace dd4hep::sim;
92 using PropertyMask = dd4hep::detail::ReferenceBitMask<int>;
93 
94 #define HEPEvtShort 1
95 #define HEPEvtLong 2
96 
97 // Local declarations in anaonymous namespace
98 namespace {
99  class Geant4EventReaderHepEvtShort : public Geant4EventReaderHepEvt {
100  public:
102  explicit Geant4EventReaderHepEvtShort(const std::string& nam) : Geant4EventReaderHepEvt(nam,HEPEvtShort) {}
104  virtual ~Geant4EventReaderHepEvtShort() {}
105  };
106  class Geant4EventReaderHepEvtLong : public Geant4EventReaderHepEvt {
107  public:
109  explicit Geant4EventReaderHepEvtLong(const std::string& nam) : Geant4EventReaderHepEvt(nam,HEPEvtLong) {}
111  virtual ~Geant4EventReaderHepEvtLong() {}
112  };
113 }
114 
115 // Factory entry
116 DECLARE_GEANT4_EVENT_READER(Geant4EventReaderHepEvtShort)
117 // Factory entry
118 DECLARE_GEANT4_EVENT_READER(Geant4EventReaderHepEvtLong)
119 
120 
121 Geant4EventReaderHepEvt::Geant4EventReaderHepEvt(const std::string& nam, int format)
123 : Geant4EventReader(nam), m_input(), m_format(format)
124 {
125  // Now open the input file:
126  m_input.open(nam.c_str(), std::ifstream::in);
127  if ( !m_input.good() ) {
128  except("Geant4EventReaderHepEvt","+++ Failed to open input stream: %s Error:%s",
129  nam.c_str(), ::strerror(errno));
130  }
131 }
132 
134 Geant4EventReaderHepEvt::~Geant4EventReaderHepEvt() {
135  m_input.close();
136 }
137 
140 Geant4EventReaderHepEvt::moveToEvent(int event_number) {
141  if( m_currEvent == 0 && event_number != 0 ) {
142  printout(INFO,"EventReaderHepEvt::moveToEvent","Skipping the first %d events ", event_number );
143  printout(INFO,"EventReaderHepEvt::moveToEvent","Event number before skipping: %d", m_currEvent );
144  while ( m_currEvent < event_number ) {
145  std::vector<Particle*> particles;
146  Vertices vertices ;
147  EventReaderStatus sc = readParticles(m_currEvent,vertices,particles);
148  for_each(vertices.begin(), vertices.end(), detail::deleteObject<Vertex>);
149  for_each(particles.begin(), particles.end(), detail::deleteObject<Particle>);
150  if ( sc != EVENT_READER_OK ) return sc;
151  //Current event is increased in readParticles already!
152  // ++m_currEvent;
153  }
154  }
155  printout(INFO,"EventReaderHepEvt::moveToEvent","Event number after skipping: %d", m_currEvent );
156  return EVENT_READER_OK;
157 }
158 
161 Geant4EventReaderHepEvt::readParticles(int /* event_number */,
162  Vertices& vertices,
163  std::vector<Particle*>& particles) {
164 
165 
166  // First check the input file status
167  if ( m_input.eof() ) {
168  return EVENT_READER_EOF;
169  }
170  else if ( !m_input.good() ) {
171  return EVENT_READER_IO_ERROR;
172  }
173  //static const double c_light = 299.792;// mm/ns
174  //
175  // Read the event, check for errors
176  //
177  unsigned NHEP(0); // number of entries
178  m_input >> NHEP;
179 
180  if( m_input.eof() ) { return EVENT_READER_EOF; }
181 
182 
183  //check loop variable read from input file and chack that is reasonable
184  // should fix coverity issue: "Using tainted variable NHEP as a loop boundary."
185 
186  if( NHEP > 1e6 ){
187  printout(ERROR,"EventReaderHepEvt::readParticles","Cannot read in more than million particles, but %d requested", NHEP );
188  return EVENT_READER_EOF;
189  }
190 
191  //
192  // Loop over particles
193  int ISTHEP(0); // status code
194  int IDHEP(0); // PDG code
195  int JMOHEP1(0); // first mother
196  int JMOHEP2(0); // last mother
197  int JDAHEP1(0); // first daughter
198  int JDAHEP2(0); // last daughter
199  double PHEP1(0); // px in GeV/c
200  double PHEP2(0); // py in GeV/c
201  double PHEP3(0); // pz in GeV/c
202  double PHEP4(0); // energy in GeV
203  double PHEP5(0); // mass in GeV/c**2
204  double VHEP1(0); // x vertex position in mm
205  double VHEP2(0); // y vertex position in mm
206  double VHEP3(0); // z vertex position in mm
207  double VHEP4(0); // production time in mm/c
208 
209  std::vector<int> daughter1;
210  std::vector<int> daughter2;
211 
212  for( unsigned IHEP=0; IHEP<NHEP; IHEP++ ) {
213  if ( m_format == HEPEvtShort )
214  m_input >> ISTHEP >> IDHEP >> JDAHEP1 >> JDAHEP2
215  >> PHEP1 >> PHEP2 >> PHEP3 >> PHEP5;
216  else
217  m_input >> ISTHEP >> IDHEP
218  >> JMOHEP1 >> JMOHEP2
219  >> JDAHEP1 >> JDAHEP2
220  >> PHEP1 >> PHEP2 >> PHEP3
221  >> PHEP4 >> PHEP5
222  >> VHEP1 >> VHEP2 >> VHEP3
223  >> VHEP4;
224 
225  if(m_input.eof())
226  return EVENT_READER_EOF;
227 
228  if(! m_input.good())
229  return EVENT_READER_IO_ERROR;
230 
231  //
232  // create a MCParticle and fill it from stdhep info
233  Particle* p = new Particle(IHEP);
234  PropertyMask status(p->status);
235  //
236  // PDGID
237  p->pdgID = IDHEP;
238  p->charge = 0;
239  //
240  // Momentum vector
241  p->pex = p->psx = PHEP1*CLHEP::GeV;
242  p->pey = p->psy = PHEP2*CLHEP::GeV;
243  p->pez = p->psz = PHEP3*CLHEP::GeV;
244  //
245  // Mass
246  p->mass = PHEP5*CLHEP::GeV;
247  //
248  // Vertex
249  p->vsx = VHEP1*CLHEP::mm;
250  p->vsy = VHEP2*CLHEP::mm;
251  p->vsz = VHEP3*CLHEP::mm;
252  // endpoint (missing information in HEPEvt files)
253  p->vex = 0.0;
254  p->vey = 0.0;
255  p->vez = 0.0;
256  //
257  // Creation time (note the units [1/c_light])
258  p->time = VHEP4 * CLHEP::mm / CLHEP::c_light;
259  p->properTime = VHEP4 * CLHEP::mm / CLHEP::c_light;
260  //
261  // Generator status
262  // Simulator status 0 until simulator acts on it
263  p->status = 0;
264  if ( ISTHEP == 0 ) status.set(G4PARTICLE_GEN_EMPTY);
265  else if ( ISTHEP == 1 ) status.set(G4PARTICLE_GEN_STABLE);
266  else if ( ISTHEP == 2 ) status.set(G4PARTICLE_GEN_DECAYED);
267  else if ( ISTHEP == 3 ) status.set(G4PARTICLE_GEN_DOCUMENTATION);
268  else status.set(G4PARTICLE_GEN_DOCUMENTATION);
269  // RAW Generator status
270  p->genStatus = ISTHEP&G4PARTICLE_GEN_STATUS_MASK;
271 
272  //
273  // Keep daughters information for later
274  daughter1.emplace_back(JDAHEP1);
275  daughter2.emplace_back(JDAHEP2);
276  //
277  // Add the particle to the collection vector
278  particles.emplace_back(p);
279  //
280  }// End loop over particles
281 
282  //
283  // Now make a second loop over the particles, checking the daughter
284  // information. This is not always consistent with parent
285  // information, and this utility assumes all parents listed are
286  // parents and all daughters listed are daughters
287  //
288  for( unsigned IHEP=0; IHEP<NHEP; IHEP++ ) {
289  struct ParticleHandler {
290  Particle* m_part;
291  ParticleHandler(Particle* p) : m_part(p) {}
292  void addParent(const Particle* p) {
293  m_part->parents.insert(p->id);
294  }
295  void addDaughter(const Particle* p) {
296  if(m_part->daughters.find(p->id) == m_part->daughters.end()) {
297  m_part->daughters.insert(p->id);
298  }
299  }
300  Particle* findParent(const Particle* p) {
301  return m_part->parents.find(p->id)==m_part->parents.end() ? 0 : m_part;
302  }
303  };
304 
305  //
306  // Get the MCParticle
307  //
308  Particle* mcp = particles[IHEP];
309  ParticleHandler theParticle(mcp);
310  //
311  // Get the daughter information, discarding extra information
312  // sometimes stored in daughter variables.
313  //
314  int fd = daughter1[IHEP] - 1;
315  int ld = daughter2[IHEP] - 1;
316 
317  //
318  // As with the parents, look for range, 2 discreet or 1 discreet
319  // daughter.
320  if( (fd > -1) && (fd < int(particles.size())) && (ld > -1) && (ld < int(particles.size())) ) {
321  if(ld >= fd) {
322  for(int id=fd;id<ld+1;id++) {
323  //
324  // Get the daughter, and see if it already lists this particle as
325  // a parent. If not, add this particle as a parent
326  //
327  ParticleHandler part(particles[id]);
328  if ( !part.findParent(mcp) ) part.addParent(mcp);
329  theParticle.addDaughter(particles[id]);
330  }
331  }
332  else { // Same logic, discrete cases
333  ParticleHandler part_fd(particles[fd]);
334  if ( !part_fd.findParent(mcp) ) part_fd.addParent(mcp);
335  theParticle.addDaughter(particles[fd]);
336 
337  ParticleHandler part_ld(particles[ld]);
338  if ( !part_ld.findParent(mcp) ) part_ld.addParent(mcp);
339  theParticle.addDaughter(particles[ld]);
340  }
341  }
342  else if(fd > -1 && fd < int(particles.size())) {
343  ParticleHandler part(particles[fd]);
344  if ( !part.findParent(mcp) ) part.addParent(mcp);
345  }
346  else if(ld > -1 && ld < int(particles.size())) {
347  ParticleHandler part(particles[ld]);
348  if ( !part.findParent(mcp) ) part.addParent(mcp);
349  }
350  } // End second loop over particles
351 
352 
353  //fg: we simply add all particles without parents as with their own vertex.
354  // This might include the incoming beam particles, e.g. in
355  // the case of lcio files written with Whizard2, which is slightly odd,
356  // however should be treated correctly in Geant4InputHandling.cpp.
357  // We no longer make an attempt to identify the incoming particles
358  // based on the generator status, as this varies widely with different
359  // generators.
360 
361  for( std::size_t i=0; i < particles.size(); ++i ) {
362  Geant4ParticleHandle p(particles[i]);
363  if ( p->parents.size() == 0 ) {
364  Geant4Vertex* vtx = new Geant4Vertex ;
365  vertices.emplace_back( vtx );
366  vtx->x = p->vsx;
367  vtx->y = p->vsy;
368  vtx->z = p->vsz;
369  vtx->time = p->time;
370 
371  vtx->out.insert(p->id) ;
372  }
373  }
374 
375  ++m_currEvent;
376  return EVENT_READER_OK;
377 }
378 
HEPEvtLong
#define HEPEvtLong
Definition: Geant4EventReaderHepEvt.cpp:95
dd4hep::sim::Geant4EventReaderHepEvt::readParticles
virtual EventReaderStatus readParticles(int event_number, Vertices &vertices, std::vector< Particle * > &particles)
Read an event and fill a vector of MCParticles.
dd4hep::sim::Geant4EventReader::EVENT_READER_EOF
@ EVENT_READER_EOF
Definition: Geant4InputAction.h:75
dd4hep::sim::G4PARTICLE_GEN_STABLE
@ G4PARTICLE_GEN_STABLE
Definition: Geant4Particle.h:71
dd4hep::sim::Geant4Vertex::y
double y
Definition: Geant4Vertex.h:53
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
DECLARE_GEANT4_EVENT_READER
#define DECLARE_GEANT4_EVENT_READER(name)
Plugin defintion to create event reader objects.
Definition: Factories.h:237
dd4hep::sim::Geant4EventReaderHepEvt::~Geant4EventReaderHepEvt
virtual ~Geant4EventReaderHepEvt()
Default destructor.
HEPEvtShort
#define HEPEvtShort
Definition: Geant4EventReaderHepEvt.cpp:94
dd4hep::sim::G4PARTICLE_GEN_STATUS_MASK
@ G4PARTICLE_GEN_STATUS_MASK
Definition: Geant4Particle.h:83
dd4hep::sim::Geant4EventReaderHepEvt::skipEvent
virtual EventReaderStatus skipEvent()
Skip event. To be implemented for sequential sources.
Definition: Geant4EventReaderHepEvt.cpp:63
dd4hep::sim::Geant4EventReaderHepEvt::m_input
std::ifstream m_input
Definition: Geant4EventReaderHepEvt.cpp:50
dd4hep::sim::G4PARTICLE_GEN_DOCUMENTATION
@ G4PARTICLE_GEN_DOCUMENTATION
Definition: Geant4Particle.h:73
dd4hep::sim::Geant4EventReader::m_currEvent
int m_currEvent
Current event number.
Definition: Geant4InputAction.h:83
dd4hep::sim::Geant4EventReaderHepEvt::m_format
int m_format
Definition: Geant4EventReaderHepEvt.cpp:51
dd4hep::sim::G4PARTICLE_GEN_EMPTY
@ G4PARTICLE_GEN_EMPTY
Definition: Geant4Particle.h:70
dd4hep::sim::Geant4EventReader::Particle
Geant4Particle Particle
Definition: Geant4InputAction.h:64
dd4hep::sim::Geant4Vertex::x
double x
Vertex data.
Definition: Geant4Vertex.h:53
dd4hep::sim::G4PARTICLE_GEN_DECAYED
@ G4PARTICLE_GEN_DECAYED
Definition: Geant4Particle.h:72
Geant4EventReaderHepEvt
Plugin to read HepEvt ASCII files.
dd4hep::sim::Geant4EventReaderHepEvt::Geant4EventReaderHepEvt
Geant4EventReaderHepEvt(const std::string &nam, int format)
Initializing constructor.
dd4hep::sim::Geant4EventReader
Basic geant4 event reader class. This interface/base-class must be implemented by concrete readers.
Definition: Geant4InputAction.h:60
dd4hep::sim::Geant4Vertex::out
Particles out
The list of outgoing particles.
Definition: Geant4Vertex.h:55
Factories.h
dd4hep::sim::Geant4EventReader::EVENT_READER_IO_ERROR
@ EVENT_READER_IO_ERROR
Definition: Geant4InputAction.h:74
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::Geant4EventReader::EVENT_READER_OK
@ EVENT_READER_OK
Definition: Geant4InputAction.h:70
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::Geant4Vertex::time
double time
Definition: Geant4Vertex.h:53
dd4hep::sim::Geant4Particle::parents
Particles parents
The list of parents of this MC particle.
Definition: Geant4Particle.h:138
dd4hep::sim::Geant4Vertex::z
double z
Definition: Geant4Vertex.h:53
dd4hep::sim::Geant4Vertex
Data structure to store the MC vertex information.
Definition: Geant4Vertex.h:45
dd4hep::sim::Geant4EventReader::Vertices
std::vector< Vertex * > Vertices
Definition: Geant4InputAction.h:66
Printout.h
dd4hep::sim::Geant4EventReaderHepEvt::moveToEvent
virtual EventReaderStatus moveToEvent(int event_number)
Move to the indicated event number.