CMS 3D CMS Logo

Generator.cc
Go to the documentation of this file.
4 
6 
8 
9 #include "HepPDT/ParticleID.hh"
10 
11 #include "G4Event.hh"
12 
13 #include "G4HEPEvtParticle.hh"
14 #include "G4Log.hh"
15 #include "G4ParticleDefinition.hh"
16 #include "G4PhysicalConstants.hh"
17 #include "G4SystemOfUnits.hh"
18 #include "G4UnitsTable.hh"
19 
20 #include <sstream>
21 
22 using namespace edm;
23 
25  : fPCuts(p.getParameter<bool>("ApplyPCuts")),
26  fPtransCut(p.getParameter<bool>("ApplyPtransCut")),
27  fEtaCuts(p.getParameter<bool>("ApplyEtaCuts")),
28  fPhiCuts(p.getParameter<bool>("ApplyPhiCuts")),
29  fFinalState(p.getParameter<bool>("FixOfFinalStateRadiation")),
30  theMinPhiCut(p.getParameter<double>("MinPhiCut")), // in radians (CMS standard)
31  theMaxPhiCut(p.getParameter<double>("MaxPhiCut")),
32  theMinEtaCut(p.getParameter<double>("MinEtaCut")),
33  theMaxEtaCut(p.getParameter<double>("MaxEtaCut")),
34  theMinPCut(p.getParameter<double>("MinPCut")), // in GeV (CMS standard)
35  theMaxPCut(p.getParameter<double>("MaxPCut")),
36  theEtaCutForHector(p.getParameter<double>("EtaCutForHector")),
37  verbose(p.getUntrackedParameter<int>("Verbosity", 0)),
38  fLumiFilter(nullptr),
39  evt_(nullptr),
40  vtx_(nullptr),
41  weight_(0),
42  Z_lmin(0),
43  Z_lmax(0),
44  Z_hector(0),
45  pdgFilterSel(false),
46  fPDGFilter(false) {
47  bool lumi = p.getParameter<bool>("ApplyLumiMonitorCuts");
48  if (lumi) {
50  }
51 
52  double theRDecLenCut = p.getParameter<double>("RDecLenCut") * cm;
53  theDecRCut2 = theRDecLenCut * theRDecLenCut;
54 
56 
57  double theDecLenCut = p.getParameter<double>("LDecLenCut") * cm;
58 
60 
61  pdgFilter.resize(0);
62  if (p.exists("PDGselection")) {
63  pdgFilterSel = (p.getParameter<edm::ParameterSet>("PDGselection")).getParameter<bool>("PDGfilterSel");
64  pdgFilter = (p.getParameter<edm::ParameterSet>("PDGselection")).getParameter<std::vector<int>>("PDGfilter");
65  if (!pdgFilter.empty()) {
66  fPDGFilter = true;
67  std::stringstream ss;
68  ss << "SimG4Core/Generator: ";
69  if (pdgFilterSel) {
70  ss << " Selecting only PDG ID = ";
71  } else {
72  ss << " Filtering out PDG ID = ";
73  }
74  for (unsigned int ii = 0; ii < pdgFilter.size(); ++ii) {
75  ss << pdgFilter[ii] << " ";
76  }
77  edm::LogVerbatim("SimG4CoreGenerator") << ss.str();
78  }
79  }
80 
81  if (fEtaCuts) {
82  Z_lmax = theRDecLenCut * ((1 - exp(-2 * theMaxEtaCut)) / (2 * exp(-theMaxEtaCut)));
83  Z_lmin = theRDecLenCut * ((1 - exp(-2 * theMinEtaCut)) / (2 * exp(-theMinEtaCut)));
84  }
85 
86  Z_hector = theRDecLenCut * ((1 - exp(-2 * theEtaCutForHector)) / (2 * exp(-theEtaCutForHector)));
87 
88  edm::LogVerbatim("SimG4CoreGenerator") << "SimG4Core/Generator: Rdecaycut= " << theRDecLenCut / cm
89  << " cm; Zdecaycut= " << theDecLenCut / cm << "Z_min= " << Z_lmin / cm
90  << " cm; Z_max= " << Z_lmax << " cm; Z_hector = " << Z_hector << " cm\n"
91  << " ApplyCuts: " << fFiductialCuts
92  << " PCuts: " << fPCuts << " PtransCut: " << fPtransCut
93  << " EtaCut: " << fEtaCuts << " PhiCut: " << fPhiCuts
94  << " LumiMonitorCut: " << lumi << " FinalState: " << fFinalState;
95  if (lumi) {
97  }
98 }
99 
101 
102 void Generator::HepMC2G4(const HepMC::GenEvent *evt_orig, G4Event *g4evt) {
103  HepMC::GenEvent *evt = new HepMC::GenEvent(*evt_orig);
104 
105  if (*(evt->vertices_begin()) == nullptr) {
106  std::stringstream ss;
107  ss << "SimG4Core/Generator: in event " << g4evt->GetEventID() << " Corrupted Event - GenEvent with no vertex \n";
108  throw SimG4Exception(ss.str());
109  }
110 
111  if (!evt->weights().empty()) {
112  weight_ = evt->weights()[0];
113  for (unsigned int iw = 1; iw < evt->weights().size(); ++iw) {
114  // terminate if the vector of weights contains a zero-weight
115  if (evt->weights()[iw] <= 0)
116  break;
117  weight_ *= evt->weights()[iw];
118  }
119  }
120 
121  if (vtx_ != nullptr) {
122  delete vtx_;
123  }
124  vtx_ = new math::XYZTLorentzVector((*(evt->vertices_begin()))->position().x(),
125  (*(evt->vertices_begin()))->position().y(),
126  (*(evt->vertices_begin()))->position().z(),
127  (*(evt->vertices_begin()))->position().t());
128 
129  edm::LogVerbatim("SimG4CoreGenerator") << "Generator: primary Vertex = (" << vtx_->x() << ", " << vtx_->y() << ", "
130  << vtx_->z() << ")";
131 
132  unsigned int ng4vtx = 0;
133 
134  for (HepMC::GenEvent::vertex_const_iterator vitr = evt->vertices_begin(); vitr != evt->vertices_end(); ++vitr) {
135  // loop for vertex, is it a real vertex?
136  // Set qvtx to true for any particles that should be propagated by GEANT,
137  // i.e., status 1 particles or status 2 particles that decay outside the
138  // beampipe.
139  G4bool qvtx = false;
140  HepMC::GenVertex::particle_iterator pitr;
141  for (pitr = (*vitr)->particles_begin(HepMC::children); pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
142  // For purposes of this function, the status is defined as follows:
143  // 1: particles are not decayed by generator
144  // 2: particles are decayed by generator but need to be propagated by
145  // GEANT 3: particles are decayed by generator but do not need to be
146  // propagated by GEANT
147  int status = (*pitr)->status();
148  int pdg = (*pitr)->pdg_id();
149  if (status > 3 && isExotic(pdg) && (!(isExoticNonDetectable(pdg)))) {
150  // In Pythia 8, there are many status codes besides 1, 2, 3.
151  // By setting the status to 2 for exotic particles, they will be
152  // checked: if its decay vertex is outside the beampipe, it will be
153  // propagated by GEANT. Some Standard Model particles, e.g., K0, cannot
154  // be propagated by GEANT, so do not change their status code.
155  status = 2;
156  }
157  if(std::abs(pdg) == 9900015) status = 3;
158 
159  // Particles which are not decayed by generator
160  if (status == 1) {
161  // filter out unwanted particles and vertices
162  if (fPDGFilter && !pdgFilterSel && IsInTheFilterList(pdg)) {
163  continue;
164  }
165 
166  qvtx = true;
167  if (verbose > 2)
168  LogDebug("SimG4CoreGenerator") << "GenVertex barcode = " << (*vitr)->barcode()
169  << " selected for GenParticle barcode = " << (*pitr)->barcode();
170  break;
171  }
172  // The selection is made considering if the partcile with status = 2
173  // have the end_vertex with a radius greater than the radius of beampipe
174  // cylinder (no requirement on the Z of the vertex is applyed).
175  else if (status == 2) {
176  if ((*pitr)->end_vertex() != nullptr) {
177  double xx = (*pitr)->end_vertex()->position().x();
178  double yy = (*pitr)->end_vertex()->position().y();
179  double r_dd = xx * xx + yy * yy;
180  if (r_dd > theDecRCut2) {
181  qvtx = true;
182  if (verbose > 2)
183  LogDebug("SimG4CoreGenerator")
184  << "GenVertex barcode = " << (*vitr)->barcode()
185  << " selected for GenParticle barcode = " << (*pitr)->barcode() << " radius = " << std::sqrt(r_dd);
186  break;
187  }
188  } else {
189  // particles with status 2 without end_vertex are
190  // equivalent to stable
191  qvtx = true;
192  break;
193  }
194  }
195  }
196 
197  // if this vertex is inside fiductial volume inside the beam pipe
198  // and has no long-lived secondary the vertex is not saved
199  if (!qvtx) {
200  continue;
201  }
202 
203  double x1 = (*vitr)->position().x() * mm;
204  double y1 = (*vitr)->position().y() * mm;
205  double z1 = (*vitr)->position().z() * mm;
206  double t1 = (*vitr)->position().t() * mm / c_light;
207 
208  G4PrimaryVertex *g4vtx = new G4PrimaryVertex(x1, y1, z1, t1);
209 
210  for (pitr = (*vitr)->particles_begin(HepMC::children); pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
211  int status = (*pitr)->status();
212  int pdg = (*pitr)->pdg_id();
213  bool hasDecayVertex = (nullptr != (*pitr)->end_vertex()) ? true : false;
214 
215  // Filter on allowed particle species if required
216  if (fPDGFilter) {
217  bool isInTheList = IsInTheFilterList(pdg);
218  if ((!pdgFilterSel && isInTheList) || (pdgFilterSel && !isInTheList)) {
219  edm::LogVerbatim("SimG4CoreGenerator")
220  << " Skiped GenParticle barcode= " << (*pitr)->barcode() << " PDGid= " << pdg << " status= " << status
221  << " isExotic: " << isExotic(pdg) << " isExoticNotDet: " << isExoticNonDetectable(pdg)
222  << " isInTheList: " << isInTheList << " hasDecayVertex: " << hasDecayVertex;
223  continue;
224  }
225  }
226 
227  edm::LogVerbatim("SimG4CoreGenerator")
228  << "Generator: pdg= " << pdg << " status= " << status << " hasPreDefinedDecay: " << hasDecayVertex
229  << " isExotic: " << isExotic(pdg) << " isExoticNotDet: " << isExoticNonDetectable(pdg)
230  << " isInTheList: " << IsInTheFilterList(pdg) << "\n (x,y,z,t): (" << x1 << "," << y1 << "," << z1
231  << "," << t1 << ")";
232 
233  if (status > 3 && isExotic(pdg) && (!(isExoticNonDetectable(pdg)))) {
234  status = hasDecayVertex ? 2 : 1;
235  }
236 
237  // this particle has predefined decay but has no vertex
238  if (2 == status && !hasDecayVertex) {
239  edm::LogWarning("SimG4CoreGenerator: in event ") << g4evt->GetEventID() << " a particle "
240  << " pdgid= " << pdg
241  << " has status=2 but has no decay vertex, so will be fully "
242  "tracked by Geant4";
243  status = 1;
244  }
245  // heavy neutrino should not be propagated by Geant4
246  if(std::abs(pdg) == 9900015) status = 3;
247 
248  double x2 = x1;
249  double y2 = y1;
250  double z2 = z1;
251  double decay_length = 0.0;
252  if (2 == status) {
253  x2 = (*pitr)->end_vertex()->position().x();
254  y2 = (*pitr)->end_vertex()->position().y();
255  z2 = (*pitr)->end_vertex()->position().z();
256  decay_length = std::sqrt((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1) + (z2 - z1) * (z2 - z1));
257  }
258 
259  bool toBeAdded = (fFiductialCuts) ? false : true;
260 
261  double px = (*pitr)->momentum().px();
262  double py = (*pitr)->momentum().py();
263  double pz = (*pitr)->momentum().pz();
264  double ptot = std::sqrt(px * px + py * py + pz * pz);
265  math::XYZTLorentzVector p(px, py, pz, (*pitr)->momentum().e());
266 
267  double ximpact = x1;
268  double yimpact = y1;
269  double zimpact = z1;
270 
271  // protection against numerical problems for extremely low momenta
272  // compute impact point at transition to Hector
273  const double minTan = 1.e-20;
274  if (std::abs(z1) < Z_hector && std::abs(pz) >= minTan * ptot) {
275  if (pz > 0.0) {
276  zimpact = Z_hector;
277  } else {
278  zimpact = -Z_hector;
279  }
280  double del = (zimpact - z1) / pz;
281  ximpact += del * px;
282  yimpact += del * py;
283  }
284  double rimpact2 = ximpact * ximpact + yimpact * yimpact;
285 
286  if (verbose > 2)
287  LogDebug("SimG4CoreGenerator") << "Processing GenParticle barcode= " << (*pitr)->barcode() << " pdg= " << pdg
288  << " status= " << (*pitr)->status() << " st= " << status
289  << " rimpact(cm)= " << std::sqrt(rimpact2) / cm
290  << " zimpact(cm)= " << zimpact / cm << " ptot(GeV)= " << ptot
291  << " pz(GeV)= " << pz;
292 
293  // Particles of status 1 trasnported along the beam pipe for forward
294  // detectors (HECTOR) always pass to Geant4 without cuts
295  if (1 == status && std::abs(zimpact) >= Z_hector && rimpact2 <= theDecRCut2) {
296  toBeAdded = true;
297  if (verbose > 2)
298  LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*pitr)->barcode() << " passed case 3";
299  } else {
300  // Standard case: particles not decayed by the generator
301  if (1 == status && (std::abs(zimpact) < Z_hector || rimpact2 > theDecRCut2)) {
302  // Ptot cut for all particles
303  if (fPCuts && (ptot < theMinPCut || ptot > theMaxPCut)) {
304  continue;
305  }
306  // phi cut is applied mainly for particle gun
307  if (fPhiCuts) {
308  double phi = p.phi();
309  if (phi < theMinPhiCut || phi > theMaxPhiCut) {
310  continue;
311  }
312  }
313  // eta cut is applied if position of the decay
314  // is within vacuum chamber and limited in Z
315  if (fEtaCuts) {
316  // eta cut
317  double xi = x1;
318  double yi = y1;
319  double zi = z1;
320 
321  // can be propagated along Z
322  if (std::abs(pz) >= minTan * ptot) {
323  if ((zi >= Z_lmax) & (pz < 0.0)) {
324  zi = Z_lmax;
325  } else if ((zi <= Z_lmin) & (pz > 0.0)) {
326  zi = Z_lmin;
327  } else {
328  if (pz > 0) {
329  zi = Z_lmax;
330  } else {
331  zi = Z_lmin;
332  }
333  }
334  double del = (zi - z1) / pz;
335  xi += del * px;
336  yi += del * py;
337  }
338  // check eta cut
339  if ((zi >= Z_lmin) & (zi <= Z_lmax) & (xi * xi + yi * yi < theDecRCut2)) {
340  continue;
341  }
342  }
343  if (fLumiFilter && !fLumiFilter->isGoodForLumiMonitor(*pitr)) {
344  continue;
345  }
346  toBeAdded = true;
347  if (verbose > 2)
348  LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*pitr)->barcode() << " passed case 1";
349 
350  // Decay chain outside the fiducial cylinder defined by theRDecLenCut
351  // are used for Geant4 tracking with predefined decay channel
352  // In the case of decay in vacuum particle is not tracked by Geant4
353  } else if (2 == status && x2 * x2 + y2 * y2 >= theDecRCut2 && std::abs(z2) < Z_hector) {
354  toBeAdded = true;
355  if (verbose > 2)
356  LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*pitr)->barcode() << " passed case 2"
357  << " decay_length(cm)= " << decay_length / cm;
358  }
359  }
360  if (toBeAdded) {
361  G4PrimaryParticle *g4prim = new G4PrimaryParticle(pdg, px * GeV, py * GeV, pz * GeV);
362 
363  if (g4prim->GetG4code() != nullptr) {
364  g4prim->SetMass(g4prim->GetG4code()->GetPDGMass());
365  double charge = g4prim->GetG4code()->GetPDGCharge();
366 
367  // apply Pt cut
368  if (fPtransCut && 0.0 != charge && px * px + py * py < theMinPtCut2) {
369  delete g4prim;
370  continue;
371  }
372  g4prim->SetCharge(charge);
373  }
374 
375  // V.I. do not use SetWeight but the same code
376  // value of the code compute inside TrackWithHistory
377  // g4prim->SetWeight( 10000*(*vpitr)->barcode() ) ;
378  setGenId(g4prim, (*pitr)->barcode());
379 
380  if (2 == status) {
381  particleAssignDaughters(g4prim, (HepMC::GenParticle *)*pitr, decay_length);
382  }
383  if (verbose > 1)
384  g4prim->Print();
385  g4vtx->SetPrimary(g4prim);
386  edm::LogVerbatim("SimG4CoreGenerator") << "Generator: new particle pdg= " << pdg << " Ptot(GeV/c)= " << ptot
387  << " Pt= " << std::sqrt(px * px + py * py);
388  }
389  }
390 
391  if (verbose > 1)
392  g4vtx->Print();
393  g4evt->AddPrimaryVertex(g4vtx);
394  ++ng4vtx;
395  }
396 
397  // Add a protection for completely empty events (produced by LHCTransport):
398  // add a dummy vertex with no particle attached to it
399  if (ng4vtx == 0) {
400  G4PrimaryVertex *g4vtx = new G4PrimaryVertex(0.0, 0.0, 0.0, 0.0);
401  if (verbose > 1)
402  g4vtx->Print();
403  g4evt->AddPrimaryVertex(g4vtx);
404  }
405 
406  delete evt;
407 }
408 
409 void Generator::particleAssignDaughters(G4PrimaryParticle *g4p, HepMC::GenParticle *vp, double decaylength) {
410  if (verbose > 1) {
411  LogDebug("SimG4CoreGenerator") << "Special case of long decay length \n"
412  << "Assign daughters with to mother with decaylength=" << decaylength / cm << " cm";
413  }
414  math::XYZTLorentzVector p(vp->momentum().px(), vp->momentum().py(), vp->momentum().pz(), vp->momentum().e());
415 
416  // defined preassigned decay time
417  double proper_time = decaylength / (p.Beta() * p.Gamma() * c_light);
418  g4p->SetProperTime(proper_time);
419 
420  if (verbose > 2) {
421  LogDebug("SimG4CoreGenerator") << " px= " << p.px() << " py= " << p.py() << " pz= " << p.pz() << " e= " << p.e()
422  << " beta= " << p.Beta() << " gamma= " << p.Gamma()
423  << " Proper time= " << proper_time / ns << " ns";
424  }
425 
426  // the particle will decay after the same length if it
427  // has not interacted before
428  double x1 = vp->end_vertex()->position().x();
429  double y1 = vp->end_vertex()->position().y();
430  double z1 = vp->end_vertex()->position().z();
431 
432  for (HepMC::GenVertex::particle_iterator vpdec = vp->end_vertex()->particles_begin(HepMC::children);
433  vpdec != vp->end_vertex()->particles_end(HepMC::children);
434  ++vpdec) {
435  // transform decay products such that in the rest frame of mother
437  (*vpdec)->momentum().px(), (*vpdec)->momentum().py(), (*vpdec)->momentum().pz(), (*vpdec)->momentum().e());
438 
439  // children should only be taken into account once
440  G4PrimaryParticle *g4daught =
441  new G4PrimaryParticle((*vpdec)->pdg_id(), pdec.x() * GeV, pdec.y() * GeV, pdec.z() * GeV);
442 
443  if (g4daught->GetG4code() != nullptr) {
444  g4daught->SetMass(g4daught->GetG4code()->GetPDGMass());
445  g4daught->SetCharge(g4daught->GetG4code()->GetPDGCharge());
446  }
447 
448  // V.I. do not use SetWeight but the same code
449  // value of the code compute inside TrackWithHistory
450  setGenId(g4daught, (*vpdec)->barcode());
451 
452  if (verbose > 2)
453  LogDebug("SimG4CoreGenerator") << "Assigning a " << (*vpdec)->pdg_id() << " as daughter of a " << vp->pdg_id();
454 
455  bool checkStatus = fFinalState ? (*vpdec)->status() > 3
456  : ((*vpdec)->status() == 23 && std::abs(vp->pdg_id()) == 1000015);
457  if (((*vpdec)->status() == 2 || checkStatus) &&
458  (*vpdec)->end_vertex() != nullptr) {
459  double x2 = (*vpdec)->end_vertex()->position().x();
460  double y2 = (*vpdec)->end_vertex()->position().y();
461  double z2 = (*vpdec)->end_vertex()->position().z();
462  double dd = std::sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2) + (z1 - z2) * (z1 - z2));
463  particleAssignDaughters(g4daught, *vpdec, dd);
464  }
465  (*vpdec)->set_status(1000 + (*vpdec)->status());
466  g4p->SetDaughter(g4daught);
467 
468  if (verbose > 1)
469  g4daught->Print();
470  }
471 }
472 
473 // Used for non-beam particles
474 bool Generator::particlePassesPrimaryCuts(const G4ThreeVector &p) const {
475  bool flag = true;
476  double ptot = p.mag();
477  if (fPCuts && (ptot < theMinPCut * GeV || ptot > theMaxPCut * GeV)) {
478  flag = false;
479  }
480  if (fEtaCuts && flag) {
481  double pz = p.z();
482  if (ptot < pz + 1.e-10) {
483  flag = false;
484 
485  } else {
486  double eta = 0.5 * G4Log((ptot + pz) / (ptot - pz));
487  if (eta < theMinEtaCut || eta > theMaxEtaCut) {
488  flag = false;
489  }
490  }
491  }
492  if (fPhiCuts && flag) {
493  double phi = p.phi();
494  if (phi < theMinPhiCut || phi > theMaxPhiCut) {
495  flag = false;
496  }
497  }
498 
499  if (verbose > 2)
500  LogDebug("SimG4CoreGenerator") << "Generator ptot(GeV)= " << ptot / GeV << " eta= " << p.eta()
501  << " phi= " << p.phi() << " Flag= " << flag;
502 
503  return flag;
504 }
505 
506 bool Generator::isExotic(int pdgcode) const {
507  int pdgid = std::abs(pdgcode);
508  return ((pdgid >= 1000000 && pdgid < 4000000 && pdgid != 3000022) || // SUSY, R-hadron, and technicolor particles
509  pdgid == 17 || // 4th generation lepton
510  pdgid == 34 || // W-prime
511  pdgid == 37) // charged Higgs
512  ? true
513  : false;
514 }
515 
516 bool Generator::isExoticNonDetectable(int pdgcode) const {
517  int pdgid = std::abs(pdgcode);
518  HepPDT::ParticleID pid(pdgcode);
519  int charge = pid.threeCharge();
520  return ((charge == 0) && (pdgid >= 1000000 && pdgid < 1000040)) // SUSY
521  ? true
522  : false;
523 }
524 
525 bool Generator::IsInTheFilterList(int pdgcode) const {
526  int pdgid = std::abs(pdgcode);
527  for (auto &pdg : pdgFilter) {
528  if (std::abs(pdg) == pdgid) {
529  return true;
530  }
531  }
532  return false;
533 }
534 
535 void Generator::nonBeamEvent2G4(const HepMC::GenEvent *evt, G4Event *g4evt) {
536  int i = 0;
537  for (HepMC::GenEvent::particle_const_iterator it = evt->particles_begin(); it != evt->particles_end(); ++it) {
538  ++i;
539  HepMC::GenParticle *gp = (*it);
540  int g_status = gp->status();
541  // storing only particle with status == 1
542  if (g_status == 1) {
543  int g_id = gp->pdg_id();
544  G4PrimaryParticle *g4p =
545  new G4PrimaryParticle(g_id, gp->momentum().px() * GeV, gp->momentum().py() * GeV, gp->momentum().pz() * GeV);
546  if (g4p->GetG4code() != nullptr) {
547  g4p->SetMass(g4p->GetG4code()->GetPDGMass());
548  g4p->SetCharge(g4p->GetG4code()->GetPDGCharge());
549  }
550  setGenId(g4p, i);
551  if (particlePassesPrimaryCuts(g4p->GetMomentum())) {
552  G4PrimaryVertex *v = new G4PrimaryVertex(gp->production_vertex()->position().x() * mm,
553  gp->production_vertex()->position().y() * mm,
554  gp->production_vertex()->position().z() * mm,
555  gp->production_vertex()->position().t() * mm / c_light);
556  v->SetPrimary(g4p);
557  g4evt->AddPrimaryVertex(v);
558  if (verbose > 0) {
559  v->Print();
560  }
561 
562  } else {
563  delete g4p;
564  }
565  }
566  } // end loop on HepMC particles
567 }
#define LogDebug(id)
double theDecRCut2
Definition: Generator.h:56
bool fPhiCuts
Definition: Generator.h:46
T getParameter(std::string const &) const
void setGenId(G4PrimaryParticle *p, int id) const
Definition: Generator.h:40
const double GeV
Definition: MathUtil.h:16
double theMinPtCut2
Definition: Generator.h:54
double theEtaCutForHector
Definition: Generator.h:57
double Z_lmax
Definition: Generator.h:64
void HepMC2G4(const HepMC::GenEvent *g, G4Event *e)
Definition: Generator.cc:102
double theMaxPhiCut
Definition: Generator.h:50
Generator(const edm::ParameterSet &p)
Definition: Generator.cc:24
bool fPtransCut
Definition: Generator.h:44
#define nullptr
bool exists(std::string const &parameterName) const
checks if a parameter exists
double theDecLenCut
Definition: Generator.h:58
double theMinEtaCut
Definition: Generator.h:51
bool isExoticNonDetectable(int pdgcode) const
Definition: Generator.cc:516
bool pdgFilterSel
Definition: Generator.h:66
bool fFinalState
Definition: Generator.h:48
void particleAssignDaughters(G4PrimaryParticle *p, HepMC::GenParticle *hp, double length)
Definition: Generator.cc:409
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double theMinPCut
Definition: Generator.h:53
LumiMonitorFilter * fLumiFilter
Definition: Generator.h:60
bool fFiductialCuts
Definition: Generator.h:47
T sqrt(T t)
Definition: SSEVec.h:18
bool fPDGFilter
Definition: Generator.h:67
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double weight_
Definition: Generator.h:63
double theMaxEtaCut
Definition: Generator.h:52
ii
Definition: cuy.py:590
double theMaxPCut
Definition: Generator.h:55
bool isExotic(int pdgcode) const
Definition: Generator.cc:506
void Describe() const
bool isGoodForLumiMonitor(const HepMC::GenParticle *) const
math::XYZTLorentzVector * vtx_
Definition: Generator.h:62
void nonBeamEvent2G4(const HepMC::GenEvent *g, G4Event *e)
Definition: Generator.cc:535
double Z_lmin
Definition: Generator.h:64
virtual ~Generator()
Definition: Generator.cc:100
HLT enums.
bool particlePassesPrimaryCuts(const G4ThreeVector &p) const
Definition: Generator.cc:474
std::vector< int > pdgFilter
Definition: Generator.h:65
bool IsInTheFilterList(int pdgcode) const
Definition: Generator.cc:525
bool fEtaCuts
Definition: Generator.h:45
double Z_hector
Definition: Generator.h:64
bool fPCuts
Definition: Generator.h:43