9 #include "HepPDT/ParticleID.hh" 13 #include "G4HEPEvtParticle.hh" 14 #include "G4ParticleDefinition.hh" 15 #include "G4UnitsTable.hh" 16 #include "G4SystemOfUnits.hh" 17 #include "G4PhysicalConstants.hh" 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")),
30 theMaxPhiCut(p.getParameter<double>(
"MaxPhiCut")),
31 theMinEtaCut(p.getParameter<double>(
"MinEtaCut")),
32 theMaxEtaCut(p.getParameter<double>(
"MaxEtaCut")),
33 theMinPCut(p.getParameter<double>(
"MinPCut")),
34 theMaxPCut(p.getParameter<double>(
"MaxPCut")),
35 theEtaCutForHector(p.getParameter<double>(
"EtaCutForHector")),
36 verbose(p.getUntrackedParameter<
int>(
"Verbosity",0)),
50 double theRDecLenCut = p.
getParameter<
double>(
"RDecLenCut")*cm;
60 if ( p.
exists(
"PDGselection") ) {
62 getParameter<bool>(
"PDGfilterSel");
64 getParameter<std::vector< int > >(
"PDGfilter");
65 if(!pdgFilter.empty()) {
68 ss <<
"SimG4Core/Generator: ";
70 ss <<
" Selecting only PDG ID = ";
72 ss <<
" Filtering out PDG ID = ";
74 for (
unsigned int ii = 0;
ii < pdgFilter.size(); ++
ii) {
75 ss << pdgFilter[
ii] <<
" ";
90 <<
"SimG4Core/Generator: Rdecaycut= " << theRDecLenCut/cm
91 <<
" cm; Zdecaycut= " << theDecLenCut/cm
93 <<
" cm; Z_hector = " <<
Z_hector <<
" cm\n" 94 <<
" ApplyCuts: " << fFiductialCuts
98 <<
" LumiMonitorCut: " <<
lumi;
110 HepMC::GenEvent *evt=
new HepMC::GenEvent(*evt_orig);
112 if ( *(evt->vertices_begin()) ==
nullptr ) {
113 std::stringstream ss;
114 ss <<
"SimG4Core/Generator: in event " << g4evt->GetEventID()
115 <<
" Corrupted Event - GenEvent with no vertex \n";
119 if (!evt->weights().empty()) {
122 for (
unsigned int iw=1; iw<evt->weights().size(); ++iw) {
125 if ( evt->weights()[iw] <= 0 )
break;
126 weight_ *= evt->weights()[iw] ;
130 if (
vtx_ !=
nullptr) {
delete vtx_; }
132 (*(evt->vertices_begin()))->position().y(),
133 (*(evt->vertices_begin()))->position().z(),
134 (*(evt->vertices_begin()))->position().t());
137 <<
"Generator: primary Vertex = (" 138 <<
vtx_->x() <<
", " <<
vtx_->y() <<
", " <<
vtx_->z() <<
")";
140 unsigned int ng4vtx = 0;
142 for(HepMC::GenEvent::vertex_const_iterator vitr= evt->vertices_begin();
143 vitr != evt->vertices_end(); ++vitr ) {
150 HepMC::GenVertex::particle_iterator pitr;
158 int status = (*pitr)->status();
159 int pdg = (*pitr)->pdg_id();
179 <<
"GenVertex barcode = " << (*vitr)->barcode()
180 <<
" selected for GenParticle barcode = " << (*pitr)->barcode();
186 else if (status == 2) {
187 if ( (*pitr)->end_vertex() !=
nullptr ) {
188 double xx = (*pitr)->end_vertex()->position().x();
189 double yy = (*pitr)->end_vertex()->position().y();
190 double r_dd = xx*xx+yy*
yy;
194 <<
"GenVertex barcode = " << (*vitr)->barcode()
195 <<
" selected for GenParticle barcode = " 196 << (*pitr)->barcode() <<
" radius = " <<
std::sqrt(r_dd);
214 double x1 = (*vitr)->position().x()*mm;
215 double y1 = (*vitr)->position().y()*mm;
216 double z1 = (*vitr)->position().z()*mm;
217 double t1 = (*vitr)->position().t()*mm/c_light;
219 G4PrimaryVertex* g4vtx =
new G4PrimaryVertex(x1, y1, z1, t1);
224 int status = (*pitr)->status();
225 int pdg = (*pitr)->pdg_id();
226 bool hasDecayVertex = (
nullptr != (*pitr)->end_vertex()) ?
true :
false;
233 <<
" Skiped GenParticle barcode= " << (*pitr)->barcode()
235 <<
" status= " << status
238 <<
" isInTheList: " << isInTheList
239 <<
" hasDecayVertex: " << hasDecayVertex;
245 <<
"Generator: pdg= " << pdg
246 <<
" status= " << status
247 <<
" hasPreDefinedDecay: " << hasDecayVertex
251 <<
"\n (x,y,z,t): (" << x1 <<
","<< y1 <<
"," << z1 <<
"," << t1 <<
")";
254 status = hasDecayVertex ? 2 : 1;
258 if (2 == status && !hasDecayVertex) {
260 << g4evt->GetEventID()
263 <<
" has status=2 but has no decay vertex, so will be fully tracked by Geant4";
270 double decay_length = 0.0;
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));
280 double px = (*pitr)->momentum().px();
281 double py = (*pitr)->momentum().py();
282 double pz = (*pitr)->momentum().pz();
283 double ptot =
std::sqrt(px*px + py*py + pz*pz);
292 const double minTan = 1.e-20;
294 if(pz > 0.0) { zimpact =
Z_hector; }
296 double del = (zimpact - z1)/pz;
300 double rimpact2 = ximpact*ximpact + yimpact*yimpact;
303 <<
"Processing GenParticle barcode= " << (*pitr)->barcode()
305 <<
" status= " << (*pitr)->status()
307 <<
" rimpact(cm)= " <<
std::sqrt(rimpact2)/cm
308 <<
" zimpact(cm)= " << zimpact/cm
309 <<
" ptot(GeV)= " << ptot
310 <<
" pz(GeV)= " << pz;
318 <<
"GenParticle barcode = " << (*pitr)->barcode()
332 double phi =
p.phi();
348 if((zi >=
Z_lmax) & (pz < 0.0)) {
350 }
else if((zi <=
Z_lmin) & (pz > 0.0)) {
353 if(pz > 0) { zi =
Z_lmax; }
356 double del = (zi - z1)/pz;
370 <<
"GenParticle barcode = " << (*pitr)->barcode()
376 }
else if(2 == status &&
381 <<
"GenParticle barcode = " << (*pitr)->barcode()
383 <<
" decay_length(cm)= " << decay_length/cm;
387 G4PrimaryParticle* g4prim=
388 new G4PrimaryParticle(pdg, px*
GeV, py*GeV, pz*GeV);
390 if ( g4prim->GetG4code() !=
nullptr ){
391 g4prim->SetMass( g4prim->GetG4code()->GetPDGMass() );
392 double charge = g4prim->GetG4code()->GetPDGCharge();
400 g4prim->SetCharge(charge);
406 setGenId( g4prim, (*pitr)->barcode() );
412 if (
verbose > 1 ) g4prim->Print();
413 g4vtx->SetPrimary(g4prim);
415 <<
"Generator: new particle pdg= " << pdg <<
" Ptot(GeV/c)= " << ptot
420 if (
verbose > 1 ) g4vtx->Print();
421 g4evt->AddPrimaryVertex(g4vtx);
428 G4PrimaryVertex* g4vtx =
new G4PrimaryVertex(0.0, 0.0, 0.0, 0.0);
429 if (
verbose > 1 ) g4vtx->Print();
430 g4evt->AddPrimaryVertex(g4vtx);
442 <<
"Special case of long decay length \n" 443 <<
"Assign daughters with to mother with decaylength=" 444 << decaylength/cm <<
" cm";
447 vp->momentum().pz(), vp->momentum().e());
450 double proper_time = decaylength/(
p.Beta()*
p.Gamma()*c_light);
451 g4p->SetProperTime(proper_time);
454 LogDebug(
"SimG4CoreGenerator") <<
" px= "<<
p.px()
458 <<
" beta= "<<
p.Beta()
459 <<
" gamma= " <<
p.Gamma()
460 <<
" Proper time= " <<proper_time/ns <<
" ns";
465 double x1 = vp->end_vertex()->position().x();
466 double y1 = vp->end_vertex()->position().y();
467 double z1 = vp->end_vertex()->position().z();
469 for (HepMC::GenVertex::particle_iterator
475 (*vpdec)->momentum().py(),
476 (*vpdec)->momentum().pz(),
477 (*vpdec)->momentum().e());
480 G4PrimaryParticle * g4daught=
481 new G4PrimaryParticle((*vpdec)->pdg_id(), pdec.x()*
GeV,
482 pdec.y()*
GeV, pdec.z()*
GeV);
484 if ( g4daught->GetG4code() !=
nullptr )
486 g4daught->SetMass( g4daught->GetG4code()->GetPDGMass() ) ;
487 g4daught->SetCharge( g4daught->GetG4code()->GetPDGCharge() ) ;
492 setGenId( g4daught, (*vpdec)->barcode() );
495 <<
"Assigning a "<< (*vpdec)->pdg_id()
496 <<
" as daughter of a " << vp->pdg_id();
497 if ( (*vpdec)->status() == 2 && (*vpdec)->end_vertex() !=
nullptr)
499 double x2 = (*vpdec)->end_vertex()->position().x();
500 double y2 = (*vpdec)->end_vertex()->position().y();
501 double z2 = (*vpdec)->end_vertex()->position().z();
502 double dd =
std::sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2));
505 (*vpdec)->set_status(1000+(*vpdec)->status());
506 g4p->SetDaughter(g4daught);
508 if (
verbose > 1 ) g4daught->Print();
516 double ptot = p.mag();
522 if(ptot < pz + 1.
e-10) {
526 double eta = 0.5*G4Log((ptot + pz)/(ptot - pz));
533 double phi = p.phi();
540 <<
"Generator ptot(GeV)= " << ptot/GeV <<
" eta= " << p.eta()
541 <<
" phi= " << p.phi() <<
" Flag= " <<
flag;
549 return ((pdgid >= 1000000 && pdgid < 4000000) ||
560 int charge = pid.threeCharge();
561 return ((charge==0) && (pdgid >= 1000000 && pdgid < 1000040))
575 for(HepMC::GenEvent::particle_const_iterator it = evt->particles_begin();
576 it != evt->particles_end(); ++it ) {
579 int g_status = gp->status();
582 int g_id = gp->pdg_id();
583 G4PrimaryParticle * g4p =
584 new G4PrimaryParticle(g_id,gp->momentum().px()*
GeV,
585 gp->momentum().py()*
GeV,
586 gp->momentum().pz()*
GeV);
587 if (g4p->GetG4code() !=
nullptr) {
588 g4p->SetMass(g4p->GetG4code()->GetPDGMass());
589 g4p->SetCharge(g4p->GetG4code()->GetPDGCharge());
593 G4PrimaryVertex *
v =
594 new G4PrimaryVertex(gp->production_vertex()->position().x()*mm,
595 gp->production_vertex()->position().y()*mm,
596 gp->production_vertex()->position().z()*mm,
597 gp->production_vertex()->position().t()*mm/c_light);
599 g4evt->AddPrimaryVertex(v);
600 if(
verbose > 0) { v->Print(); }
T getParameter(std::string const &) const
void setGenId(G4PrimaryParticle *p, int id) const
double theEtaCutForHector
void HepMC2G4(const HepMC::GenEvent *g, G4Event *e)
Generator(const edm::ParameterSet &p)
bool exists(std::string const ¶meterName) const
checks if a parameter exists
bool isExoticNonDetectable(int pdgcode) const
void particleAssignDaughters(G4PrimaryParticle *p, HepMC::GenParticle *hp, double length)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
LumiMonitorFilter * fLumiFilter
Abs< T >::type abs(const T &t)
bool isExotic(int pdgcode) const
bool isGoodForLumiMonitor(const HepMC::GenParticle *) const
math::XYZTLorentzVector * vtx_
void nonBeamEvent2G4(const HepMC::GenEvent *g, G4Event *e)
bool particlePassesPrimaryCuts(const G4ThreeVector &p) const
std::vector< int > pdgFilter
bool IsInTheFilterList(int pdgcode) const