13 #include <Geant4/G4ParticleTable.hh>
14 #include <Geant4/G4ParticleDefinition.hh>
16 #include <TLorentzVector.h>
20 #include <gsl/gsl_randist.h>
31 _vertex_func_x(Uniform),
32 _vertex_func_y(Uniform),
33 _vertex_func_z(Uniform),
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),
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.));
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;
197 double mmuon = 105.64e-3;
198 double melectron = 0.511e-3;
203 else if(
decay1.compare(
"mu+")==0 ||
decay1.compare(
"mu-")==0)
207 cout <<
"Do not recognize the decay type " <<
decay1 <<
" will assume e+" << endl;
216 else if(
decay2.compare(
"mu+")==0 ||
decay2.compare(
"mu-")==0)
220 cout <<
"Do not recognize the decay type " <<
decay2 <<
" will assume e-" << endl;
232 cout <<
"PHG4ParticleGeneratorVectorMeson::InitRun started." << endl;
234 trand =
new TRandom3();
236 cout <<
Name() <<
" random seed: " << iseed << endl;
237 trand->SetSeed(iseed);
241 cout <<
Name() <<
" histrand random seed: " << iseed << endl;
242 gRandom->SetSeed(iseed);
245 fsin =
new TF1(
"fsin",
"sin(x)",0,M_PI);
249 frap->SetParameter(0,1.0);
250 frap->SetParameter(1,0.0);
251 frap->SetParameter(2,0.8749);
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);
259 ineve = findNode::getClass<PHG4InEvent>(topNode,
"PHG4INEVENT");
271 cout <<
"PHG4ParticleGeneratorVectorMeson::InitRun endeded." << endl;
280 if (!
ineve) cout<<
" G4InEvent not found "<<endl;
290 pt =
fpt->GetRandom();
301 y =
frap->GetRandom();
316 double mt = sqrt(pow(mnow,2) + pow(pt,2));
317 double eta = TMath::ASinH(TMath::SinH(y)*mt/pt);
322 vm.SetPtEtaPhiM(pt,eta,phi,mnow);
341 for (std::map<unsigned int, std::string>::iterator it=decay1_names.begin(); it !=decay1_names.end() ; ++it){
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);
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;
358 cout <<
PHWHERE <<
"Decay particles && vertex info can't be found !!" << endl;
375 else if (decay_id==0)
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);
391 double th1 =
fsin->GetRandom();
397 double px1 = p1*sin(th1)*cos(phi1);
398 double py1 = p1*sin(th1)*sin(phi1);
399 double pz1 = p1*cos(th1);
401 v1.SetPxPyPzE(px1,py1,pz1,E1);
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);
413 TLorentzVector v2 = vm - v1;
422 vector<PHG4Particle *>::const_iterator iter;
437 cout << endl <<
"Output some sanity check info from PHG4ParticleGeneratorVectorMeson:" << endl;
439 cout <<
" Vertex for this event (X,Y,Z) is (" <<
vtx_x <<
", " <<
vtx_y <<
", " <<
vtx_z <<
")" << endl;
442 cout <<
" Decay particle 1:"
446 <<
" eta " << v1.PseudoRapidity()
447 <<
" phi " << v1.Phi()
448 <<
" theta " << v1.Theta()
450 <<
" mass " << v1.M()
454 cout <<
" Decay particle 2:"
458 <<
" eta " << v2.PseudoRapidity()
459 <<
" phi " << v2.Phi()
460 <<
" theta " << v2.Theta()
462 <<
" mass " << v2.M()
467 cout <<
" Vector meson input kinematics: mass " << vm.M()
471 <<
" eta " << vm.PseudoRapidity()
472 <<
" y " << vm.Rapidity()
479 TLorentzVector vreco = v1 + v2;
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()
506 PHG4ParticleGeneratorVectorMeson::smearvtx(
const double position,
const double width, FUNCTION dist)
const
508 double res = position;
513 else if (dist ==
Gaus)
int verbosity
The verbosity level. 0 means not verbose at all.
virtual const std::string Name() const
Returns the name of this module.
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)
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)
gsl_rng * RandomGenerator
double _vertex_size_width
void set_vertex_size_parameters(const double mean, const double width)
set the dimensions of the distribution of particles about the vertex
void set_mass(const double mass)
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)
int process_event(PHCompositeNode *topNode)
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_width(const double width)
void set_rapidity_range(const double y_min, const double y_max)
void set_mom_range(const double mom_min, const double mom_max)
FUNCTION _vertex_size_func_r
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
int InitRun(PHCompositeNode *topNode)
PHG4ParticleGeneratorVectorMeson(const std::string &name="PGUN")
FUNCTION
supported function distributions
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 &)