Class Reference for E1039 Core & Analysis Software
PHG4ParticleGeneratorVectorMeson.cc
Go to the documentation of this file.
2 #include "PHG4Particlev1.h"
3 
4 #include "PHG4InEvent.h"
5 
6 #include <phool/getClass.h>
7 #include <phool/recoConsts.h>
8 
10 #include <phool/PHIODataNode.h>
11 #include <phool/PHRandomSeed.h>
12 
13 #include <Geant4/G4ParticleTable.hh>
14 #include <Geant4/G4ParticleDefinition.hh>
15 
16 #include <TLorentzVector.h>
17 #include <TF1.h>
18 #include <TRandom3.h>
19 
20 #include <gsl/gsl_randist.h>
21 
22 using namespace std;
23 
26  decay1_names(),
27  decay2_names(),
28  decay_vtx_offset_x(),
29  decay_vtx_offset_y(),
30  decay_vtx_offset_z(),
31  _vertex_func_x(Uniform),
32  _vertex_func_y(Uniform),
33  _vertex_func_z(Uniform),
34  _t0(0.0),
35  _vertex_x(0.0),
36  _vertex_y(0.0),
37  _vertex_z(0.0),
38  _vertex_width_x(0.0),
39  _vertex_width_y(0.0),
40  _vertex_width_z(0.0),
41  _vertex_offset_x(0.0),
42  _vertex_offset_y(0.0),
43  _vertex_offset_z(0.0),
44  _vertex_size_func_r(Uniform),
45  _vertex_size_mean(0.0),
46  _vertex_size_width(0.0),
47  read_vtx_from_hepmc(true),
48  y_min(0.),
49  y_max(0.),
50  eta_min(-1.0),
51  eta_max(1.0),
52  mom_min(0.0),
53  mom_max(10.0),
54  pt_min(4.),
55  pt_max(4.),
56  mass(9.46),
57  width(54.02e-6),
58  m1(0.511e-3),
59  m2(0.511e-3),
60  _histrand_init(0),
61  decay1("e+"),
62  decay2("e-"),
63  fsin(NULL),
64  frap(NULL),
65  fpt(NULL),
66  trand(NULL),
67  ineve(NULL)
68 {
69 
70  // From PDG:
71  // Upsilon 1S has mass 9.4603, width 54.02 keV
72  // Upsilon 2S has mass 10.0233, width 31.98 keV
73  // Upsilon 3S has mass 10.3552, width 20.32 keV
74  return;
75 }
76 
77 void
78 PHG4ParticleGeneratorVectorMeson::add_decay_particles(const std::string &name1, const std::string &name2, const unsigned int decay_id)
79 {
80  decay1_names.insert(std::pair<unsigned int, std::string>(decay_id, name1));
81  decay2_names.insert(std::pair<unsigned int, std::string>(decay_id, name2));
82  decay_vtx_offset_x.insert(std::pair<unsigned int, double>(decay_id, 0.));
83  decay_vtx_offset_y.insert(std::pair<unsigned int, double>(decay_id, 0.));
84  decay_vtx_offset_z.insert(std::pair<unsigned int, double>(decay_id, 0.));
85  return;
86 }
87 
88 void
89 PHG4ParticleGeneratorVectorMeson::set_decay_vertex_offset(double dx, double dy, double dz, const unsigned int decay_id)
90 {
91  decay_vtx_offset_x.find(decay_id)->second = dx;
92  decay_vtx_offset_y.find(decay_id)->second = dy;
93  decay_vtx_offset_z.find(decay_id)->second = dz;
94  return;
95 }
96 
97 void
98 PHG4ParticleGeneratorVectorMeson::set_eta_range(const double min, const double max)
99 {
100  eta_min = min;
101  eta_max = max;
102  return;
103 }
104 
105 
106 void
107 PHG4ParticleGeneratorVectorMeson::set_rapidity_range(const double min, const double max)
108 {
109  y_min = min;
110  y_max = max;
111  return;
112 }
113 
114 
115 void
116 PHG4ParticleGeneratorVectorMeson::set_mom_range(const double min, const double max)
117 {
118  mom_min = min;
119  mom_max = max;
120  return;
121 }
122 
123 void
124 PHG4ParticleGeneratorVectorMeson::set_pt_range(const double min, const double max)
125 {
126  pt_min = min;
127  pt_max = max;
128  return;
129 }
130 
131 void
133  _vertex_func_x = x;
134  _vertex_func_y = y;
135  _vertex_func_z = z;
136  return;
137 }
138 
139 
140 void
141 PHG4ParticleGeneratorVectorMeson::set_vertex_distribution_mean(const double x, const double y, const double z) {
142  _vertex_x = x;
143  _vertex_y = y;
144  _vertex_z = z;
145  return;
146 }
147 
148 
149 void
150 PHG4ParticleGeneratorVectorMeson::set_vertex_distribution_width(const double x, 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 
157 void
158 PHG4ParticleGeneratorVectorMeson::set_existing_vertex_offset_vector(const double x, 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 
165 void
168  return;
169 }
170 
171 void
172 PHG4ParticleGeneratorVectorMeson::set_vertex_size_parameters(const double mean, const double width) {
173  _vertex_size_mean = mean;
175  return;
176 }
177 
178 
179 void
181 {
182  mass = mass_in;
183  return;
184 }
185 
186 void
188 {
189  width = width_in;
190  return;
191 }
192 
193 
194 void
195 PHG4ParticleGeneratorVectorMeson::set_decay_types(const std::string &name1, const std::string &name2)
196 {
197  double mmuon = 105.64e-3;
198  double melectron = 0.511e-3;
199 
200  decay1 = name1;
201  if(decay1.compare("e+")==0 || decay1.compare("e-")==0)
202  m1 = melectron;
203  else if(decay1.compare("mu+")==0 || decay1.compare("mu-")==0)
204  m1 = mmuon;
205  else
206  {
207  cout << "Do not recognize the decay type " << decay1 << " will assume e+" << endl;
208  decay1 = "e+";
209  m1 = melectron;
210  }
211 
212 
213  decay2 = name2;
214  if(decay2.compare("e+")==0 || decay2.compare("e-")==0)
215  m2 = melectron;
216  else if(decay2.compare("mu+")==0 || decay2.compare("mu-")==0)
217  m2 = mmuon;
218  else
219  {
220  cout << "Do not recognize the decay type " << decay2 << " will assume e-" << endl;
221  decay2 = "e-";
222  m2 = melectron;
223  }
224 
225  return;
226 }
227 
228 
229 int
231 {
232  cout << "PHG4ParticleGeneratorVectorMeson::InitRun started." << endl;
233 
234  trand = new TRandom3();
235  unsigned int iseed = PHRandomSeed(); // fixed seed handles in PHRandomSeed()
236  cout << Name() << " random seed: " << iseed << endl;
237  trand->SetSeed(iseed);
238  if (_histrand_init)
239  {
240  iseed = PHRandomSeed();
241  cout << Name() << " histrand random seed: " << iseed << endl;
242  gRandom->SetSeed(iseed);
243  }
244 
245  fsin = new TF1("fsin","sin(x)",0,M_PI);
246 
247  // From a fit to Pythia rapidity distribution for Upsilon(1S)
248  frap = new TF1("frap","gaus(0)",y_min,y_max);
249  frap->SetParameter(0,1.0);
250  frap->SetParameter(1,0.0);
251  frap->SetParameter(2,0.8749);
252 
253  // The dN/dPT distribution is described by:
254  fpt = new TF1("fpt","2.0*3.14159*x*[0]*pow((1 + x*x/(4*[1]) ),-[2])",pt_min,pt_max);
255  fpt->SetParameter(0,72.1);
256  fpt->SetParameter(1,26.516);
257  fpt->SetParameter(2,10.6834);
258 
259  ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
260  if (!ineve) {
261  PHNodeIterator iter( topNode );
262  PHCompositeNode *dstNode;
263  dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
264 
265  ineve = new PHG4InEvent();
266  PHIODataNode<PHObject> *newNode = new PHIODataNode<PHObject>(ineve, "PHG4INEVENT", "PHObject");
267  dstNode->addNode(newNode);
268  }
269 
270 
271  cout << "PHG4ParticleGeneratorVectorMeson::InitRun endeded." << endl;
272  return 0;
273 }
274 
275 
276 
277 int
279 {
280  if (!ineve) cout<<" G4InEvent not found "<<endl;
281 
282  // Generate a new set of vectors for the vector meson for each event
283  // These are the momentum and direction vectors for the pre-decay vector meson
284 
285  // taken randomly from a fitted pT distribution to Pythia Upsilons
286 
287  double pt = 0.0;
288  if(pt_max !=pt_min)
289  {
290  pt = fpt->GetRandom();
291  }
292  else
293  {
294  pt = pt_min;
295  }
296  // taken randomly from a fitted rapidity distribution to Pythia Upsilons
297 
298  double y = 0.0;
299  if(y_max != y_min)
300  {
301  y = frap->GetRandom();
302  }
303  else
304  {
305  y = y_min;
306  }
307  // 0 and 2*M_PI identical, so use gsl_rng_uniform which excludes 1.0
308  double phi = (2.0*M_PI)*gsl_rng_uniform(RandomGenerator);
309 
310  // The mass of the meson is taken from a Breit-Wigner lineshape
311 
312  double mnow = trand->BreitWigner(mass,width);
313 
314  // Get the pseudorapidity, eta, from the rapidity, mass and pt
315 
316  double mt = sqrt(pow(mnow,2) + pow(pt,2));
317  double eta = TMath::ASinH(TMath::SinH(y)*mt/pt);
318 
319  // Put it in a TLorentzVector
320 
321  TLorentzVector vm;
322  vm.SetPtEtaPhiM(pt,eta,phi,mnow);
323 
324  int vtxindex = -9;
325 
326  if (!ReuseExistingVertex(topNode))
327  {
328  // If not reusing existing vertex Randomly generate vertex position in z
329 
330  // mean width
334  }
338 
339 
340 
341  for (std::map<unsigned int, std::string>::iterator it=decay1_names.begin(); it !=decay1_names.end() ; ++it){
342 
343  unsigned int decay_id = it->first;
344  std::string decay1_name = it->second;
345  std::string decay2_name;
346  std::map<unsigned int, std::string>::iterator jt = decay2_names.find(decay_id);
347  std::map<unsigned int, double>::iterator xt = decay_vtx_offset_x.find(decay_id);
348  std::map<unsigned int, double>::iterator yt = decay_vtx_offset_y.find(decay_id);
349  std::map<unsigned int, double>::iterator zt = decay_vtx_offset_z.find(decay_id);
350 
351 
352  if (jt != decay2_names.end() && xt != decay_vtx_offset_x.end() && yt != decay_vtx_offset_y.end() && zt != decay_vtx_offset_z.end()){
353  decay2_name = jt->second;
354  set_decay_types(decay1_name, decay2_name);
355  set_existing_vertex_offset_vector(xt->second, yt->second, zt->second);
356 
357  } else{
358  cout << PHWHERE << "Decay particles && vertex info can't be found !!" << endl;
359  exit(1);
360  }
361 
362  // 3D Randomized vertex
363  if ((_vertex_size_width > 0.0)||(_vertex_size_mean != 0.0)) {
364  _vertex_size_mean = sqrt(pow(vtx_x,2)+pow(vtx_y,2)+pow(vtx_z,2));
366  double x = 0.0;
367  double y = 0.0;
368  double z = 0.0;
369  gsl_ran_dir_3d(RandomGenerator,&x,&y,&z);
370  x *= r;
371  y *= r;
372  z *= r;
373  vtxindex = ineve->AddVtx(vtx_x+x,vtx_y+y,vtx_z+z,t0);
374  }
375  else if (decay_id==0)
376  {
377  vtxindex = ineve->AddVtx(vtx_x,vtx_y,vtx_z,t0);
378  }
379 
380 
381  // Now decay it
382  // Get the decay energy and momentum in the frame of the vector meson - this correctly handles decay particles of any mass.
383 
384  double E1 = (pow(mnow,2) - pow(m2,2) + pow(m1,2)) / (2.0 * mnow);
385  double p1 = sqrt( ( pow(mnow,2) - pow(m1+m2,2) )*( pow(mnow,2) - pow(m1-m2,2) ) ) / (2.0 * mnow);
386 
387  // In the frame of the vector meson, get a random theta and phi angle for particle 1
388  // Assume angular distribution in the frame of the decaying meson that is uniform in phi and goes as sin(theta) in theta
389  // particle 2 has particle 1 momentum reflected through the origin
390 
391  double th1 = fsin->GetRandom();
392  // 0 and 2*M_PI identical, so use gsl_rng_uniform which excludes 1.0
393  double phi1 = 2.0*M_PI*gsl_rng_uniform(RandomGenerator);
394 
395  // Put particle 1 into a TLorentzVector
396 
397  double px1 = p1*sin(th1)*cos(phi1);
398  double py1 = p1*sin(th1)*sin(phi1);
399  double pz1 = p1*cos(th1);
400  TLorentzVector v1;
401  v1.SetPxPyPzE(px1,py1,pz1,E1);
402 
403  // now boost the decay product v1 into the lab using a vector consisting of the beta values of the vector meson
404  // where p/E is v/c if we use GeV/c for p and GeV for E
405 
406  double betax = vm.Px()/vm.E();
407  double betay = vm.Py()/vm.E();
408  double betaz = vm.Pz()/vm.E();
409  v1.Boost(betax,betay,betaz);
410 
411  // The second decay product's lab vector is the difference between the original meson and the boosted decay product 1
412 
413  TLorentzVector v2 = vm - v1;
414 
415  // Add the boosted decay particles to the particle list for the event
416 
417  AddParticle(decay1_name,v1.Px(),v1.Py(),v1.Pz());
418  AddParticle(decay2_name,v2.Px(),v2.Py(),v2.Pz());
419 
420  // Now output the list of boosted decay particles to the node tree
421 
422  vector<PHG4Particle *>::const_iterator iter;
423  for (iter = particlelist.begin(); iter != particlelist.end(); ++iter)
424  {
425  PHG4Particle *particle = new PHG4Particlev1(*iter);
426  SetParticleId(particle,ineve);
427  ineve->AddParticle(vtxindex, particle);
428  if(embedflag!=0) { ineve->AddEmbeddedParticle(particle,embedflag); }
429  }
430  // List what has been put into ineve for this event
431 
432  if(verbosity > 0)
433  {
434  ineve->identify();
435 
436  // Print some check output
437  cout << endl << "Output some sanity check info from PHG4ParticleGeneratorVectorMeson:" << endl;
438 
439  cout << " Vertex for this event (X,Y,Z) is (" << vtx_x << ", " << vtx_y << ", " << vtx_z << ")" << endl;
440  // Print the decay particle kinematics
441 
442  cout << " Decay particle 1:"
443  << " px " << v1.Px()
444  << " py " << v1.Py()
445  << " pz " << v1.Pz()
446  << " eta " << v1.PseudoRapidity()
447  << " phi " << v1.Phi()
448  << " theta " << v1.Theta()
449  << " pT " << v1.Pt()
450  << " mass " << v1.M()
451  << " E " << v1.E()
452  << endl;
453 
454  cout << " Decay particle 2:"
455  << " px " << v2.Px()
456  << " py " << v2.Py()
457  << " pz " << v2.Pz()
458  << " eta " << v2.PseudoRapidity()
459  << " phi " << v2.Phi()
460  << " theta " << v2.Theta()
461  << " pT " << v2.Pt()
462  << " mass " << v2.M()
463  << " E " << v2.E()
464  << endl;
465 
466  // Print the input vector meson kinematics
467  cout << " Vector meson input kinematics: mass " << vm.M()
468  << " px " << vm.Px()
469  << " py " << vm.Py()
470  << " pz " << vm.Pz()
471  << " eta " << vm.PseudoRapidity()
472  << " y " << vm.Rapidity()
473  << " pt " << vm.Pt()
474  << " E " << vm.E()
475  << endl;
476 
477  // Now, as a check, reconstruct the mass from the particle 1 and 2 kinematics
478 
479  TLorentzVector vreco = v1 + v2;
480 
481  cout << " Reco'd vector meson kinematics: mass " << vreco.M()
482  << " px " << vreco.Px()
483  << " py " << vreco.Py()
484  << " pz " << vreco.Pz()
485  << " eta " << vreco.PseudoRapidity()
486  << " y " << vreco.Rapidity()
487  << " pt " << vreco.Pt()
488  << " E " << vreco.E()
489  << endl;
490 
491 
492  }
493  } // decay particles
494 
495  // Reset particlelist for the next event
496  while(particlelist.begin() != particlelist.end())
497  {
498  delete particlelist.back();
499  particlelist.pop_back();
500  }
501 
502  return 0;
503 }
504 
505 double
506 PHG4ParticleGeneratorVectorMeson::smearvtx(const double position, const double width, FUNCTION dist) const
507 {
508  double res = position;
509  if (dist == Uniform)
510  {
511  res = (position-width) + 2*gsl_rng_uniform_pos(RandomGenerator)*width;
512  }
513  else if (dist == Gaus)
514  {
515  res = position + gsl_ran_gaussian(RandomGenerator,width);
516  }
517  return res;
518 }
519 
#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 *)
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 int ReuseExistingVertex(PHCompositeNode *topNode)
std::vector< PHG4Particle * > particlelist
virtual void AddParticle(const std::string &particle, const double x, const double y, const double z)
void SetParticleId(PHG4Particle *particle, PHG4InEvent *ineve)
void set_vertex_size_parameters(const double mean, const double width)
set the dimensions of the distribution of particles about the vertex
void set_existing_vertex_offset_vector(const double x, const double y, const double z)
set an offset vector from the existing vertex
void set_decay_types(const std::string &decay1, const std::string &decay2)
void set_vertex_size_function(FUNCTION r)
set the distribution function of particles about the vertex
void set_decay_vertex_offset(double dx, double dy, double dz, const unsigned int decay_id)
void set_eta_range(const double eta_min, const double eta_max)
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_rapidity_range(const double y_min, const double y_max)
void set_mom_range(const double mom_min, const double mom_max)
void set_pt_range(const double pt_min, const double pt_max)
void add_decay_particles(const std::string &name1, const std::string &name2, const unsigned int decay_id)
interface for adding particles by name
PHG4ParticleGeneratorVectorMeson(const std::string &name="PGUN")
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
PHNode * findFirst(const std::string &, const std::string &)
#define PHWHERE
Definition: phool.h:23