Class Reference for E1039 Core & Analysis Software
PHHepMCGenHelper.cc
Go to the documentation of this file.
1 // $Id: $
2 
11 #include "PHHepMCGenHelper.h"
12 
13 #include "PHHepMCGenEvent.h"
14 #include "PHHepMCGenEventMap.h"
15 
17 #include <phool/PHCompositeNode.h>
18 #include <phool/PHDataNode.h>
19 #include <phool/PHRandomSeed.h>
20 #include <phool/getClass.h>
21 
22 #include <HepMC/GenEvent.h>
23 #include <HepMC/IO_GenEvent.h>
25 
26 #include <gsl/gsl_const.h>
27 #include <gsl/gsl_randist.h>
28 #include <gsl/gsl_rng.h>
29 
30 #include <cassert>
31 #include <climits>
32 #include <iostream>
33 
34 using namespace std;
35 
37  : _vertex_func_x(Gaus)
38  , _vertex_func_y(Gaus)
39  , _vertex_func_z(Gaus)
40  , _vertex_func_t(Gaus)
41  , _vertex_x(0)
42  , _vertex_y(0)
43  , _vertex_z(0)
44  , _vertex_t(0)
45  , _vertex_width_x(0)
46  , _vertex_width_y(0)
47  , _vertex_width_z(0)
48  , _vertex_width_t(0)
49  , _embedding_id(0)
50  , _reuse_vertex(false)
51  , _reuse_vertex_embedding_id(numeric_limits<int>::min())
52  , _geneventmap(nullptr)
53  , _legacy_vertexgenerator(false)
54  ,_paratio(0)
55 {
56 
57  RandomGenerator = gsl_rng_alloc(gsl_rng_mt19937);
58  unsigned int seed = PHRandomSeed(); // fixed seed is handled in this function
59  gsl_rng_set(RandomGenerator, seed );
60  _vertexGen = new SQPrimaryVertexGen(); //Abi
61 }
62 
64 {
65  gsl_rng_free(RandomGenerator);
66  delete _vertexGen;
67 }
68 
71 {
72  PHNodeIterator iter(topNode);
73  PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
74  if (!dstNode)
75  {
76  cout << PHWHERE << "DST Node missing doing nothing" << endl;
78  }
79 
80  _geneventmap = findNode::getClass<PHHepMCGenEventMap>(dstNode, "PHHepMCGenEventMap");
81  if (!_geneventmap)
82  {
83  _geneventmap = new PHHepMCGenEventMap();
84  PHIODataNode<PHObject> *newmapnode = new PHIODataNode<PHObject>(_geneventmap, "PHHepMCGenEventMap", "PHObject");
85  dstNode->addNode(newmapnode);
86  }
87 
88  assert(_geneventmap);
89 
91  _vertexGen->InitRun(topNode);//abi
92 
94 }
95 
98 {
99  assert(evt);
100  assert(_geneventmap);
101 
102  PHHepMCGenEvent *genevent = _geneventmap->insert_event(_embedding_id);
103 
104  genevent->addEvent(evt);
105  move_vertex(genevent);
106 
107  return genevent;
108 }
109 
111 {
112  assert(genevent);
113 
114  assert(_vertex_width_x >= 0);
115 
117  if( _legacy_vertexgenerator)
118  {
119  TVector3 vtx = _vertexGen->generateVertex();
120  _paratio = _vertexGen->getPARatio();
121  genevent->moveVertex(vtx.X(), vtx.Y(), vtx.Z(), _vertex_func_t);
122  }
123 
124  if (_reuse_vertex)
125  {
126  assert(_geneventmap);
127 
128  PHHepMCGenEvent *vtx_evt =
129  _geneventmap->get(_reuse_vertex_embedding_id);
130 
131  if (!vtx_evt)
132  {
133  cout << "PHHepMCGenHelper::move_vertex - Fatal Error - the requested source subevent with embedding ID "
134  << _reuse_vertex_embedding_id << " does not exist. Current HepMCEventMap:";
135  _geneventmap->identify();
136  exit(11);
137  }
138 
139 
140  genevent->moveVertex(
141  vtx_evt->get_collision_vertex().x(),
142  vtx_evt->get_collision_vertex().y(),
143  vtx_evt->get_collision_vertex().z(),
144  vtx_evt->get_collision_vertex().t());
145  }
146 
147 
148  genevent->moveVertex(
149  (smear(_vertex_x, _vertex_width_x, _vertex_func_x)),
150  (smear(_vertex_y, _vertex_width_y, _vertex_func_y)),
151  (smear(_vertex_z, _vertex_width_z, _vertex_func_z)),
152  (smear(_vertex_t, _vertex_width_t, _vertex_func_t)));
153 
154 
155 }
156 
158 {
159  _vertex_func_x = x;
160  _vertex_func_y = y;
161  _vertex_func_z = z;
162  _vertex_func_t = t;
163  return;
164 }
165 
166 void PHHepMCGenHelper::set_vertex_distribution_mean(const double x, const double y, const double z, const double t)
167 {
168  _vertex_x = x;
169  _vertex_y = y;
170  _vertex_z = z;
171  _vertex_t = t;
172  return;
173 }
174 
175 void PHHepMCGenHelper::set_vertex_distribution_width(const double x, const double y, const double z, const double t)
176 {
177  _vertex_width_x = x;
178  _vertex_width_y = y;
179  _vertex_width_z = z;
180  _vertex_width_t = t;
181  return;
182 }
183 
184 double PHHepMCGenHelper::smear(const double position,
185  const double width,
186  VTXFUNC dist) const
187 {
188  assert(width >= 0);
189 
190  double res = position;
191 
192  if (width == 0)
193  return res;
194 
195  if (dist == Uniform)
196  {
197  res = (position - width) + 2 * gsl_rng_uniform_pos(RandomGenerator) * width;
198  }
199  else if (dist == Gaus)
200  {
201  res = position + gsl_ran_gaussian(RandomGenerator, width);
202  }
203  else
204  {
205  cout << "PHHepMCGenHelper::smear - FATAL Error - unknown vertex function " << dist << endl;
206  exit(10);
207  }
208  return res;
209 }
PHBoolean addNode(PHNode *)
PHHepMCGenEventMap is collection of HEPMC events input into this simulation map of embedding ID -> PH...
const PHHepMCGenEvent * get(int idkey) const
fetch event
PHHepMCGenEvent * insert_event(const int embedding_id, const PHHepMCGenEvent *event=nullptr)
insert a event with specific embedding ID
void identify(std::ostream &os=std::cout) const
const HepMC::FourVector & get_collision_vertex() const
collision vertex position in the Hall coordinate system, use PHENIX units of cm, ns
bool addEvent(HepMC::GenEvent *evt)
host an HepMC event
virtual void moveVertex(double x, double y, double z, double t=0)
move the collision vertex position in the Hall coordinate system, use PHENIX units of cm,...
void set_vertex_distribution_function(VTXFUNC x, VTXFUNC y, VTXFUNC z, VTXFUNC t)
toss a new vertex according to a Uniform or Gaus distribution
void set_vertex_distribution_width(const double x, const double y, const double z, const double t)
set the width of the vertex distribution function about the mean, use PHENIX units of cm,...
virtual ~PHHepMCGenHelper()
int create_node_tree(PHCompositeNode *topNode)
init interface nodes
VTXFUNC
supported function distributions
@ Gaus
normal distribution with sigma width set via set_vertex_distribution_width()
@ Uniform
uniform distribution with half width set via set_vertex_distribution_width()
void set_vertex_distribution_mean(const double x, const double y, const double z, const double t)
set the mean value of the vertex distribution, use PHENIX units of cm, ns
PHHepMCGenEvent * insert_event(HepMC::GenEvent *evt)
send HepMC::GenEvent to DST tree. This function takes ownership of evt
gsl_rng * RandomGenerator
void move_vertex(PHHepMCGenEvent *genevent)
move vertex according to vertex settings
double smear(const double position, const double width, VTXFUNC dist) const
PHNode * findFirst(const std::string &, const std::string &)
Class to generate the event vertex, based on the beam profile and the target+spectrometer materials g...
double getPARatio()
get the proton/neutron ratio of the piece, must be called after generateVertex
TVector3 generateVertex()
generate 3-D vertex position
int InitRun(PHCompositeNode *node)
Initialize at the begining of Run.
#define PHWHERE
Definition: phool.h:23