Class Reference for E1039 Core & Analysis Software
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 <TRandom3.h>
25 
26 #include <cassert>
27 #include <list>
28 
29 using namespace std;
30 
33 //
35 //const double mm_over_c_to_sec = 0.1 / GSL_CONST_CGS_SPEED_OF_LIGHT;
37 //const double mm_over_c_to_nanosecond = mm_over_c_to_sec * 1e9;
39 
41 class IsStateFinal
42 {
43 public:
45  bool operator()(const HepMC::GenParticle *p)
46  {
47  if (!p->end_vertex() && p->status() == 1) return 1;
48  return 0;
49  }
50 };
51 
52 static IsStateFinal isfinal;
53 
54 HepMCNodeReader::HepMCNodeReader(const std::string &name)
55  : SubsysReco(name)
56  , use_seed(0)
57  , seed(0)
58  , vertex_pos_x(0.0)
59  , vertex_pos_y(0.0)
60  , vertex_pos_z(0.0)
61  , vertex_t0(0.0)
62  , width_vx(0.0)
63  , width_vy(0.0)
64  , width_vz(0.0)
65  , _bkg_mode(false)//Abi
66  , _pxy2pz_rat(0.25)
67  , _particle_filter_on(false)
68  , _position_filter_on(false)
69  , _pos_filter_x_min(0.0)
70  , _pos_filter_x_max(0.0)
71  , _pos_filter_y_min(0.0)
72  , _pos_filter_y_max(0.0)
73  , _pos_filter_z_min(0.0)
74  , _pos_filter_z_max(0.0)
75 {
76  RandomGenerator = gsl_rng_alloc(gsl_rng_mt19937);
77  _particle_filter_pid.clear();
78  return;
79 }
80 
82 {
83  gsl_rng_free(RandomGenerator);
84 }
85 
87 {
88  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
89  if (!ineve)
90  {
91  PHNodeIterator iter(topNode);
92  PHCompositeNode *dstNode;
93  dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
94 
95  ineve = new PHG4InEvent();
96  PHIODataNode<PHObject> *newNode =
97  new PHIODataNode<PHObject>(ineve, "PHG4INEVENT", "PHObject");
98  dstNode->addNode(newNode);
99  }
100  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
101  if (use_seed)
102  {
103  phseed = seed;
104  cout << Name() << " override random seed: " << phseed << endl;
105  }
106  gsl_rng_set(RandomGenerator, phseed);
107  return 0;
108 }
109 
111 {
112  // For pile-up simulation: define GenEventMap
113  PHHepMCGenEventMap *genevtmap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
114 
115  if (!genevtmap)
116  {
117  static bool once = true;
118 
119  if (once and Verbosity())
120  {
121  once = false;
122 
123  cout << "HepMCNodeReader::process_event - No PHHepMCGenEventMap node. Do not perform HepMC->Geant4 input" << endl;
124  }
125 
127  }
128 
129  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
130  if (!ineve)
131  {
132  cout << PHWHERE << "no PHG4INEVENT node" << endl;
134  }
135 
137 
138  float worldsizex = rc->get_FloatFlag("WorldSizex");
139  float worldsizey = rc->get_FloatFlag("WorldSizey");
140  float worldsizez = rc->get_FloatFlag("WorldSizez");
141  string worldshape = rc->get_CharFlag("WorldShape");
142 
143  enum
144  {
145  ShapeG4Tubs = 0,
146  ShapeG4Box = 1
147  };
148 
149  int ishape;
150  if (worldshape == "G4Tubs")
151  {
152  ishape = ShapeG4Tubs;
153  }
154  else if (worldshape == "G4Box")
155  {
156  ishape = ShapeG4Box;
157  }
158  else
159  {
160  cout << PHWHERE << " unknown world shape " << worldshape << endl;
161  exit(1);
162  }
163 
164  // For pile-up simulation: loop over PHHepMC event map
165  // insert highest embedding ID event first, whose vertex maybe resued in PHG4ParticleGeneratorBase::ReuseExistingVertex()
166  int vtxindex = -1;
167  for (PHHepMCGenEventMap::ReverseIter iter = genevtmap->rbegin(); iter != genevtmap->rend(); ++iter)
168  {
169  PHHepMCGenEvent *genevt = iter->second;
170  assert(genevt);
171 
172  if (genevt->is_simulated())
173  {
174  if (Verbosity())
175  {
176  cout << "HepMCNodeReader::process_event - this event is already simulated. Move on: ";
177  genevt->identify();
178  }
179 
180  continue;
181  }
182 
183  HepMC::GenEvent *evt = genevt->getEvent();
184  if (!evt)
185  {
186  cout << PHWHERE << " no evt pointer under HEPMC Node found:";
187  genevt->identify();
189  }
190 
191  genevt->is_simulated(true);
192 
193  const int embed_flag = genevt->get_embedding_id();
194 
195  double xshift = vertex_pos_x + genevt->get_collision_vertex().x();
196  double yshift = vertex_pos_y + genevt->get_collision_vertex().y();
197  double zshift = vertex_pos_z + genevt->get_collision_vertex().z();
198  double tshift = vertex_t0 + genevt->get_collision_vertex().t();
199 
200  if (width_vx > 0.0)
201  xshift += smeargauss(width_vx);
202  else if (width_vx < 0.0)
203  xshift += smearflat(width_vx);
204 
205  if (width_vy > 0.0)
206  yshift += smeargauss(width_vy);
207  else if (width_vy < 0.0)
208  yshift += smearflat(width_vy);
209 
210  if (width_vz > 0.0)
211  zshift += smeargauss(width_vz);
212  else if (width_vz < 0.0)
213  zshift += smearflat(width_vz);
214 
215  std::list<HepMC::GenParticle *> finalstateparticles;
216  std::list<HepMC::GenParticle *>::const_iterator fiter;
217 
218  // units in G4 interface are GeV and CM as in PHENIX convention
219  const double mom_factor = HepMC::Units::conversion_factor(evt->momentum_unit(), HepMC::Units::GEV);
220  const double length_factor = HepMC::Units::conversion_factor(evt->length_unit(), HepMC::Units::CM);
221  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
222 
223  for (HepMC::GenEvent::vertex_iterator v = evt->vertices_begin();
224  v != evt->vertices_end();
225  ++v)
226  {
227  finalstateparticles.clear();
228  for (HepMC::GenVertex::particle_iterator p =
229  (*v)->particles_begin(HepMC::children);
230  p != (*v)->particles_end(HepMC::children); ++p)
231  {
232  if (isfinal(*p))
233  {
234  if(_particle_filter_on) {
235  if(!PassParticleFilter(*p)) continue;
236  }
237  finalstateparticles.push_back(*p);
238  }
239  }
240 
241  if (!finalstateparticles.empty())
242  {
243  double xpos = (*v)->position().x() * length_factor + xshift;
244  double ypos = (*v)->position().y() * length_factor + yshift;
245  double zpos = (*v)->position().z() * length_factor + zshift;
246  double time = (*v)->position().t() * time_factor + tshift;
247 
248  if (Verbosity() > 1)
249  {
250  cout << "Vertex : " << endl;
251  (*v)->print();
252  cout << "id: " << (*v)->barcode() << endl;
253  cout << "x: " << xpos << endl;
254  cout << "y: " << ypos << endl;
255  cout << "z: " << zpos << endl;
256  cout << "t: " << time << endl;
257  cout << "Particles" << endl;
258  }
259 
260  if (_position_filter_on &&
261  (xpos < _pos_filter_x_min || xpos > _pos_filter_x_max ||
262  ypos < _pos_filter_y_min || ypos > _pos_filter_y_max ||
263  zpos < _pos_filter_z_min || zpos > _pos_filter_z_max )) continue;
264 
265  if (ishape == ShapeG4Tubs)
266  {
267  if (sqrt(xpos * xpos + ypos * ypos) > worldsizey / 2 ||
268  fabs(zpos) > worldsizez / 2)
269  {
270  if (Verbosity()>0)
271  {
272  cout << "vertex x/y/z" << xpos << "/" << ypos << "/" << zpos
273  << " outside world volume radius/z (+-) " << worldsizex / 2
274  << "/" << worldsizez / 2 << ", dropping it and its particles"
275  << endl;
276  }
277  continue;
278  }
279  }
280  else if (ishape == ShapeG4Box)
281  {
282  if (fabs(xpos) > worldsizex / 2 || fabs(ypos) > worldsizey / 2 ||
283  fabs(zpos) > worldsizez / 2)
284  {
285  if (Verbosity()>0)
286  {
287  cout << "Vertex x/y/z " << xpos << "/" << ypos << "/" << zpos
288  << " outside world volume x/y/z (+-) " << worldsizex / 2 << "/"
289  << worldsizey / 2 << "/" << worldsizez / 2
290  << ", dropping it and its particles" << endl;
291  }
292  continue;
293  }
294  }
295  else
296  {
297  cout << PHWHERE << " shape " << ishape << " not implemented. exiting"
298  << endl;
299  exit(1);
300  }
301 
302  // For pile-up simulation: vertex position
303  vtxindex = ineve->AddVtx(xpos, ypos, zpos, time);
304  int trackid = -1;
305  for (fiter = finalstateparticles.begin();
306  fiter != finalstateparticles.end();
307  ++fiter)
308  {
309  ++trackid;
310 
311  if (verbosity > 1) (*fiter)->print();
312 
313  PHG4Particle *particle = new PHG4Particlev1();
314  particle->set_pid((*fiter)->pdg_id());
315  particle->set_px((*fiter)->momentum().px() * mom_factor);
316  particle->set_py((*fiter)->momentum().py() * mom_factor);
317  particle->set_pz((*fiter)->momentum().pz() * mom_factor);
318  particle->set_barcode((*fiter)->barcode());
319 
320 
321  if(_bkg_mode){//@ Filters for the background candidates
322 
324  double FMAG_hole_R = 2.5;//circular hole of diameter 5 cm and length 25 cm at FMAG
325  double FMAG_hole_L = 25.;
326 
327  double Transverse_pos = sqrt(pow(xpos,2)+pow(ypos,2));
328  double FMAG_X = 160.;
329  double FMAG_Y = 160.; //
330  double FMAG_Z = 503.;
331 
332  double pion_lambda_Fe = 20.42;
333  double pion_lambda_Cu = 18.51;
334  double pion_lambda_B = 48.75;
335  double pion_lambda_sh = 55.92;
336 
338  if(fabs((*fiter)->momentum().px()* mom_factor/((*fiter)->momentum().pz()* mom_factor))>_pxy2pz_rat) continue;
339  if(fabs((*fiter)->momentum().py()* mom_factor/((*fiter)->momentum().pz()* mom_factor))>_pxy2pz_rat) continue;
340  if((*fiter)->momentum().pz()<0) continue;
341 
343  if(fabs(xpos)>FMAG_X/2.) continue;
344  if(fabs(ypos)>FMAG_Y/2.) continue;
345 
346 
348  double coll_x = 7.8;// cm
349  double coll_y = 3.48;
350  double totMom =sqrt(pow((*fiter)->momentum().px() * mom_factor,2)+pow((*fiter)->momentum().py() * mom_factor,2)+pow((*fiter)->momentum().pz() * mom_factor,2));
351 
353 
354  if(genevt->get_collision_vertex().z()>-664.&&genevt->get_collision_vertex().z()<-542. && (fabs(xpos)>coll_x/2.||fabs(ypos)>coll_y/2.) )
355  {
356  if(zpos<-542.)
357  {
358  if (gRandom->Uniform(0,1) > exp (-(-genevt->get_collision_vertex().z()+zpos)/(pion_lambda_Cu))) continue;
359  }
360  else
361  {
362  if (gRandom->Uniform(0,1) > exp (-(-genevt->get_collision_vertex().z()-542.)/(pion_lambda_Cu))) continue;
363  }
364  }
365 
366 
368 
369  if(genevt->get_collision_vertex().z()>-171.5&& genevt->get_collision_vertex().z()<-0.001)
370  {
371  if(zpos<-0.001)
372  {
373  if (gRandom->Uniform(0,1) > exp (-(-genevt->get_collision_vertex().z()+zpos)/pion_lambda_sh)) continue;
374  }
375  else
376  {
377  if (gRandom->Uniform(0,1) > exp (-(171.5)/pion_lambda_sh)) continue;
378  }
379 
380  }
381 
382 
383  if(genevt->get_collision_vertex().z()<-171.5 && Transverse_pos<7.62 )
384  {
385  if(zpos>-171.5 && zpos<-0.001)
386  {
387  if (gRandom->Uniform(0,1) > exp (-(171.5+zpos)/pion_lambda_sh)) continue;
388  }
389  if(zpos>-.001)
390  {
391  if (gRandom->Uniform(0,1) > exp (-(171.5)/pion_lambda_sh)) continue;
392  }
393 
394  }
395 
396 
397 
399  if (genevt->get_collision_vertex().z()>=0.)
400  {
401  if (gRandom->Uniform(0,1)> exp (-(zpos-genevt->get_collision_vertex().z())/pion_lambda_Fe)) continue;//decay length cut
402  }
403 
404 
405  if (genevt->get_collision_vertex().z()<0. && zpos>0.)
406  {
407  if (Transverse_pos> FMAG_hole_R)
408  {
409  if (gRandom->Uniform(0,1) > exp (-zpos/pion_lambda_Fe)) continue;
410  }
411  else
412  {
413  if (gRandom->Uniform(0,1) > exp (-(zpos-FMAG_hole_L)/pion_lambda_Fe)) continue;
414  }
415  }
416 
417 
418  //Filter to get decay muons having enough energy to pass trough FMAG
419  if (zpos<0. && totMom/(FMAG_Z)<0.01) continue;
420  if (zpos>0. && totMom/(FMAG_Z-zpos)<0.01) continue;
421  }//@
422 
423  ineve->AddParticle(vtxindex, particle);
424 
425  if (embed_flag != 0) ineve->AddEmbeddedParticle(particle, embed_flag);
426  }
427  } // if (!finalstateparticles.empty())
428 
429  } // for (HepMC::GenEvent::vertex_iterator v = evt->vertices_begin();
430 
431  } // For pile-up simulation: loop end for PHHepMC event map
432 
433 
434  if (verbosity > 0) ineve->identify();
435 
437 }
438 
439 double HepMCNodeReader::smeargauss(const double width)
440 {
441  if (width == 0) return 0;
442  return gsl_ran_gaussian(RandomGenerator, width);
443 }
444 
446 void HepMCNodeReader::enable_position_filter(const double x_min, const double x_max, const double y_min, const double y_max, const double z_min, const double z_max)
447 {
448  _position_filter_on = true;
449  _pos_filter_x_min = x_min;
450  _pos_filter_x_max = x_max;
451  _pos_filter_y_min = y_min;
452  _pos_filter_y_max = y_max;
453  _pos_filter_z_min = z_min;
454  _pos_filter_z_max = z_max;
455 }
456 
457 bool HepMCNodeReader::PassParticleFilter(HepMC::GenParticle* p) {
458  for(int pid : _particle_filter_pid) {
459  if(p->pdg_id() == pid) return true;
460  }
461  return false;
462 }
463 
464 double HepMCNodeReader::smearflat(const double width)
465 {
466  if (width == 0) return 0;
467  return 2.0 * width * (gsl_rng_uniform_pos(RandomGenerator) - 0.5);
468 }
469 
470 void HepMCNodeReader::VertexPosition(const double v_x, const double v_y,
471  const double v_z)
472 {
473  cout << "HepMCNodeReader::VertexPosition - WARNING - this function is depreciated. "
474  << "HepMCNodeReader::VertexPosition() move all HEPMC subevents to a new vertex location. "
475  << "This also leads to a different vertex is used for HepMC subevent in Geant4 than that recorded in the HepMCEvent Node."
476  << "Recommendation: the vertex shifts are better controlled for individually HEPMC subevents in Fun4AllHepMCInputManagers and event generators."
477  << endl;
478 
479  vertex_pos_x = v_x;
480  vertex_pos_y = v_y;
481  vertex_pos_z = v_z;
482  return;
483 }
484 
485 void HepMCNodeReader::SmearVertex(const double s_x, const double s_y,
486  const double s_z)
487 {
488  cout << "HepMCNodeReader::SmearVertex - WARNING - this function is depreciated. "
489  << "HepMCNodeReader::SmearVertex() smear each HEPMC subevents to a new vertex location. "
490  << "This also leads to a different vertex is used for HepMC subevent in Geant4 than that recorded in the HepMCEvent Node."
491  << "Recommendation: the vertex smears are better controlled for individually HEPMC subevents in Fun4AllHepMCInputManagers and event generators."
492  << endl;
493 
494  width_vx = s_x;
495  width_vy = s_y;
496  width_vz = s_z;
497  return;
498 }
499 
500 void HepMCNodeReader::Embed(const int i)
501 {
502  cout << "HepMCNodeReader::Embed - WARNING - this function is depreciated. "
503  << "Embedding IDs are controlled for individually HEPMC subevents in Fun4AllHepMCInputManagers and event generators."
504  << endl;
505 
506  return;
507 }
static IsStateFinal isfinal
int verbosity
The verbosity level. 0 means not verbose at all.
Definition: Fun4AllBase.h:75
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
virtual int Verbosity() const
Gets the verbosity of this module.
Definition: Fun4AllBase.h:64
void VertexPosition(const double v_x, const double v_y, const double v_z)
void enable_position_filter(const double x_min, const double x_max, const double y_min, const double y_max, const double z_min, const double z_max)
Enable and define a position filter. The x,y,z limits are in cm.
void SmearVertex(const double s_x, const double s_y, const double s_z)
int Init(PHCompositeNode *topNode)
bool PassParticleFilter(HepMC::GenParticle *p)
void Embed(const int i=1)
virtual ~HepMCNodeReader()
int process_event(PHCompositeNode *topNode)
HepMCNodeReader(const std::string &name="HEPMCREADER")
PHBoolean addNode(PHNode *)
virtual const std::string get_CharFlag(const std::string &flag) const
Definition: PHFlag.cc:13
virtual float get_FloatFlag(const std::string &name) const
Definition: PHFlag.cc:83
int AddParticle(const int vtxid, PHG4Particle *particle)
Definition: PHG4InEvent.cc:63
int AddVtx(const double x, const double y, const double z, const double t)
Definition: PHG4InEvent.cc:38
virtual void identify(std::ostream &os=std::cout) const
Definition: PHG4InEvent.cc:134
void AddEmbeddedParticle(PHG4Particle *particle, int flag)
Definition: PHG4InEvent.h:25
virtual void set_pid(const int i)
Definition: PHG4Particle.h:33
virtual void set_barcode(const int barcode)
Definition: PHG4Particle.h:39
virtual void set_py(const double x)
Definition: PHG4Particle.h:35
virtual void set_px(const double x)
Definition: PHG4Particle.h:34
virtual void set_pz(const double x)
Definition: PHG4Particle.h:36
PHHepMCGenEventMap is collection of HEPMC events input into this simulation map of embedding ID -> PH...
ConstReverseIter rend() const
ConstReverseIter rbegin() const
iterator from lowest ID to highest, i.e. signal to background
std::map< int, PHHepMCGenEvent * >::reverse_iterator ReverseIter
virtual HepMC::GenEvent * getEvent()
const HepMC::FourVector & get_collision_vertex() const
collision vertex position in the Hall coordinate system, use PHENIX units of cm, ns
int get_embedding_id() const
virtual void identify(std::ostream &os=std::cout) const
bool is_simulated() const
whether this event has been processed in Geant4 simulation
PHNode * findFirst(const std::string &, const std::string &)
static recoConsts * instance()
Definition: recoConsts.cc:7
#define PHWHERE
Definition: phool.h:23