9 #include "HepPDT/ParticleID.hh" 12 #include "G4HEPEvtParticle.hh" 14 #include "G4ParticleDefinition.hh" 15 #include "G4PhysicalConstants.hh" 16 #include "G4SystemOfUnits.hh" 17 #include "G4UnitsTable.hh" 24 : fPCuts(
p.getParameter<
bool>(
"ApplyPCuts")),
25 fPtransCut(
p.getParameter<
bool>(
"ApplyPtransCut")),
26 fEtaCuts(
p.getParameter<
bool>(
"ApplyEtaCuts")),
27 fPhiCuts(
p.getParameter<
bool>(
"ApplyPhiCuts")),
28 theMinPhiCut(
p.getParameter<double>(
"MinPhiCut")),
29 theMaxPhiCut(
p.getParameter<double>(
"MaxPhiCut")),
30 theMinEtaCut(
p.getParameter<double>(
"MinEtaCut")),
31 theMaxEtaCut(
p.getParameter<double>(
"MaxEtaCut")),
32 theMinPCut(
p.getParameter<double>(
"MinPCut")),
33 theMaxPCut(
p.getParameter<double>(
"MaxPCut")),
34 theEtaCutForHector(
p.getParameter<double>(
"EtaCutForHector")),
35 verbose(
p.getUntrackedParameter<
int>(
"Verbosity", 0)),
45 bool lumi =
p.getParameter<
bool>(
"ApplyLumiMonitorCuts");
50 double theRDecLenCut =
p.getParameter<
double>(
"RDecLenCut") * CLHEP::cm;
55 double theDecLenCut =
p.getParameter<
double>(
"LDecLenCut") * CLHEP::cm;
62 if (
p.exists(
"PDGselection")) {
68 ss <<
"SimG4Core/Generator: ";
70 ss <<
" Selecting only PDG ID = ";
72 ss <<
" Filtering out PDG ID = ";
88 edm::LogVerbatim(
"SimG4CoreGenerator") <<
"SimG4Core/Generator: Rdecaycut= " << theRDecLenCut / CLHEP::cm
90 <<
"Z_min= " <<
Z_lmin / CLHEP::cm <<
" cm; Z_max= " <<
Z_lmax / CLHEP::cm
94 <<
" Z_hector = " <<
Z_hector / CLHEP::cm <<
" cm\n" 98 <<
" LumiMonitorCut: " <<
lumi;
115 if (*(evt->vertices_begin()) ==
nullptr) {
116 std::stringstream
ss;
117 ss <<
"SimG4Core/Generator: in event " << g4evt->GetEventID() <<
" Corrupted Event - GenEvent with no vertex \n";
121 if (!evt->weights().empty()) {
123 for (
unsigned int iw = 1; iw < evt->weights().size(); ++iw) {
125 if (evt->weights()[iw] <= 0)
131 if (
vtx_ !=
nullptr) {
135 (*(evt->vertices_begin()))->position().y(),
136 (*(evt->vertices_begin()))->position().z(),
137 (*(evt->vertices_begin()))->position().t());
142 unsigned int ng4vtx = 0;
143 unsigned int ng4par = 0;
145 for (HepMC::GenEvent::vertex_const_iterator vitr = evt->vertices_begin(); vitr != evt->vertices_end(); ++vitr) {
151 HepMC::GenVertex::particle_iterator pitr;
157 int status = (*pitr)->status();
158 int pdg = (*pitr)->pdg_id();
180 LogDebug(
"SimG4CoreGenerator") <<
"GenVertex barcode = " << (*vitr)->barcode()
181 <<
" selected for GenParticle barcode = " << (*pitr)->barcode();
188 if ((*pitr)->end_vertex() !=
nullptr) {
189 double xx = (*pitr)->end_vertex()->position().x();
190 double yy = (*pitr)->end_vertex()->position().y();
196 <<
"GenVertex barcode = " << (*vitr)->barcode()
197 <<
" selected for GenParticle barcode = " << (*pitr)->barcode() <<
" radius = " <<
std::sqrt(r_dd);
215 double x1 = (*vitr)->position().x() * CLHEP::mm;
216 double y1 = (*vitr)->position().y() * CLHEP::mm;
217 double z1 = (*vitr)->position().z() * CLHEP::mm;
218 double t1 = (*vitr)->position().t() * CLHEP::mm / CLHEP::c_light;
220 G4PrimaryVertex *g4vtx =
new G4PrimaryVertex(
x1,
y1, z1,
t1);
223 int status = (*pitr)->status();
224 int pdg = (*pitr)->pdg_id();
225 bool hasDecayVertex = (
nullptr != (*pitr)->end_vertex());
233 <<
" Skiped GenParticle barcode= " << (*pitr)->barcode() <<
" PDGid= " <<
pdg <<
" status= " <<
status 235 <<
" isInTheList: " << isInTheList <<
" hasDecayVertex: " << hasDecayVertex;
242 <<
"Generator: pdg= " <<
pdg <<
" status= " <<
status <<
" hasPreDefinedDecay: " << hasDecayVertex
249 status = hasDecayVertex ? 2 : 1;
256 if (2 ==
status && !hasDecayVertex) {
258 << g4evt->GetEventID() <<
" a particle " 259 <<
" pdgid= " <<
pdg <<
" has status=2 but has no decay vertex, so will be fully tracked by Geant4";
266 double decay_length = 0.0;
268 x2 = (*pitr)->end_vertex()->position().x();
269 y2 = (*pitr)->end_vertex()->position().y();
270 z2 = (*pitr)->end_vertex()->position().z();
276 double px = (*pitr)->momentum().px();
277 double py = (*pitr)->momentum().py();
278 double pz = (*pitr)->momentum().pz();
288 const double minTan = 1.e-20;
291 double del = (zimpact - z1) / pz;
295 double rimpact2 = ximpact * ximpact + yimpact * yimpact;
298 LogDebug(
"SimG4CoreGenerator") <<
"Processing GenParticle barcode= " << (*pitr)->barcode() <<
" pdg= " <<
pdg 299 <<
" status= " << (*pitr)->status() <<
" st= " <<
status 300 <<
" rimpact(cm)= " <<
std::sqrt(rimpact2) / cm
301 <<
" zimpact(cm)= " << zimpact / cm <<
" ptot(GeV)= " << ptot
302 <<
" pz(GeV)= " << pz;
311 <<
"GenParticle barcode = " << (*pitr)->barcode() <<
" very forward; to be added: " << toBeAdded;
322 double phi =
p.phi();
336 if (
std::abs(pz) >= minTan * ptot) {
337 if ((zi >=
Z_lmax) & (pz < 0.0)) {
339 }
else if ((zi <=
Z_lmin) & (pz > 0.0)) {
348 double del = (zi - z1) / pz;
363 <<
"GenParticle barcode = " << (*pitr)->barcode() <<
" passed case 1";
371 edm::LogVerbatim(
"SimG4CoreGenerator") <<
"GenParticle barcode = " << (*pitr)->barcode() <<
" passed case 2" 372 <<
" decay_length(cm)= " << decay_length / CLHEP::cm;
376 G4PrimaryParticle *g4prim =
new G4PrimaryParticle(
pdg,
px * GeV,
py * GeV, pz * GeV);
378 if (g4prim->GetG4code() !=
nullptr) {
379 g4prim->SetMass(g4prim->GetG4code()->GetPDGMass());
380 double charge = g4prim->GetG4code()->GetPDGCharge();
387 g4prim->SetCharge(
charge);
393 setGenId(g4prim, (*pitr)->barcode());
402 g4vtx->SetPrimary(g4prim);
403 edm::LogVerbatim(
"SimG4CoreGenerator") <<
" " << ng4par <<
". new Geant4 particle pdg= " <<
pdg 405 <<
" status= " <<
status <<
"; dir= " << g4prim->GetMomentumDirection();
411 g4evt->AddPrimaryVertex(g4vtx);
418 G4PrimaryVertex *g4vtx =
new G4PrimaryVertex(0.0, 0.0, 0.0, 0.0);
422 g4evt->AddPrimaryVertex(g4vtx);
425 edm::LogVerbatim(
"SimG4CoreGenerator") <<
"The list of Geant4 primaries includes " << ng4par <<
" particles in " 426 << ng4vtx <<
" vertex";
432 LogDebug(
"SimG4CoreGenerator") <<
"Special case of long decay length \n" 433 <<
"Assign daughters with to mother with decaylength=" << decaylength / cm <<
" cm";
438 double proper_time = decaylength / (
p.Beta() *
p.Gamma() * c_light);
439 g4p->SetProperTime(proper_time);
442 LogDebug(
"SimG4CoreGenerator") <<
" px= " <<
p.px() <<
" py= " <<
p.py() <<
" pz= " <<
p.pz() <<
" e= " <<
p.e()
443 <<
" beta= " <<
p.Beta() <<
" gamma= " <<
p.Gamma()
444 <<
" Proper time= " << proper_time / ns <<
" ns";
449 double x1 = vp->end_vertex()->position().x();
450 double y1 = vp->end_vertex()->position().y();
451 double z1 = vp->end_vertex()->position().z();
453 for (HepMC::GenVertex::particle_iterator vpdec = vp->end_vertex()->particles_begin(
HepMC::children);
458 (*vpdec)->momentum().px(), (*vpdec)->momentum().py(), (*vpdec)->momentum().pz(), (*vpdec)->momentum().e());
461 G4PrimaryParticle *g4daught =
462 new G4PrimaryParticle((*vpdec)->pdg_id(), pdec.x() * CLHEP::GeV, pdec.y() * CLHEP::GeV, pdec.z() * CLHEP::GeV);
464 if (g4daught->GetG4code() !=
nullptr) {
465 g4daught->SetMass(g4daught->GetG4code()->GetPDGMass());
466 g4daught->SetCharge(g4daught->GetG4code()->GetPDGCharge());
471 setGenId(g4daught, (*vpdec)->barcode());
473 int status = (*vpdec)->status();
475 LogDebug(
"SimG4CoreGenerator::::particleAssignDaughters")
476 <<
"Assigning a " << (*vpdec)->pdg_id() <<
" as daughter of a " << vp->pdg_id() <<
" status=" <<
status;
479 (*vpdec)->end_vertex() !=
nullptr) {
480 double x2 = (*vpdec)->end_vertex()->position().x();
481 double y2 = (*vpdec)->end_vertex()->position().y();
482 double z2 = (*vpdec)->end_vertex()->position().z();
486 (*vpdec)->set_status(1000 +
status);
487 g4p->SetDaughter(g4daught);
497 double ptot =
p.mag();
498 if (
fPCuts && (ptot < theMinPCut * CLHEP::GeV || ptot >
theMaxPCut * CLHEP::GeV)) {
503 if (ptot < pz + 1.
e-10) {
507 double eta = 0.5 * G4Log((ptot + pz) / (ptot - pz));
514 double phi =
p.phi();
521 LogDebug(
"SimG4CoreGenerator") <<
"Generator ptot(GeV)= " << ptot / GeV <<
" eta= " <<
p.eta()
522 <<
" phi= " <<
p.phi() <<
" Flag= " <<
flag;
538 int charge = pid.threeCharge();
553 int i = g4evt->GetNumberOfPrimaryVertex();
554 for (HepMC::GenEvent::particle_const_iterator
it = evt->particles_begin();
it != evt->particles_end(); ++
it) {
559 if (
gp->status() != 1)
562 int pdg =
gp->pdg_id();
563 G4PrimaryParticle *g4p =
new G4PrimaryParticle(
564 pdg,
gp->momentum().px() * CLHEP::GeV,
gp->momentum().py() * CLHEP::GeV,
gp->momentum().pz() * CLHEP::GeV);
565 if (g4p->GetG4code() !=
nullptr) {
566 g4p->SetMass(g4p->GetG4code()->GetPDGMass());
567 g4p->SetCharge(g4p->GetG4code()->GetPDGCharge());
570 G4PrimaryVertex *
v =
new G4PrimaryVertex(
gp->production_vertex()->position().x() * CLHEP::mm,
571 gp->production_vertex()->position().y() * CLHEP::mm,
572 gp->production_vertex()->position().z() * CLHEP::mm,
573 gp->production_vertex()->position().t() * CLHEP::mm / CLHEP::c_light);
575 g4evt->AddPrimaryVertex(
v);
Log< level::Info, true > LogVerbatim
bool isExoticNonDetectable(int pdgcode) const
void nonCentralEvent2G4(const HepMC::GenEvent *g, G4Event *e)
double theEtaCutForHector
void HepMC2G4(const HepMC::GenEvent *g, G4Event *e)
Generator(const edm::ParameterSet &p)
bool isExotic(int pdgcode) const
bool isGoodForLumiMonitor(const HepMC::GenParticle *) const
void particleAssignDaughters(G4PrimaryParticle *p, HepMC::GenParticle *hp, double length)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
LumiMonitorFilter * fLumiFilter
bool particlePassesPrimaryCuts(const G4ThreeVector &p) const
Abs< T >::type abs(const T &t)
math::XYZTLorentzVector * vtx_
void setGenId(G4PrimaryParticle *p, int id) const
std::vector< int > pdgFilter
Log< level::Warning, false > LogWarning
bool IsInTheFilterList(int pdgcode) const