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