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