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