Class Reference for E1039 Core & Analysis Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
HepMCNodeReader.cc
Go to the documentation of this file.
1 #include "HepMCNodeReader.h"
2 #include "PHG4InEvent.h"
3 #include "PHG4Particlev1.h"
4 #include "PHG4Particlev2.h"
5 
6 #include <Geant4/G4ParticleDefinition.hh>
7 #include <Geant4/G4ParticleTable.hh>
8 
10 
13 
14 #include <phool/PHRandomSeed.h>
15 #include <phool/getClass.h>
16 #include <phool/recoConsts.h>
17 
18 #include <HepMC/GenEvent.h>
19 
20 #include <gsl/gsl_const.h>
21 #include <gsl/gsl_randist.h>
22 #include <gsl/gsl_rng.h>
23 
24 #include <cassert>
25 #include <list>
26 
27 using namespace std;
28 
31 //
33 //const double mm_over_c_to_sec = 0.1 / GSL_CONST_CGS_SPEED_OF_LIGHT;
35 //const double mm_over_c_to_nanosecond = mm_over_c_to_sec * 1e9;
37 
39 class IsStateFinal
40 {
41  public:
43  bool operator()(const HepMC::GenParticle *p)
44  {
45  if (!p->end_vertex() && p->status() == 1) return 1;
46  return 0;
47  }
48 };
49 
50 static IsStateFinal isfinal;
51 
52 HepMCNodeReader::HepMCNodeReader(const std::string &name)
53  : SubsysReco(name)
54  , use_seed(0)
55  , seed(0)
56  , vertex_pos_x(0.0)
57  , vertex_pos_y(0.0)
58  , vertex_pos_z(0.0)
59  , vertex_t0(0.0)
60  , width_vx(0.0)
61  , width_vy(0.0)
62  , width_vz(0.0)
63  , _particle_filter_on(false)
64 {
65  RandomGenerator = gsl_rng_alloc(gsl_rng_mt19937);
66  _particle_filter_pid.clear();
67  return;
68 }
69 
71 {
72  gsl_rng_free(RandomGenerator);
73 }
74 
76 {
77  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
78  if (!ineve)
79  {
80  PHNodeIterator iter(topNode);
81  PHCompositeNode *dstNode;
82  dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
83 
84  ineve = new PHG4InEvent();
85  PHDataNode<PHObject> *newNode =
86  new PHDataNode<PHObject>(ineve, "PHG4INEVENT", "PHObject");
87  dstNode->addNode(newNode);
88  }
89  unsigned int phseed = PHRandomSeed(); // fixed seed is handled in this funtcion, need to call it to preserve random numbder order even if we override it
90  if (use_seed)
91  {
92  phseed = seed;
93  cout << Name() << " override random seed: " << phseed << endl;
94  }
95  gsl_rng_set(RandomGenerator, phseed);
96  return 0;
97 }
98 
100 {
101  // For pile-up simulation: define GenEventMap
102  PHHepMCGenEventMap *genevtmap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
103 
104  if (!genevtmap)
105  {
106  static bool once = true;
107 
108  if (once and Verbosity())
109  {
110  once = false;
111 
112  cout << "HepMCNodeReader::process_event - No PHHepMCGenEventMap node. Do not perform HepMC->Geant4 input" << endl;
113  }
114 
116  }
117 
118  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
119  if (!ineve)
120  {
121  cout << PHWHERE << "no PHG4INEVENT node" << endl;
123  }
124 
126 
127  float worldsizex = rc->get_FloatFlag("WorldSizex");
128  float worldsizey = rc->get_FloatFlag("WorldSizey");
129  float worldsizez = rc->get_FloatFlag("WorldSizez");
130  string worldshape = rc->get_CharFlag("WorldShape");
131 
132  enum
133  {
134  ShapeG4Tubs = 0,
135  ShapeG4Box = 1
136  };
137 
138  int ishape;
139  if (worldshape == "G4Tubs")
140  {
141  ishape = ShapeG4Tubs;
142  }
143  else if (worldshape == "G4Box")
144  {
145  ishape = ShapeG4Box;
146  }
147  else
148  {
149  cout << PHWHERE << " unknown world shape " << worldshape << endl;
150  exit(1);
151  }
152 
153  // For pile-up simulation: loop over PHHepMC event map
154  // insert highest embedding ID event first, whose vertex maybe resued in PHG4ParticleGeneratorBase::ReuseExistingVertex()
155  int vtxindex = -1;
156  for (PHHepMCGenEventMap::ReverseIter iter = genevtmap->rbegin(); iter != genevtmap->rend(); ++iter)
157  {
158  PHHepMCGenEvent *genevt = iter->second;
159  assert(genevt);
160 
161  if (genevt->is_simulated())
162  {
163  if (Verbosity())
164  {
165  cout << "HepMCNodeReader::process_event - this event is already simulated. Move on: ";
166  genevt->identify();
167  }
168 
169  continue;
170  }
171 
172  HepMC::GenEvent *evt = genevt->getEvent();
173  if (!evt)
174  {
175  cout << PHWHERE << " no evt pointer under HEPMC Node found:";
176  genevt->identify();
178  }
179 
180  genevt->is_simulated(true);
181 
182  const int embed_flag = genevt->get_embedding_id();
183 
184  double xshift = vertex_pos_x + genevt->get_collision_vertex().x();
185  double yshift = vertex_pos_y + genevt->get_collision_vertex().y();
186  double zshift = vertex_pos_z + genevt->get_collision_vertex().z();
187  double tshift = vertex_t0 + genevt->get_collision_vertex().t();
188 
189  if (width_vx > 0.0)
190  xshift += smeargauss(width_vx);
191  else if (width_vx < 0.0)
192  xshift += smearflat(width_vx);
193 
194  if (width_vy > 0.0)
195  yshift += smeargauss(width_vy);
196  else if (width_vy < 0.0)
197  yshift += smearflat(width_vy);
198 
199  if (width_vz > 0.0)
200  zshift += smeargauss(width_vz);
201  else if (width_vz < 0.0)
202  zshift += smearflat(width_vz);
203 
204  std::list<HepMC::GenParticle *> finalstateparticles;
205  std::list<HepMC::GenParticle *>::const_iterator fiter;
206 
207  // units in G4 interface are GeV and CM as in PHENIX convention
208  const double mom_factor = HepMC::Units::conversion_factor(evt->momentum_unit(), HepMC::Units::GEV);
209  const double length_factor = HepMC::Units::conversion_factor(evt->length_unit(), HepMC::Units::CM);
210  const double time_factor = HepMC::Units::conversion_factor(evt->length_unit(), HepMC::Units::CM) / GSL_CONST_CGS_SPEED_OF_LIGHT * 1e9; // from length_unit()/c to ns
211 
212  for (HepMC::GenEvent::vertex_iterator v = evt->vertices_begin();
213  v != evt->vertices_end();
214  ++v)
215  {
216  finalstateparticles.clear();
217  for (HepMC::GenVertex::particle_iterator p =
218  (*v)->particles_begin(HepMC::children);
219  p != (*v)->particles_end(HepMC::children); ++p)
220  {
221  if (isfinal(*p))
222  {
223  if(_particle_filter_on) {
224  if(!PassParticleFilter(*p)) continue;
225  }
226  finalstateparticles.push_back(*p);
227  }
228  }
229 
230  if (!finalstateparticles.empty())
231  {
232  double xpos = (*v)->position().x() * length_factor + xshift;
233  double ypos = (*v)->position().y() * length_factor + yshift;
234  double zpos = (*v)->position().z() * length_factor + zshift;
235  double time = (*v)->position().t() * time_factor + tshift;
236 
237  if (Verbosity() > 1)
238  {
239  cout << "Vertex : " << endl;
240  (*v)->print();
241  cout << "id: " << (*v)->barcode() << endl;
242  cout << "x: " << xpos << endl;
243  cout << "y: " << ypos << endl;
244  cout << "z: " << zpos << endl;
245  cout << "t: " << time << endl;
246  cout << "Particles" << endl;
247  }
248 
249  if (ishape == ShapeG4Tubs)
250  {
251  if (sqrt(xpos * xpos + ypos * ypos) > worldsizey / 2 ||
252  fabs(zpos) > worldsizez / 2)
253  {
254  cout << "vertex x/y/z" << xpos << "/" << ypos << "/" << zpos
255  << " outside world volume radius/z (+-) " << worldsizex / 2
256  << "/" << worldsizez / 2 << ", dropping it and its particles"
257  << endl;
258  continue;
259  }
260  }
261  else if (ishape == ShapeG4Box)
262  {
263  if (fabs(xpos) > worldsizex / 2 || fabs(ypos) > worldsizey / 2 ||
264  fabs(zpos) > worldsizez / 2)
265  {
266  cout << "Vertex x/y/z " << xpos << "/" << ypos << "/" << zpos
267  << " outside world volume x/y/z (+-) " << worldsizex / 2 << "/"
268  << worldsizey / 2 << "/" << worldsizez / 2
269  << ", dropping it and its particles" << endl;
270  continue;
271  }
272  }
273  else
274  {
275  cout << PHWHERE << " shape " << ishape << " not implemented. exiting"
276  << endl;
277  exit(1);
278  }
279 
280  // For pile-up simulation: vertex position
281  vtxindex = ineve->AddVtx(xpos, ypos, zpos, time);
282  int trackid = -1;
283  for (fiter = finalstateparticles.begin();
284  fiter != finalstateparticles.end();
285  ++fiter)
286  {
287  ++trackid;
288 
289  if (verbosity > 1) (*fiter)->print();
290 
291  PHG4Particle *particle = new PHG4Particlev1();
292  particle->set_pid((*fiter)->pdg_id());
293  particle->set_px((*fiter)->momentum().px() * mom_factor);
294  particle->set_py((*fiter)->momentum().py() * mom_factor);
295  particle->set_pz((*fiter)->momentum().pz() * mom_factor);
296  particle->set_barcode((*fiter)->barcode());
297 
298  ineve->AddParticle(vtxindex, particle);
299 
300  if (embed_flag != 0) ineve->AddEmbeddedParticle(particle, embed_flag);
301  }
302  } // if (!finalstateparticles.empty())
303 
304  } // for (HepMC::GenEvent::vertex_iterator v = evt->vertices_begin();
305 
306  } // For pile-up simulation: loop end for PHHepMC event map
307 
308  if (verbosity > 0) ineve->identify();
309 
311 }
312 
313 double HepMCNodeReader::smeargauss(const double width)
314 {
315  if (width == 0) return 0;
316  return gsl_ran_gaussian(RandomGenerator, width);
317 }
318 
319 bool HepMCNodeReader::PassParticleFilter(HepMC::GenParticle* p) {
320  for(int pid : _particle_filter_pid) {
321  if(p->pdg_id() == pid) return true;
322  }
323  return false;
324 }
325 
326 double HepMCNodeReader::smearflat(const double width)
327 {
328  if (width == 0) return 0;
329  return 2.0 * width * (gsl_rng_uniform_pos(RandomGenerator) - 0.5);
330 }
331 
332 void HepMCNodeReader::VertexPosition(const double v_x, const double v_y,
333  const double v_z)
334 {
335  cout << "HepMCNodeReader::VertexPosition - WARNING - this function is depreciated. "
336  << "HepMCNodeReader::VertexPosition() move all HEPMC subevents to a new vertex location. "
337  << "This also leads to a different vertex is used for HepMC subevent in Geant4 than that recorded in the HepMCEvent Node."
338  << "Recommendation: the vertex shifts are better controlled for individually HEPMC subevents in Fun4AllHepMCInputManagers and event generators."
339  << endl;
340 
341  vertex_pos_x = v_x;
342  vertex_pos_y = v_y;
343  vertex_pos_z = v_z;
344  return;
345 }
346 
347 void HepMCNodeReader::SmearVertex(const double s_x, const double s_y,
348  const double s_z)
349 {
350  cout << "HepMCNodeReader::SmearVertex - WARNING - this function is depreciated. "
351  << "HepMCNodeReader::SmearVertex() smear each HEPMC subevents to a new vertex location. "
352  << "This also leads to a different vertex is used for HepMC subevent in Geant4 than that recorded in the HepMCEvent Node."
353  << "Recommendation: the vertex smears are better controlled for individually HEPMC subevents in Fun4AllHepMCInputManagers and event generators."
354  << endl;
355 
356  width_vx = s_x;
357  width_vy = s_y;
358  width_vz = s_z;
359  return;
360 }
361 
362 void HepMCNodeReader::Embed(const int i)
363 {
364  cout << "HepMCNodeReader::Embed - WARNING - this function is depreciated. "
365  << "Embedding IDs are controlled for individually HEPMC subevents in Fun4AllHepMCInputManagers and event generators."
366  << endl;
367 
368  return;
369 }
int get_embedding_id() const
PHHepMCGenEventMap is collection of HEPMC events input into this simulation map of embedding ID -&gt; PH...
int verbosity
The verbosity level. 0 means not verbose at all.
Definition: Fun4AllBase.h:75
int Init(PHCompositeNode *topNode)
PHNode * findFirst(const std::string &, const std::string &)
HepMCNodeReader(const std::string &name="HEPMCREADER")
virtual ~HepMCNodeReader()
#define PHWHERE
Definition: phool.h:23
virtual void set_pz(const double x)
Definition: PHG4Particle.h:36
void AddEmbeddedParticle(PHG4Particle *particle, int flag)
Definition: PHG4InEvent.h:25
ConstReverseIter rbegin() const
iterator from lowest ID to highest, i.e. signal to background
PHBoolean addNode(PHNode *)
int AddParticle(const int vtxid, PHG4Particle *particle)
Definition: PHG4InEvent.cc:63
virtual HepMC::GenEvent * getEvent()
virtual void set_px(const double x)
Definition: PHG4Particle.h:34
void VertexPosition(const double v_x, const double v_y, const double v_z)
virtual void set_py(const double x)
Definition: PHG4Particle.h:35
std::map< int, PHHepMCGenEvent * >::reverse_iterator ReverseIter
static recoConsts * instance()
Definition: recoConsts.cc:7
virtual float get_FloatFlag(const std::string &name) const
Definition: PHFlag.cc:83
ConstReverseIter rend() const
const HepMC::FourVector & get_collision_vertex() const
collision vertex position in the Hall coordinate system, use PHENIX units of cm, ns ...
void SmearVertex(const double s_x, const double s_y, const double s_z)
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
bool is_simulated() const
whether this event has been processed in Geant4 simulation
bool PassParticleFilter(HepMC::GenParticle *p)
virtual int Verbosity() const
Gets the verbosity of this module.
Definition: Fun4AllBase.h:64
virtual void set_barcode(const int barcode)
Definition: PHG4Particle.h:39
virtual const std::string get_CharFlag(const std::string &flag) const
Definition: PHFlag.cc:13
static IsStateFinal isfinal
virtual void set_pid(const int i)
Definition: PHG4Particle.h:33
int AddVtx(const double x, const double y, const double z, const double t)
Definition: PHG4InEvent.cc:38
void Embed(const int i=1)
int process_event(PHCompositeNode *topNode)
virtual void identify(std::ostream &os=std::cout) const
Definition: PHG4InEvent.cc:134
virtual void identify(std::ostream &os=std::cout) const