Class Reference for E1039 Core & Analysis Software
PHG4ParticleGeneratorBase.cc
Go to the documentation of this file.
2 #include "PHG4Particlev1.h"
3 
4 #include "PHG4InEvent.h"
6 #include "PHG4VtxPoint.h"
7 
10 
11 #include <phool/getClass.h>
12 #include <phool/recoConsts.h>
13 
14 #include <phool/PHCompositeNode.h>
15 #include <phool/PHIODataNode.h>
16 #include <phool/PHRandomSeed.h>
17 
18 #include <Geant4/G4ParticleDefinition.hh>
19 #include <Geant4/G4ParticleTable.hh>
20 
21 #include <gsl/gsl_rng.h>
22 #include <cassert>
23 
24 using namespace std;
25 
27  : SubsysReco(name)
28  , embedflag(0)
29  , reuse_existing_vertex(0)
30  , vtx_x(0)
31  , vtx_y(0)
32  , vtx_z(0)
33  , t0(0)
34 {
35  RandomGenerator = gsl_rng_alloc(gsl_rng_mt19937);
36  seed = PHRandomSeed(); // fixed seed is handled in this funtcion
37  gsl_rng_set(RandomGenerator, seed);
38  return;
39 }
40 
42 {
43  while (particlelist.begin() != particlelist.end())
44  {
45  delete particlelist.back();
46  particlelist.pop_back();
47  }
48  gsl_rng_free(RandomGenerator);
49  return;
50 }
51 
52 int PHG4ParticleGeneratorBase::get_pdgcode(const std::string &name) const
53 {
54  G4ParticleTable *particleTable = G4ParticleTable::GetParticleTable();
55  G4String particleName = name;
56  G4ParticleDefinition *particledef = particleTable->FindParticle(particleName);
57  if (particledef)
58  {
59  return particledef->GetPDGEncoding();
60  }
61  return 0;
62 }
63 
64 string
66 {
67  G4ParticleTable *particleTable = G4ParticleTable::GetParticleTable();
68  G4ParticleDefinition *particledef = particleTable->FindParticle(pdgcode);
69  if (particledef)
70  {
71  return particledef->GetParticleName();
72  }
73  // if we cannot find the particle definition we'll make it ia geantino
74  return "geantino";
75 }
76 
77 double
78 PHG4ParticleGeneratorBase::get_mass(const int pdgcode) const
79 {
80  G4ParticleTable *particleTable = G4ParticleTable::GetParticleTable();
81  G4ParticleDefinition *particledef = particleTable->FindParticle(get_pdgname(pdgcode));
82  if (particledef)
83  {
84  return particledef->GetPDGMass();
85  }
86  return 0;
87 }
88 
89 void PHG4ParticleGeneratorBase::set_name(const std::string &particle)
90 {
92  particlelist[0]->set_name(particle);
93  particlelist[0]->set_pid(get_pdgcode(particle));
94  return;
95 }
96 
98 {
100  particlelist[0]->set_pid(pid);
101 }
102 
103 void PHG4ParticleGeneratorBase::set_mom(const double x, const double y, const double z)
104 {
106  particlelist[0]->set_px(x);
107  particlelist[0]->set_py(y);
108  particlelist[0]->set_pz(z);
109  return;
110 }
111 
112 void PHG4ParticleGeneratorBase::set_vtx(const double x, const double y, const double z)
113 {
114  vtx_x = x;
115  vtx_y = y;
116  vtx_z = z;
117  return;
118 }
119 
121 {
122  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
123  if (!ineve)
124  {
125  PHNodeIterator iter(topNode);
126  PHCompositeNode *dstNode;
127  dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
128 
129  ineve = new PHG4InEvent();
130  PHIODataNode<PHObject> *newNode = new PHIODataNode<PHObject>(ineve, "PHG4INEVENT", "PHObject");
131  dstNode->addNode(newNode);
132  }
133  return 0;
134 }
135 
137 {
138  cout << PHWHERE << " " << Name() << " using empty process_event" << endl;
139  return 0;
140 }
141 
142 void PHG4ParticleGeneratorBase::Print(const std::string &what) const
143 {
144  vector<PHG4Particle *>::const_iterator iter;
145  int i = 0;
146  for (iter = particlelist.begin(); iter != particlelist.end(); ++iter)
147  {
148  cout << "particle " << i << endl;
149  (*iter)->identify();
150  i++;
151  }
152 }
153 
154 void PHG4ParticleGeneratorBase::AddParticle(const std::string &particle, const double x, const double y, const double z)
155 {
156  PHG4Particle *part = new PHG4Particlev1(particle, get_pdgcode(particle), x, y, z);
157  particlelist.push_back(part);
158 }
159 
160 void PHG4ParticleGeneratorBase::AddParticle(const int pid, const double x, const double y, const double z)
161 {
162  PHG4Particle *particle = new PHG4Particlev1();
163  particle->set_pid(pid);
164  particle->set_px(x);
165  particle->set_py(y);
166  particle->set_pz(z);
167  particlelist.push_back(particle);
168 }
169 
171 {
172  if (!particlelist.size())
173  {
174  PHG4Particle *part = new PHG4Particlev1();
175  particlelist.push_back(part);
176  }
177  return;
178 }
179 
181 {
182  if ((particle->get_name()).size() == 0) // no size -> empty name string
183  {
184  particle->set_name(get_pdgname(particle->get_pid()));
185  }
186  if (particle->get_pid() == 0)
187  {
188  particle->set_pid(get_pdgcode(particle->get_name()));
189  }
190  if (embedflag)
191  {
192  ineve->AddEmbeddedParticle(particle, embedflag);
193  }
194  return;
195 }
196 
197 void PHG4ParticleGeneratorBase::set_seed(const unsigned int iseed)
198 {
199  seed = iseed;
200  cout << Name() << " random seed: " << seed << endl;
201  gsl_rng_set(RandomGenerator, seed);
202 }
203 
205 {
207  {
208  return 0;
209  }
210 
211  PHHepMCGenEventMap *genevtmap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
212 
213  if (genevtmap)
214  {
215  // use the highest priority HepMC subevent's vertex, ie that with highest embedding ID
216 
218  genevtmap->rbegin();
219 
220  if (iter != genevtmap->rend())
221  {
222  const PHHepMCGenEvent * hepmc_evt = iter->second;
223 
224  assert(hepmc_evt);
225 
226  const HepMC::FourVector &vtx = hepmc_evt->get_collision_vertex();
227 
228  set_vtx(vtx.x(), vtx.y(), vtx.z());
229 
230  if (verbosity > 0) {
231  cout <<"PHG4ParticleGeneratorBase::ReuseExistingVertex - reuse PHHepMCGenEventMap vertex "
232  << vtx.x()<<", "<< vtx.y()<<", "<< vtx.z()<<" cm. Source event:"
233  <<endl;
234  hepmc_evt->identify();
235  }
236 
237  return 1;
238  }
239  }
240 
241  // try PHG4INEVENT
242  PHG4InEvent *ineve = findNode::getClass<PHG4InEvent>(topNode, "PHG4INEVENT");
243 
244  if (ineve->GetNVtx() > 0)
245  {
246  std::pair<std::map<int, PHG4VtxPoint *>::const_iterator,
247  std::map<int, PHG4VtxPoint *>::const_iterator>
248  range = ineve->GetVertices();
249  std::map<int, PHG4VtxPoint *>::const_iterator iter = range.first;
250  PHG4VtxPoint *vtx = iter->second;
251 
252  if (!vtx)
253  {
254  cout << PHWHERE << "::Error - PHG4SimpleEventGenerator expects an existing vertex in PHG4InEvent, but none exists" << endl;
255  exit(1);
256  }
257  if (verbosity > 0)
258  {
259  cout << PHWHERE << "::Info - use this primary vertex from PHG4InEvent:" << endl;
260  vtx->identify();
261  }
262 
263  set_vtx(vtx->get_x(), vtx->get_y(), vtx->get_z());
264  return 1;
265 
266  } // if (_ineve->GetNVtx() > 0) {
267  // try the next option for getting primary vertex
268 
269  PHG4TruthInfoContainer *truthInfoList = //
270  findNode::getClass<PHG4TruthInfoContainer>(topNode,
271  "G4TruthInfo");
272  if (truthInfoList)
273  {
274  // embed to vertexes as set by primary vertex from the truth container, e.g. when embedding into DST
275 
276  PHG4VtxPoint *vtx = truthInfoList->GetPrimaryVtx(1);
277 
278  if (!vtx)
279  {
280  truthInfoList->identify();
281  cout << PHWHERE << "::Error - PHG4SimpleEventGenerator expects an existing truth vertex in PHG4TruthInfoContainer, but none exists"
282  << endl;
283  exit(1);
284  }
285 
286  if (verbosity > 0)
287  {
288  cout << PHWHERE << "::Info - use this primary vertex from PHG4TruthInfoContainer:" << endl;
289  vtx->identify();
290  }
291 
292  set_vtx(vtx->get_x(), vtx->get_y(), vtx->get_z());
293  return 1;
294  }
295 
296  // I am out of options.....
297 
298  cout << PHWHERE << "::Error - PHG4SimpleEventGenerator expects an existing truth vertex, but none exists"
299  << endl;
300  exit(1);
301 
302  return 0;
303 }
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 GetNVtx() const
Definition: PHG4InEvent.h:36
std::pair< std::map< int, PHG4VtxPoint * >::const_iterator, std::map< int, PHG4VtxPoint * >::const_iterator > GetVertices() const
Definition: PHG4InEvent.cc:105
void AddEmbeddedParticle(PHG4Particle *particle, int flag)
Definition: PHG4InEvent.h:25
int get_pdgcode(const std::string &name) const
virtual int ReuseExistingVertex(PHCompositeNode *topNode)
std::vector< PHG4Particle * > particlelist
virtual void set_mom(const double x, const double y, const double z)
std::string get_pdgname(const int pdgcode) const
virtual void AddParticle(const std::string &particle, const double x, const double y, const double z)
virtual void set_vtx(const double x, const double y, const double z)
PHG4ParticleGeneratorBase(const std::string &name="GENERATORBASE")
void set_seed(const unsigned int iseed)
virtual int InitRun(PHCompositeNode *topNode)
double get_mass(const int pdgcode) const
void SetParticleId(PHG4Particle *particle, PHG4InEvent *ineve)
virtual int process_event(PHCompositeNode *topNode)
virtual void set_name(const std::string &particle="proton")
virtual void set_pid(const int pid)
virtual void Print(const std::string &what="ALL") const
virtual int get_pid() const
Definition: PHG4Particle.h:14
virtual void set_pid(const int i)
Definition: PHG4Particle.h:33
virtual std::string get_name() const
Definition: PHG4Particle.h:15
virtual void set_py(const double x)
Definition: PHG4Particle.h:35
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 identify(std::ostream &os=std::cout) const
PHG4VtxPoint * GetPrimaryVtx(const int vtxid)
virtual double get_y() const
Definition: PHG4VtxPoint.h:21
virtual double get_x() const
Definition: PHG4VtxPoint.h:20
virtual void identify(std::ostream &os=std::cout) const
Definition: PHG4VtxPoint.cc:9
virtual double get_z() const
Definition: PHG4VtxPoint.h:22
PHHepMCGenEventMap is collection of HEPMC events input into this simulation map of embedding ID -> PH...
ConstReverseIter rend() const
ConstReverseIter rbegin() const
iterator from lowest ID to highest, i.e. signal to background
std::map< int, PHHepMCGenEvent * >::const_reverse_iterator ConstReverseIter
const HepMC::FourVector & get_collision_vertex() const
collision vertex position in the Hall coordinate system, use PHENIX units of cm, ns
virtual void identify(std::ostream &os=std::cout) const
PHNode * findFirst(const std::string &, const std::string &)
#define PHWHERE
Definition: phool.h:23