Class Reference for E1039 Core & Analysis Software
PHG4PSScanGenerator.cc
Go to the documentation of this file.
1 #include "PHG4PSScanGenerator.h"
2 
3 #include "PHG4Particlev2.h"
4 #include "PHG4InEvent.h"
5 #include "PHG4VtxPoint.h"
7 
9 #include <phool/getClass.h>
10 
11 #include <phool/PHCompositeNode.h>
12 #include <phool/PHIODataNode.h>
13 
14 #include <Geant4/G4ParticleTable.hh>
15 #include <Geant4/G4ParticleDefinition.hh>
16 
17 #include <gsl/gsl_randist.h>
18 
19 #include <cstdlib>
20 #include <cmath>
21 #include <cassert>
22 
23 using namespace std;
24 
27  _particle_codes(),
28  _particle_names(),
29  _vertex_func_x(Uniform),
30  _vertex_func_y(Uniform),
31  _vertex_func_z(Uniform),
32  _t0(0.0),
33  _vertex_x(0.0),
34  _vertex_y(0.0),
35  _vertex_z(0.0),
36  _vertex_width_x(0.0),
37  _vertex_width_y(0.0),
38  _vertex_width_z(0.0),
39  _vertex_offset_x(0.0),
40  _vertex_offset_y(0.0),
41  _vertex_offset_z(0.0),
42  _vertex_size_func_r(Uniform),
43  _vertex_size_mean(0.0),
44  _vertex_size_width(0.0),
45  _eta_min(-1.25),
46  _eta_max(1.25),
47  _phi_min(-M_PI),
48  _phi_max(M_PI),
49  _pt_min(0.0),
50  _pt_max(10.0),
51  _pt_gaus_width(0.0),
52  _p_min(NAN),
53  _p_max(NAN),
54  _p_gaus_width(NAN),
55  _px_min(NAN), _px_max(NAN),
56  _py_min(NAN), _py_max(NAN),
57  _pz_min(NAN), _pz_max(NAN),
58  _ineve(NULL)
59 {
60  _event = 0;
61  return;
62 }
63 
64 void PHG4PSScanGenerator::add_particles(const std::string &name,
65  const unsigned int num) {
66  _particle_names.push_back(std::make_pair(name, num));
67  return;
68 }
69 
70 void PHG4PSScanGenerator::add_particles(const int pid, const unsigned int num) {
71  _particle_codes.push_back(std::make_pair(pid, num));
72  return;
73 }
74 
75 void PHG4PSScanGenerator::set_t0(const double t0) {
76  _t0 = t0;
77  return;
78 }
79 
80 void PHG4PSScanGenerator::set_eta_range(const double min, const double max) {
81  if (min > max) {
82  cout << "not setting eta bc etamin " << min << " > etamax: " << max << endl;
83  return;
84  }
85  _eta_min = min;
86  _eta_max = max;
87  return;
88 }
89 
90 void PHG4PSScanGenerator::set_phi_range(const double min, const double max) {
91  if (min > max) {
92  cout << "not setting phi bc phimin " << min << " > phimax: " << max << endl;
93  return;
94  }
95  _phi_min = min;
96  _phi_max = max;
97  return;
98 }
99 
100 void PHG4PSScanGenerator::set_pt_range(const double min, const double max,
101  const double pt_gaus_width) {
102  if (min > max) {
103  cout << "not setting pt bc ptmin " << min << " > ptmax: " << max << endl;
104  return;
105  }
106  assert(pt_gaus_width >= 0);
107 
108  _pt_min = min;
109  _pt_max = max;
110  _pt_gaus_width = pt_gaus_width;
111  _p_min = NAN;
112  _p_max = NAN;
113  _p_gaus_width = NAN;
114  return;
115 }
116 
117 void PHG4PSScanGenerator::set_p_range(const double min, const double max,
118  const double p_gaus_width) {
119  if (min > max) {
120  cout << "not setting p bc ptmin " << min << " > ptmax: " << max << endl;
121  return;
122  }
123  assert(p_gaus_width >= 0);
124  _pt_min = NAN;
125  _pt_max = NAN;
126  _pt_gaus_width = NAN;
127  _p_min = min;
128  _p_max = max;
129  _p_gaus_width = p_gaus_width;
130  return;
131 }
132 
134  FUNCTION y, FUNCTION z) {
135  _vertex_func_x = x;
136  _vertex_func_y = y;
137  _vertex_func_z = z;
138  return;
139 }
140 
142  const double y, const double z) {
143  _vertex_x = x;
144  _vertex_y = y;
145  _vertex_z = z;
146  return;
147 }
148 
150  const double y, const double z) {
151  _vertex_width_x = x;
152  _vertex_width_y = y;
153  _vertex_width_z = z;
154  return;
155 }
156 
158  const double y, const double z) {
159  _vertex_offset_x = x;
160  _vertex_offset_y = y;
161  _vertex_offset_z = z;
162  return;
163 }
164 
166  _vertex_size_func_r = r;
167  return;
168 }
169 
171  const double width) {
172  _vertex_size_mean = mean;
173  _vertex_size_width = width;
174  return;
175 }
176 
178 
179  if ((_vertex_func_x != Uniform) && (_vertex_func_x != Gaus)) {
180  cout << PHWHERE
181  << "::Error - unknown vertex distribution function requested" << endl;
183  }
184  if ((_vertex_func_y != Uniform) && (_vertex_func_y != Gaus)) {
185  cout << PHWHERE
186  << "::Error - unknown vertex distribution function requested" << endl;
188  }
189  if ((_vertex_func_z != Uniform) && (_vertex_func_z != Gaus)) {
190  cout << PHWHERE
191  << "::Error - unknown vertex distribution function requested" << endl;
193  }
194 
195  _ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
196  if (!_ineve) {
197  PHNodeIterator iter(topNode);
198  PHCompositeNode *dstNode;
199  dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode",
200  "DST"));
201 
202  _ineve = new PHG4InEvent();
203  PHDataNode<PHObject> *newNode = new PHDataNode<PHObject>(_ineve,
204  "PHG4INEVENT", "PHObject");
205  dstNode->addNode(newNode);
206  }
207 
208  if (verbosity > 0) {
209  cout
210  << "================ PHG4PSScanGenerator::InitRun() ======================"
211  << endl;
212  cout << " Random seed = " << get_seed() << endl;
213  cout << " Particles:" << endl;
214  for (unsigned int i = 0; i < _particle_codes.size(); ++i) {
215  cout << " " << _particle_codes[i].first << ", count = "
216  << _particle_codes[i].second << endl;
217  }
218  for (unsigned int i = 0; i < _particle_names.size(); ++i) {
219  cout << " " << _particle_names[i].first << ", count = "
220  << _particle_names[i].second << endl;
221  }
223  cout
224  << " Vertex Distribution: Set to reuse a previously generated sim vertex"
225  << endl;
226  cout << " Vertex offset vector (x,y,z) = (" << _vertex_offset_x << ","
227  << _vertex_offset_y << "," << _vertex_offset_z << ")" << endl;
228  } else {
229  cout << " Vertex Distribution Function (x,y,z) = (";
230  if (_vertex_func_x == Uniform)
231  cout << "Uniform,";
232  else if (_vertex_func_x == Gaus)
233  cout << "Gaus,";
234  if (_vertex_func_y == Uniform)
235  cout << "Uniform,";
236  else if (_vertex_func_y == Gaus)
237  cout << "Gaus,";
238  if (_vertex_func_z == Uniform)
239  cout << "Uniform";
240  else if (_vertex_func_z == Gaus)
241  cout << "Gaus";
242  cout << ")" << endl;
243  cout << " Vertex mean (x,y,z) = (" << _vertex_x << "," << _vertex_y << ","
244  << _vertex_z << ")" << endl;
245  cout << " Vertex width (x,y,z) = (" << _vertex_width_x << ","
246  << _vertex_width_y << "," << _vertex_width_z << ")" << endl;
247  }
248  cout << " Vertex size function (r) = (";
249  if (_vertex_size_func_r == Uniform)
250  cout << "Uniform";
251  else if (_vertex_size_func_r == Gaus)
252  cout << "Gaus";
253  cout << ")" << endl;
254  cout << " Vertex size (mean) = (" << _vertex_size_mean << ")" << endl;
255  cout << " Vertex size (width) = (" << _vertex_size_width << ")" << endl;
256  cout << " Eta range = " << _eta_min << " - " << _eta_max << endl;
257  cout << " Phi range = " << _phi_min << " - " << _phi_max << endl;
258  cout << " pT range = " << _pt_min << " - " << _pt_max << endl;
259  cout << " t0 = " << _t0 << endl;
260  cout
261  << "==========================================================================="
262  << endl;
263  }
264 
265  // the definition table should be filled now, so convert codes into names
266  for (unsigned int i = 0; i < _particle_codes.size(); ++i) {
267  int pdgcode = _particle_codes[i].first;
268  unsigned int count = _particle_codes[i].second;
269  string pdgname = get_pdgname(pdgcode);
270  _particle_names.push_back(make_pair(pdgname, count));
271  }
272 
274 }
275 
277 
278  if (verbosity > 0) {
279  cout
280  << "====================== PHG4PSScanGenerator::process_event() ====================="
281  << endl;
282  cout << "PHG4PSScanGenerator::process_event - reuse_existing_vertex = "
283  << reuse_existing_vertex << endl;
284  }
285 
286  // vtx_x, vtx_y and vtx_z are doubles from the base class
287  // common methods modify those, please no private copies
288  // at some point we might rely on them being up to date
289  if (!ReuseExistingVertex(topNode)) {
290  // generate a new vertex point
291  vtx_x = smearvtx(_vertex_x, _vertex_width_x, _vertex_func_x);
292  vtx_y = smearvtx(_vertex_y, _vertex_width_y, _vertex_func_y);
293  vtx_z = smearvtx(_vertex_z, _vertex_width_z, _vertex_func_z);
294  }
295 
296  vtx_x += _vertex_offset_x;
297  vtx_y += _vertex_offset_y;
298  vtx_z += _vertex_offset_z;
299 
300  if (verbosity > 0) {
301  cout << "PHG4PSScanGenerator::process_event - vertex center"
302  << reuse_existing_vertex << vtx_x << ", " << vtx_y << ", " << vtx_z
303  << " cm" << endl;
304  }
305 
306  int vtxindex = -1;
307  int trackid = -1;
308  for (unsigned int i = 0; i < _particle_names.size(); ++i) {
309 
310  string pdgname = _particle_names[i].first;
311  int pdgcode = get_pdgcode(pdgname);
312  unsigned int nparticles = _particle_names[i].second;
313 
314  for (unsigned int j = 0; j < nparticles; ++j) {
315 
316  if ((_vertex_size_width > 0.0) || (_vertex_size_mean != 0.0)) {
317 
318  double r = smearvtx(_vertex_size_mean, _vertex_size_width,
319  _vertex_size_func_r);
320 
321  double x = 0.0;
322  double y = 0.0;
323  double z = 0.0;
324  gsl_ran_dir_3d(RandomGenerator, &x, &y, &z);
325  x *= r;
326  y *= r;
327  z *= r;
328 
329  vtxindex = _ineve->AddVtx(vtx_x + x, vtx_y + y, vtx_z + z, _t0);
330  } else if ((i == 0) && (j == 0)) {
331  vtxindex = _ineve->AddVtx(vtx_x, vtx_y, vtx_z, _t0);
332  }
333 
334  ++trackid;
335 
336  double px = std::numeric_limits<double>::max();
337  double py = std::numeric_limits<double>::max();
338  double pz = std::numeric_limits<double>::max();
339 
340  if (std::isnan(_px_min)) {
341  double eta = (_eta_max - _eta_min)
342  * gsl_rng_uniform_pos(RandomGenerator) + _eta_min;
343  double phi = (_phi_max - _phi_min)
344  * gsl_rng_uniform_pos(RandomGenerator) + _phi_min;
345 
346  double pt;
347  if (!std::isnan(_p_min) && !std::isnan(_p_max)
348  && !std::isnan(_p_gaus_width)) {
349  pt = ((_p_max - _p_min) * gsl_rng_uniform_pos(RandomGenerator)
350  + _p_min + gsl_ran_gaussian(RandomGenerator, _p_gaus_width))
351  / cosh(eta);
352  } else if (!std::isnan(_pt_min) && !std::isnan(_pt_max)
353  && !std::isnan(_pt_gaus_width)) {
354  pt = (_pt_max - _pt_min) * gsl_rng_uniform_pos(RandomGenerator)
355  + _pt_min + gsl_ran_gaussian(RandomGenerator, _pt_gaus_width);
356  } else {
357  cout << PHWHERE
358  << "Error: neither a p range or pt range was specified" << endl;
359  exit(-1);
360  }
361 
362  px = pt * cos(phi);
363  py = pt * sin(phi);
364  pz = pt * sinh(eta);
365  } else {
366  int npx = (_px_max - _px_min)/_px_step + 1;
367  int npy = (_py_max - _py_min)/_py_step + 1;
368  int npz = (_pz_max - _pz_min)/_pz_step + 1;
369 
370  int ibin = _event % (npx*npy*npz);
371 
372  int ipz = ibin / (npx*npy);
373  int ipy = ( ibin % (npx*npy) ) / npx;
374  int ipx = ( ibin % (npx*npy) ) % npx;
375 
376  px = _px_step*ipx + _px_min;
377  py = _py_step*ipy + _py_min;
378  pz = _pz_step*ipz + _pz_min;
379 
380  std::cout << _event << ": " << ipx << ", " << ipy << ", " << ipz << std::endl;
381 
382  if(px > _px_max
383  or py > _py_max
384  or pz > _pz_max)
386  }
387 
388  double m = get_mass(pdgcode);
389  double e = sqrt(px * px + py * py + pz * pz + m * m);
390 
391  PHG4Particle *particle = new PHG4Particlev2();
392  particle->set_track_id(trackid);
393  particle->set_vtx_id(vtxindex);
394  particle->set_parent_id(0);
395  particle->set_name(pdgname);
396  particle->set_pid(pdgcode);
397  particle->set_px(px);
398  particle->set_py(py);
399  particle->set_pz(pz);
400  particle->set_e(e);
401 
402  _ineve->AddParticle(vtxindex, particle);
403  if (embedflag != 0)
404  _ineve->AddEmbeddedParticle(particle, embedflag);
405  }
406  }
407 
408  if (verbosity > 0) {
409  _ineve->identify();
410  cout
411  << "======================================================================================"
412  << endl;
413  }
414 
415  _event++;
416 
418 }
419 
421  const double x_max, const double y_min, const double y_max,
422  const double z_min, const double z_max) {
423  _px_min = x_min;
424  _px_max = x_max;
425  _py_min = y_min;
426  _py_max = y_max;
427  _pz_min = z_min;
428  _pz_max = z_max;
429 }
430 
431 double PHG4PSScanGenerator::smearvtx(const double position, const double width,
432  FUNCTION dist) const {
433  double res = position;
434  if (dist == Uniform) {
435  res = (position - width) + 2 * gsl_rng_uniform_pos(RandomGenerator) * width;
436  } else if (dist == Gaus) {
437  res = position + gsl_ran_gaussian(RandomGenerator, width);
438  }
439  return res;
440 }
#define NULL
Definition: Pdb.h:9
int verbosity
The verbosity level. 0 means not verbose at all.
Definition: Fun4AllBase.h:75
PHBoolean addNode(PHNode *)
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
void set_vertex_distribution_width(const double x, const double y, const double z)
set the width of the vertex distribution function about the mean
void set_p_range(const double p_min, const double p_max, const double p_gaus_width=0)
void set_eta_range(const double eta_min, const double eta_max)
range of randomized eta values
void set_pxpypz_range(const double x_min, const double x_max, const double y_min, const double y_max, const double z_min, const double z_max)
void set_vertex_size_function(FUNCTION r)
set the distribution function of particles about the vertex
FUNCTION
supported function distributions
void set_t0(const double t0)
set the starting time for the event
void set_vertex_distribution_mean(const double x, const double y, const double z)
set the mean value of the vertex distribution
void set_existing_vertex_offset_vector(const double x, const double y, const double z)
set an offset vector from the existing vertex
int process_event(PHCompositeNode *topNode)
void add_particles(const std::string &name, const unsigned int count)
interface for adding particles by name
void set_vertex_size_parameters(const double mean, const double width)
set the dimensions of the distribution of particles about the vertex
int InitRun(PHCompositeNode *topNode)
void set_vertex_distribution_function(FUNCTION x, FUNCTION y, FUNCTION z)
toss a new vertex according to a Uniform or Gaus distribution
void set_phi_range(const double phi_min, const double phi_max)
range of randomized phi values
PHG4PSScanGenerator(const std::string &name="EVTGENERATOR")
void set_pt_range(const double pt_min, const double pt_max, const double pt_gaus_width=0)
int get_pdgcode(const std::string &name) const
virtual int ReuseExistingVertex(PHCompositeNode *topNode)
std::string get_pdgname(const int pdgcode) const
double get_mass(const int pdgcode) const
virtual void set_parent_id(const int i)
Definition: PHG4Particle.h:30
virtual void set_vtx_id(const int i)
Definition: PHG4Particle.h:29
virtual void set_pid(const int i)
Definition: PHG4Particle.h:33
virtual void set_e(const double e)
Definition: PHG4Particle.h:37
virtual void set_py(const double x)
Definition: PHG4Particle.h:35
virtual void set_track_id(const int i)
Definition: PHG4Particle.h:28
virtual void set_name(const std::string &name)
Definition: PHG4Particle.h:32
virtual void set_px(const double x)
Definition: PHG4Particle.h:34
virtual void set_pz(const double x)
Definition: PHG4Particle.h:36
PHNode * findFirst(const std::string &, const std::string &)
#define PHWHERE
Definition: phool.h:23