Class Reference for E1039 Core & Analysis Software
PHG4ParticleGeneratorD0.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 #include <gsl/gsl_const.h>
22 
23 using namespace std;
24 
27  vtx_zmin(-10.),
28  vtx_zmax(10),
29  y_min(0.),
30  y_max(0.),
31  eta_min(-1.0),
32  eta_max(1.0),
33  mom_min(0.0),
34  mom_max(10.0),
35  pt_min(4.),
36  pt_max(4.),
37  mass(1.86486),
38  m1(0.493677),
39  m2(0.13957018),
40  fsin(NULL),
41  frap(NULL),
42  fpt(NULL)
43 {
44  return;
45 }
46 
47 void
48 PHG4ParticleGeneratorD0::set_eta_range(const double min, const double max)
49 {
50  eta_min = min;
51  eta_max = max;
52  return;
53 }
54 
55 
56 void
57 PHG4ParticleGeneratorD0::set_rapidity_range(const double min, const double max)
58 {
59  y_min = min;
60  y_max = max;
61  return;
62 }
63 
64 
65 void
66 PHG4ParticleGeneratorD0::set_mom_range(const double min, const double max)
67 {
68  mom_min = min;
69  mom_max = max;
70  return;
71 }
72 
73 void
74 PHG4ParticleGeneratorD0::set_pt_range(const double min, const double max)
75 {
76  pt_min = min;
77  pt_max = max;
78  return;
79 }
80 
81 void
82 PHG4ParticleGeneratorD0::set_vtx_zrange(const double zmin, const double zmax)
83 {
84  vtx_zmin = zmin;
85  vtx_zmax = zmax;
86 
87  return;
88 }
89 
90 void
91 PHG4ParticleGeneratorD0::set_mass(const double mass_in)
92 {
93  mass = mass_in;
94  return;
95 }
96 
97 int
99 {
100  unsigned int iseed = PHRandomSeed(); // fixed seed handled in PHRandomSeed()
101  cout << Name() << " random seed: " << iseed << endl;
102  gRandom->SetSeed(iseed);
103 
104  fsin = new TF1("fsin","sin(x)",0,M_PI);
105 
106  // From a fit to Pythia rapidity distribution for Upsilon(1S)
107  frap = new TF1("frap","gaus(0)",y_min,y_max);
108  frap->SetParameter(0,1.0);
109  frap->SetParameter(1,0.0);
110  frap->SetParameter(2,0.8749);
111 
112 
113  // The dN/dPT distribution is described by:
114  fpt = new TF1("fpt","[0]*x*pow((1.0 + x*x/[1]/[1]),[2])",pt_min,pt_max);
115  fpt->SetParameter(0,1.16930e+04);
116  fpt->SetParameter(1,3.03486e+00);
117  fpt->SetParameter(2,-5.42819e+00);
118 
119 
120  return 0;
121 }
122 
123 int
125 {
126  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode,"PHG4INEVENT");
127 
128  double tau = 410.1e-15; // D0 life time in seconds
129  double cc = GSL_CONST_CGSM_SPEED_OF_LIGHT; // speed of light cm/sec
130  double ctau = tau*cc*1.0e+04; // ctau in micrometers
131 
132  // If not reusing existing vertex Randomly generate vertex position in z
133  if (! ReuseExistingVertex(topNode))
134  {
135 
136  if (vtx_zmax != vtx_zmin)
137  {
138  vtx_z = (vtx_zmax - vtx_zmin) * gsl_rng_uniform_pos(RandomGenerator) + vtx_zmin;
139  }
140  else
141  {
142  vtx_z = vtx_zmin;
143  }
144  }
145  // taken randomly from a fitted pT distribution to Pythia Upsilons
146 
147  double pt = 0.0;
148  if(pt_max !=pt_min)
149  {
150  pt = fpt->GetRandom();
151  }
152  else
153  {
154  pt = pt_min;
155  }
156 
157  // taken randomly from a fitted rapidity distribution to Pythia Upsilons
158 
159  double y = 0.0;
160  if(y_max != y_min)
161  {
162  y = frap->GetRandom();
163  }
164  else
165  {
166  y = y_min;
167  }
168  // 0 and 2*M_PI identical, so use gsl_rng_uniform which excludes 1.0
169  double phi = (2.0*M_PI)*gsl_rng_uniform(RandomGenerator);
170 
171  double mnow = mass;
172 
173  // Get the pseudorapidity, eta, from the rapidity, mass and pt
174 
175  double mt = sqrt(pow(mnow,2) + pow(pt,2));
176  double eta = TMath::ASinH(TMath::SinH(y)*mt/pt);
177 
178  // Put it in a TLorentzVector
179 
180  TLorentzVector vd0;
181  vd0.SetPtEtaPhiM(pt,eta,phi,mnow);
182 
183  double beta = vd0.Beta();
184  double gamma = vd0.Gamma();
185  double lifetime = gsl_ran_exponential(RandomGenerator,410.1e-03)* 1.0e-12; // life time in seconds
186  double lifepath = lifetime*gamma*beta*cc; // path in cm
187  if (verbosity > 0)
188  {
189  cout << "D0 px,py,pz: " << vd0.Px() << " " << vd0.Py() << " " << vd0.Pz() << " " << beta << " " << gamma << endl;
190  cout << " ctau = " << ctau << " " << lifetime << " " << lifepath << " " << lifepath*1.0e+04 << endl;
191  }
192 
193  vtx_x = vd0.Px()/vd0.P()*lifepath;
194  vtx_y = vd0.Py()/vd0.P()*lifepath;
195  vtx_z = vtx_z + vd0.Pz()/vd0.P()*lifepath;
196  t0 = lifetime;
197  int vtxindex = ineve->AddVtx(vtx_x,vtx_y,vtx_z,t0);
198  if (verbosity > 0)
199  {
200  cout << " XY vertex: " << sqrt(vtx_x*vtx_x+vtx_y*vtx_y) << " " << sqrt(vtx_x*vtx_x+vtx_y*vtx_y)*1.0e+04 << endl;
201  }
202 
203  // Now decay it
204  // Get the decay energy and momentum in the frame of the D0 - this correctly handles decay particles of any mass.
205 
206  double E1 = (pow(mnow,2) - pow(m2,2) + pow(m1,2)) / (2.0 * mnow);
207  double p1 = sqrt( ( pow(mnow,2) - pow(m1+m2,2) )*( pow(mnow,2) - pow(m1-m2,2) ) ) / (2.0 * mnow);
208 
209  // In the frame of the D0, get a random theta and phi angle for particle 1
210  // Assume angular distribution in the frame of the decaying D0 that is uniform in phi and goes as sin(theta) in theta
211  // particle 2 has particle 1 momentum reflected through the origin
212 
213  double th1 = fsin->GetRandom();
214  double phi1 = 2.0*M_PI*gsl_rng_uniform_pos(RandomGenerator);
215 
216  // Put particle 1 into a TLorentzVector
217 
218  double px1 = p1*sin(th1)*cos(phi1);
219  double py1 = p1*sin(th1)*sin(phi1);
220  double pz1 = p1*cos(th1);
221  TLorentzVector v1;
222  v1.SetPxPyPzE(px1,py1,pz1,E1);
223 
224  // now boost the decay product v1 into the lab using a vector consisting of the beta values of the vector meson
225  // where p/E is v/c if we use GeV/c for p and GeV for E
226 
227  double betax = vd0.Px()/vd0.E();
228  double betay = vd0.Py()/vd0.E();
229  double betaz = vd0.Pz()/vd0.E();
230  v1.Boost(betax,betay,betaz);
231 
232  // The second decay product's lab vector is the difference between the original D0 and the boosted decay product 1
233 
234  TLorentzVector v2 = vd0 - v1;
235 
236  // Add the boosted decay particles to the particle list for the event
237 
238  AddParticle(-321,v1.Px(),v1.Py(),v1.Pz()); // K-
239  AddParticle(211,v2.Px(),v2.Py(),v2.Pz()); // pi+
240 
241  // Now output the list of boosted decay particles to the node tree
242 
243  vector<PHG4Particle *>::const_iterator iter;
244  for (iter = particlelist.begin(); iter != particlelist.end(); ++iter)
245  {
246  PHG4Particle *particle = new PHG4Particlev1(*iter);
247  SetParticleId(particle,ineve);
248  ineve->AddParticle(vtxindex, particle);
249  if(embedflag!=0) { ineve->AddEmbeddedParticle(particle,embedflag); }
250  }
251 
252  // List what has been put into ineve for this event
253 
254  if(verbosity > 0)
255  {
256  ineve->identify();
257 
258  // Print some check output
259  cout << endl << "Output some sanity check info from PHG4ParticleGeneratorD0:" << endl;
260 
261  cout << "Event vertex: (" << vtx_x << ", " << vtx_y << ", " << vtx_z << ")" << endl;
262 
263  cout << "Kaon : " << v1.Pt() << " " << v1.PseudoRapidity() << " " << v1.M() << endl;
264  cout << "pion : " << v2.Pt() << " " << v2.PseudoRapidity() << " " << v2.M() << endl;
265  cout << "D0 : " << vd0.Pt() << " " << vd0.PseudoRapidity() << " " << vd0.M() << endl;
266  TLorentzVector vreco = v1 + v2;
267  cout << "reconstructed D0 : " << vreco.Pt() << " " << vreco.PseudoRapidity() << " " << vreco.M() << endl;
268 
269 /*
270  cout << " Decay particle 1:"
271  << " px " << v1.Px()
272  << " py " << v1.Py()
273  << " pz " << v1.Pz()
274  << " eta " << v1.PseudoRapidity()
275  << " phi " << v1.Phi()
276  << " theta " << v1.Theta()
277  << " pT " << v1.Pt()
278  << " mass " << v1.M()
279  << " E " << v1.E()
280  << endl;
281 
282  cout << " Decay particle 2:"
283  << " px " << v2.Px()
284  << " py " << v2.Py()
285  << " pz " << v2.Pz()
286  << " eta " << v2.PseudoRapidity()
287  << " phi " << v2.Phi()
288  << " theta " << v2.Theta()
289  << " pT " << v2.Pt()
290  << " mass " << v2.M()
291  << " E " << v2.E()
292  << endl;
293 
294  // Print the input vector meson kinematics
295  cout << " D0 input kinematics: mass " << vd0.M()
296  << " px " << vd0.Px()
297  << " py " << vd0.Py()
298  << " pz " << vd0.Pz()
299  << " eta " << vd0.PseudoRapidity()
300  << " y " << vd0.Rapidity()
301  << " pt " << vd0.Pt()
302  << " E " << vd0.E()
303  << endl;
304 
305  // Now, as a check, reconstruct the mass from the particle 1 and 2 kinematics
306 
307  TLorentzVector vreco = v1 + v2;
308 
309  cout << " Reconstructed D0 kinematics: mass " << vreco.M()
310  << " px " << vreco.Px()
311  << " py " << vreco.Py()
312  << " pz " << vreco.Pz()
313  << " eta " << vreco.PseudoRapidity()
314  << " y " << vreco.Rapidity()
315  << " pt " << vreco.Pt()
316  << " E " << vreco.E()
317  << endl;
318 */
319 
320  }
321 
322  // Reset particlelist for the next event
323  while(particlelist.begin() != particlelist.end())
324  {
325  delete particlelist.back();
326  particlelist.pop_back();
327  }
328 
329  return 0;
330 }
331 
#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
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_mom_range(const double mom_min, const double mom_max)
int process_event(PHCompositeNode *topNode)
void set_pt_range(const double pt_min, const double pt_max)
void set_eta_range(const double eta_min, const double eta_max)
PHG4ParticleGeneratorD0(const std::string &name="D0GEN")
void set_rapidity_range(const double y_min, const double y_max)
int InitRun(PHCompositeNode *topNode)
void set_vtx_zrange(const double zmin, const double zmax)
void set_mass(const double mass)