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