14 #include <Geant4/G4ParticleTable.hh>
15 #include <Geant4/G4ParticleDefinition.hh>
17 #include <gsl/gsl_randist.h>
29 _vertex_func_x(Uniform),
30 _vertex_func_y(Uniform),
31 _vertex_func_z(Uniform),
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),
55 _px_min(NAN), _px_max(NAN),
56 _py_min(NAN), _py_max(NAN),
57 _pz_min(NAN), _pz_max(NAN),
65 const unsigned int num) {
66 _particle_names.push_back(std::make_pair(name, num));
71 _particle_codes.push_back(std::make_pair(pid, num));
82 cout <<
"not setting eta bc etamin " << min <<
" > etamax: " << max << endl;
92 cout <<
"not setting phi bc phimin " << min <<
" > phimax: " << max << endl;
101 const double pt_gaus_width) {
103 cout <<
"not setting pt bc ptmin " << min <<
" > ptmax: " << max << endl;
106 assert(pt_gaus_width >= 0);
110 _pt_gaus_width = pt_gaus_width;
118 const double p_gaus_width) {
120 cout <<
"not setting p bc ptmin " << min <<
" > ptmax: " << max << endl;
123 assert(p_gaus_width >= 0);
126 _pt_gaus_width = NAN;
129 _p_gaus_width = p_gaus_width;
142 const double y,
const double z) {
150 const double y,
const double z) {
158 const double y,
const double z) {
159 _vertex_offset_x = x;
160 _vertex_offset_y = y;
161 _vertex_offset_z = z;
166 _vertex_size_func_r = r;
171 const double width) {
172 _vertex_size_mean = mean;
173 _vertex_size_width = width;
179 if ((_vertex_func_x !=
Uniform) && (_vertex_func_x !=
Gaus)) {
181 <<
"::Error - unknown vertex distribution function requested" << endl;
184 if ((_vertex_func_y !=
Uniform) && (_vertex_func_y !=
Gaus)) {
186 <<
"::Error - unknown vertex distribution function requested" << endl;
189 if ((_vertex_func_z !=
Uniform) && (_vertex_func_z !=
Gaus)) {
191 <<
"::Error - unknown vertex distribution function requested" << endl;
195 _ineve = findNode::getClass<PHG4InEvent>(topNode,
"PHG4INEVENT");
204 "PHG4INEVENT",
"PHObject");
210 <<
"================ PHG4PSScanGenerator::InitRun() ======================"
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;
218 for (
unsigned int i = 0; i < _particle_names.size(); ++i) {
219 cout <<
" " << _particle_names[i].first <<
", count = "
220 << _particle_names[i].second << endl;
224 <<
" Vertex Distribution: Set to reuse a previously generated sim vertex"
226 cout <<
" Vertex offset vector (x,y,z) = (" << _vertex_offset_x <<
","
227 << _vertex_offset_y <<
"," << _vertex_offset_z <<
")" << endl;
229 cout <<
" Vertex Distribution Function (x,y,z) = (";
232 else if (_vertex_func_x ==
Gaus)
236 else if (_vertex_func_y ==
Gaus)
240 else if (_vertex_func_z ==
Gaus)
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;
248 cout <<
" Vertex size function (r) = (";
249 if (_vertex_size_func_r ==
Uniform)
251 else if (_vertex_size_func_r ==
Gaus)
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;
261 <<
"==========================================================================="
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;
270 _particle_names.push_back(make_pair(pdgname, count));
280 <<
"====================== PHG4PSScanGenerator::process_event() ====================="
282 cout <<
"PHG4PSScanGenerator::process_event - reuse_existing_vertex = "
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);
296 vtx_x += _vertex_offset_x;
297 vtx_y += _vertex_offset_y;
298 vtx_z += _vertex_offset_z;
301 cout <<
"PHG4PSScanGenerator::process_event - vertex center"
308 for (
unsigned int i = 0; i < _particle_names.size(); ++i) {
310 string pdgname = _particle_names[i].first;
312 unsigned int nparticles = _particle_names[i].second;
314 for (
unsigned int j = 0; j < nparticles; ++j) {
316 if ((_vertex_size_width > 0.0) || (_vertex_size_mean != 0.0)) {
318 double r = smearvtx(_vertex_size_mean, _vertex_size_width,
319 _vertex_size_func_r);
330 }
else if ((i == 0) && (j == 0)) {
336 double px = std::numeric_limits<double>::max();
337 double py = std::numeric_limits<double>::max();
338 double pz = std::numeric_limits<double>::max();
340 if (std::isnan(_px_min)) {
341 double eta = (_eta_max - _eta_min)
343 double phi = (_phi_max - _phi_min)
347 if (!std::isnan(_p_min) && !std::isnan(_p_max)
348 && !std::isnan(_p_gaus_width)) {
352 }
else if (!std::isnan(_pt_min) && !std::isnan(_pt_max)
353 && !std::isnan(_pt_gaus_width)) {
358 <<
"Error: neither a p range or pt range was specified" << endl;
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;
370 int ibin = _event % (npx*npy*npz);
372 int ipz = ibin / (npx*npy);
373 int ipy = ( ibin % (npx*npy) ) / npx;
374 int ipx = ( ibin % (npx*npy) ) % npx;
376 px = _px_step*ipx + _px_min;
377 py = _py_step*ipy + _py_min;
378 pz = _pz_step*ipz + _pz_min;
380 std::cout << _event <<
": " << ipx <<
", " << ipy <<
", " << ipz << std::endl;
389 double e = sqrt(px * px + py * py + pz * pz + m * m);
411 <<
"======================================================================================"
421 const double x_max,
const double y_min,
const double y_max,
422 const double z_min,
const double z_max) {
431 double PHG4PSScanGenerator::smearvtx(
const double position,
const double width,
432 FUNCTION dist)
const {
433 double res = position;
435 res = (position - width) + 2 * gsl_rng_uniform_pos(
RandomGenerator) * width;
436 }
else if (dist ==
Gaus) {
int verbosity
The verbosity level. 0 means not verbose at all.
PHBoolean addNode(PHNode *)
int AddParticle(const int vtxid, PHG4Particle *particle)
int AddVtx(const double x, const double y, const double z, const double t)
virtual void identify(std::ostream &os=std::cout) const
void AddEmbeddedParticle(PHG4Particle *particle, int flag)
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_reuse_existing_vertex() const
int get_pdgcode(const std::string &name) const
virtual int ReuseExistingVertex(PHCompositeNode *topNode)
std::string get_pdgname(const int pdgcode) const
unsigned int get_seed() const
double get_mass(const int pdgcode) const
gsl_rng * RandomGenerator
int reuse_existing_vertex
virtual void set_parent_id(const int i)
virtual void set_vtx_id(const int i)
virtual void set_pid(const int i)
virtual void set_e(const double e)
virtual void set_py(const double x)
virtual void set_track_id(const int i)
virtual void set_name(const std::string &name)
virtual void set_px(const double x)
virtual void set_pz(const double x)
PHNode * findFirst(const std::string &, const std::string &)