CMS 3D CMS Logo

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