9 #include "HepPDT/ParticleID.hh" 13 #include "G4HEPEvtParticle.hh" 15 #include "G4ParticleDefinition.hh" 16 #include "G4PhysicalConstants.hh" 17 #include "G4SystemOfUnits.hh" 18 #include "G4UnitsTable.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 fFinalState(p.getParameter<
bool>(
"FixOfFinalStateRadiation")),
30 theMinPhiCut(p.getParameter<double>(
"MinPhiCut")),
31 theMaxPhiCut(p.getParameter<double>(
"MaxPhiCut")),
32 theMinEtaCut(p.getParameter<double>(
"MinEtaCut")),
33 theMaxEtaCut(p.getParameter<double>(
"MaxEtaCut")),
34 theMinPCut(p.getParameter<double>(
"MinPCut")),
35 theMaxPCut(p.getParameter<double>(
"MaxPCut")),
36 theEtaCutForHector(p.getParameter<double>(
"EtaCutForHector")),
37 verbose(p.getUntrackedParameter<
int>(
"Verbosity", 0)),
52 double theRDecLenCut = p.
getParameter<
double>(
"RDecLenCut") * cm;
62 if (p.
exists(
"PDGselection")) {
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] <<
" ";
88 edm::LogVerbatim(
"SimG4CoreGenerator") <<
"SimG4Core/Generator: Rdecaycut= " << theRDecLenCut / cm
89 <<
" cm; Zdecaycut= " << theDecLenCut / cm <<
"Z_min= " <<
Z_lmin / cm
90 <<
" cm; Z_max= " <<
Z_lmax <<
" cm; Z_hector = " <<
Z_hector <<
" cm\n" 91 <<
" ApplyCuts: " << fFiductialCuts
94 <<
" LumiMonitorCut: " << lumi <<
" FinalState: " <<
fFinalState;
105 if (*(evt->vertices_begin()) ==
nullptr) {
106 std::stringstream ss;
107 ss <<
"SimG4Core/Generator: in event " << g4evt->GetEventID() <<
" Corrupted Event - GenEvent with no vertex \n";
111 if (!evt->weights().empty()) {
113 for (
unsigned int iw = 1; iw < evt->weights().size(); ++iw) {
115 if (evt->weights()[iw] <= 0)
121 if (
vtx_ !=
nullptr) {
125 (*(evt->vertices_begin()))->position().y(),
126 (*(evt->vertices_begin()))->position().z(),
127 (*(evt->vertices_begin()))->position().t());
132 unsigned int ng4vtx = 0;
134 for (HepMC::GenEvent::vertex_const_iterator vitr = evt->vertices_begin(); vitr != evt->vertices_end(); ++vitr) {
140 HepMC::GenVertex::particle_iterator pitr;
147 int status = (*pitr)->status();
148 int pdg = (*pitr)->pdg_id();
157 if(
std::abs(pdg) == 9900015) status = 3;
168 LogDebug(
"SimG4CoreGenerator") <<
"GenVertex barcode = " << (*vitr)->barcode()
169 <<
" selected for GenParticle barcode = " << (*pitr)->barcode();
175 else if (status == 2) {
176 if ((*pitr)->end_vertex() !=
nullptr) {
177 double xx = (*pitr)->end_vertex()->position().x();
178 double yy = (*pitr)->end_vertex()->position().y();
179 double r_dd = xx * xx + yy *
yy;
184 <<
"GenVertex barcode = " << (*vitr)->barcode()
185 <<
" selected for GenParticle barcode = " << (*pitr)->barcode() <<
" radius = " <<
std::sqrt(r_dd);
203 double x1 = (*vitr)->position().x() * mm;
204 double y1 = (*vitr)->position().y() * mm;
205 double z1 = (*vitr)->position().z() * mm;
206 double t1 = (*vitr)->position().t() * mm / c_light;
208 G4PrimaryVertex *g4vtx =
new G4PrimaryVertex(x1, y1, z1, t1);
211 int status = (*pitr)->status();
212 int pdg = (*pitr)->pdg_id();
213 bool hasDecayVertex = (
nullptr != (*pitr)->end_vertex()) ?
true :
false;
220 <<
" Skiped GenParticle barcode= " << (*pitr)->barcode() <<
" PDGid= " << pdg <<
" status= " << status
222 <<
" isInTheList: " << isInTheList <<
" hasDecayVertex: " << hasDecayVertex;
228 <<
"Generator: pdg= " << pdg <<
" status= " << status <<
" hasPreDefinedDecay: " << hasDecayVertex
230 <<
" isInTheList: " <<
IsInTheFilterList(pdg) <<
"\n (x,y,z,t): (" << x1 <<
"," << y1 <<
"," << z1
234 status = hasDecayVertex ? 2 : 1;
238 if (2 == status && !hasDecayVertex) {
239 edm::LogWarning(
"SimG4CoreGenerator: in event ") << g4evt->GetEventID() <<
" a particle " 241 <<
" has status=2 but has no decay vertex, so will be fully " 246 if(
std::abs(pdg) == 9900015) status = 3;
251 double decay_length = 0.0;
253 x2 = (*pitr)->end_vertex()->position().x();
254 y2 = (*pitr)->end_vertex()->position().y();
255 z2 = (*pitr)->end_vertex()->position().z();
256 decay_length =
std::sqrt((x2 - x1) * (x2 - x1) + (y2 - y1) * (y2 - y1) + (z2 - z1) * (z2 - z1));
261 double px = (*pitr)->momentum().px();
262 double py = (*pitr)->momentum().py();
263 double pz = (*pitr)->momentum().pz();
264 double ptot =
std::sqrt(px * px + py * py + pz * pz);
273 const double minTan = 1.e-20;
280 double del = (zimpact - z1) / pz;
284 double rimpact2 = ximpact * ximpact + yimpact * yimpact;
287 LogDebug(
"SimG4CoreGenerator") <<
"Processing GenParticle barcode= " << (*pitr)->barcode() <<
" pdg= " << pdg
288 <<
" status= " << (*pitr)->status() <<
" st= " << status
289 <<
" rimpact(cm)= " <<
std::sqrt(rimpact2) / cm
290 <<
" zimpact(cm)= " << zimpact / cm <<
" ptot(GeV)= " << ptot
291 <<
" pz(GeV)= " << pz;
298 LogDebug(
"SimG4CoreGenerator") <<
"GenParticle barcode = " << (*pitr)->barcode() <<
" passed case 3";
308 double phi =
p.phi();
322 if (
std::abs(pz) >= minTan * ptot) {
323 if ((zi >=
Z_lmax) & (pz < 0.0)) {
325 }
else if ((zi <=
Z_lmin) & (pz > 0.0)) {
334 double del = (zi - z1) / pz;
348 LogDebug(
"SimG4CoreGenerator") <<
"GenParticle barcode = " << (*pitr)->barcode() <<
" passed case 1";
356 LogDebug(
"SimG4CoreGenerator") <<
"GenParticle barcode = " << (*pitr)->barcode() <<
" passed case 2" 357 <<
" decay_length(cm)= " << decay_length / cm;
361 G4PrimaryParticle *g4prim =
new G4PrimaryParticle(pdg, px *
GeV, py * GeV, pz * GeV);
363 if (g4prim->GetG4code() !=
nullptr) {
364 g4prim->SetMass(g4prim->GetG4code()->GetPDGMass());
365 double charge = g4prim->GetG4code()->GetPDGCharge();
372 g4prim->SetCharge(charge);
378 setGenId(g4prim, (*pitr)->barcode());
385 g4vtx->SetPrimary(g4prim);
386 edm::LogVerbatim(
"SimG4CoreGenerator") <<
"Generator: new particle pdg= " << pdg <<
" Ptot(GeV/c)= " << ptot
387 <<
" Pt= " <<
std::sqrt(px * px + py * py);
393 g4evt->AddPrimaryVertex(g4vtx);
400 G4PrimaryVertex *g4vtx =
new G4PrimaryVertex(0.0, 0.0, 0.0, 0.0);
403 g4evt->AddPrimaryVertex(g4vtx);
411 LogDebug(
"SimG4CoreGenerator") <<
"Special case of long decay length \n" 412 <<
"Assign daughters with to mother with decaylength=" << decaylength / cm <<
" cm";
417 double proper_time = decaylength / (
p.Beta() *
p.Gamma() * c_light);
418 g4p->SetProperTime(proper_time);
421 LogDebug(
"SimG4CoreGenerator") <<
" px= " <<
p.px() <<
" py= " <<
p.py() <<
" pz= " <<
p.pz() <<
" e= " <<
p.e()
422 <<
" beta= " <<
p.Beta() <<
" gamma= " <<
p.Gamma()
423 <<
" Proper time= " << proper_time / ns <<
" ns";
428 double x1 = vp->end_vertex()->position().x();
429 double y1 = vp->end_vertex()->position().y();
430 double z1 = vp->end_vertex()->position().z();
432 for (HepMC::GenVertex::particle_iterator vpdec = vp->end_vertex()->particles_begin(
HepMC::children);
437 (*vpdec)->momentum().px(), (*vpdec)->momentum().py(), (*vpdec)->momentum().pz(), (*vpdec)->momentum().e());
440 G4PrimaryParticle *g4daught =
441 new G4PrimaryParticle((*vpdec)->pdg_id(), pdec.x() *
GeV, pdec.y() *
GeV, pdec.z() *
GeV);
443 if (g4daught->GetG4code() !=
nullptr) {
444 g4daught->SetMass(g4daught->GetG4code()->GetPDGMass());
445 g4daught->SetCharge(g4daught->GetG4code()->GetPDGCharge());
450 setGenId(g4daught, (*vpdec)->barcode());
453 LogDebug(
"SimG4CoreGenerator") <<
"Assigning a " << (*vpdec)->pdg_id() <<
" as daughter of a " << vp->pdg_id();
455 bool checkStatus =
fFinalState ? (*vpdec)->status() > 3
456 : ((*vpdec)->status() == 23 &&
std::abs(vp->pdg_id()) == 1000015);
457 if (((*vpdec)->status() == 2 || checkStatus) &&
458 (*vpdec)->end_vertex() !=
nullptr) {
459 double x2 = (*vpdec)->end_vertex()->position().x();
460 double y2 = (*vpdec)->end_vertex()->position().y();
461 double z2 = (*vpdec)->end_vertex()->position().z();
462 double dd =
std::sqrt((x1 - x2) * (x1 - x2) + (y1 - y2) * (y1 - y2) + (z1 - z2) * (z1 - z2));
465 (*vpdec)->set_status(1000 + (*vpdec)->status());
466 g4p->SetDaughter(g4daught);
476 double ptot = p.mag();
482 if (ptot < pz + 1.
e-10) {
486 double eta = 0.5 * G4Log((ptot + pz) / (ptot - pz));
493 double phi = p.phi();
500 LogDebug(
"SimG4CoreGenerator") <<
"Generator ptot(GeV)= " << ptot / GeV <<
" eta= " << p.eta()
501 <<
" phi= " << p.phi() <<
" Flag= " <<
flag;
508 return ((pdgid >= 1000000 && pdgid < 4000000 && pdgid != 3000022) ||
519 int charge = pid.threeCharge();
520 return ((charge == 0) && (pdgid >= 1000000 && pdgid < 1000040))
537 for (HepMC::GenEvent::particle_const_iterator it = evt->particles_begin(); it != evt->particles_end(); ++it) {
540 int g_status = gp->status();
543 int g_id = gp->pdg_id();
544 G4PrimaryParticle *g4p =
545 new G4PrimaryParticle(g_id, gp->momentum().px() *
GeV, gp->momentum().py() *
GeV, gp->momentum().pz() *
GeV);
546 if (g4p->GetG4code() !=
nullptr) {
547 g4p->SetMass(g4p->GetG4code()->GetPDGMass());
548 g4p->SetCharge(g4p->GetG4code()->GetPDGCharge());
552 G4PrimaryVertex *
v =
new G4PrimaryVertex(gp->production_vertex()->position().x() * mm,
553 gp->production_vertex()->position().y() * mm,
554 gp->production_vertex()->position().z() * mm,
555 gp->production_vertex()->position().t() * mm / c_light);
557 g4evt->AddPrimaryVertex(v);
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