Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Process modifier for displaced SUSY simulation fix #47197

Merged
merged 1 commit into from
Feb 3, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Original file line number Diff line number Diff line change
@@ -0,0 +1,4 @@
import FWCore.ParameterSet.Config as cms

# Designed to disable a bug affecting long lived slepton decays in HepMC-G4 interface
fixLongLivedSleptonSim = cms.Modifier()
10 changes: 9 additions & 1 deletion SimG4Core/Application/python/g4SimHits_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -285,7 +285,7 @@
MinPhiCut = cms.double(-3.14159265359), ## (radians)
MaxPhiCut = cms.double(3.14159265359), ## according to CMS conventions
ApplyLumiMonitorCuts = cms.bool(False), ## primary for lumi monitors
IsSmuon = cms.bool(False),
IsSlepton = cms.bool(False),
Verbosity = cms.untracked.int32(0),
PDGselection = cms.PSet(
PDGfilterSel = cms.bool(False), ## filter out unwanted particles
Expand Down Expand Up @@ -790,3 +790,11 @@
HGCSD = dict(
HitCollection = 2)
)

##
## Fix for long-lived slepton simulation
##
from Configuration.ProcessModifiers.fixLongLivedSleptonSim_cff import fixLongLivedSleptonSim
dd4hep.toModify( g4SimHits,
Generator = dict(IsSlepton = True)
)
2 changes: 1 addition & 1 deletion SimG4Core/Generators/interface/Generator.h
Original file line number Diff line number Diff line change
Expand Up @@ -45,7 +45,7 @@ class Generator {
bool fEtaCuts;
bool fPhiCuts;
bool fFiductialCuts;
bool fSmuon;
bool fSlepton;
double theMinPhiCut;
double theMaxPhiCut;
double theMinEtaCut;
Expand Down
1 change: 1 addition & 0 deletions SimG4Core/Generators/interface/Generator3.h
Original file line number Diff line number Diff line change
Expand Up @@ -45,6 +45,7 @@ class Generator3 {
bool fEtaCuts;
bool fPhiCuts;
bool fFiductialCuts;
bool fSlepton;
double theMinPhiCut;
double theMaxPhiCut;
double theMinEtaCut;
Expand Down
6 changes: 3 additions & 3 deletions SimG4Core/Generators/src/Generator.cc
Original file line number Diff line number Diff line change
Expand Up @@ -25,7 +25,7 @@ Generator::Generator(const ParameterSet &p)
fPtransCut(p.getParameter<bool>("ApplyPtransCut")),
fEtaCuts(p.getParameter<bool>("ApplyEtaCuts")),
fPhiCuts(p.getParameter<bool>("ApplyPhiCuts")),
fSmuon(p.getParameter<bool>("IsSmuon")),
fSlepton(p.getParameter<bool>("IsSlepton")),
theMinPhiCut(p.getParameter<double>("MinPhiCut")), // in radians (CMS standard)
theMaxPhiCut(p.getParameter<double>("MaxPhiCut")),
theMinEtaCut(p.getParameter<double>("MinEtaCut")),
Expand Down Expand Up @@ -366,7 +366,7 @@ void Generator::HepMC2G4(const HepMC::GenEvent *evt_orig, G4Event *g4evt) {
// Decay chain outside the fiducial cylinder defined by theRDecLenCut
// are used for Geant4 tracking with predefined decay channel
// In the case of decay in vacuum particle is not tracked by Geant4
} else if (2 == status && x2 * x2 + y2 * y2 >= theDecRCut2 && (fSmuon || std::abs(z2) < Z_hector)) {
} else if (2 == status && x2 * x2 + y2 * y2 >= theDecRCut2 && (fSlepton || std::abs(z2) < Z_hector)) {
toBeAdded = true;
if (verbose > 1)
edm::LogVerbatim("SimG4CoreGenerator") << "GenParticle barcode = " << (*pitr)->barcode() << " passed case 2"
Expand Down Expand Up @@ -478,7 +478,7 @@ void Generator::particleAssignDaughters(G4PrimaryParticle *g4p, HepMC::GenPartic
<< "Assigning a " << (*vpdec)->pdg_id() << " as daughter of a " << vp->pdg_id() << " status=" << status;

bool isInList;
if (fSmuon) {
if (fSlepton) {
std::vector<int> fParticleList = {1000011, 1000013, 1000015, 2000011, 2000013, 2000015};
isInList = (std::find(fParticleList.begin(), fParticleList.end(), std::abs(vp->pdg_id())) != fParticleList.end());
} else {
Expand Down
12 changes: 10 additions & 2 deletions SimG4Core/Generators/src/Generator3.cc
Original file line number Diff line number Diff line change
Expand Up @@ -27,6 +27,7 @@ Generator3::Generator3(const ParameterSet &p)
fPtransCut(p.getParameter<bool>("ApplyPtransCut")),
fEtaCuts(p.getParameter<bool>("ApplyEtaCuts")),
fPhiCuts(p.getParameter<bool>("ApplyPhiCuts")),
fSlepton(p.getParameter<bool>("IsSlepton")),
theMinPhiCut(p.getParameter<double>("MinPhiCut")), // in radians (CMS standard)
theMaxPhiCut(p.getParameter<double>("MaxPhiCut")),
theMinEtaCut(p.getParameter<double>("MinEtaCut")),
Expand Down Expand Up @@ -367,7 +368,7 @@ void Generator3::HepMC2G4(const HepMC3::GenEvent *evt_orig, G4Event *g4evt) {
// Decay chain outside the fiducial cylinder defined by theRDecLenCut
// are used for Geant4 tracking with predefined decay channel
// In the case of decay in vacuum particle is not tracked by Geant4
} else if (2 == status && x2 * x2 + y2 * y2 >= theDecRCut2 && std::abs(z2) < Z_hector) {
} else if (2 == status && x2 * x2 + y2 * y2 >= theDecRCut2 && (fSlepton || std::abs(z2) < Z_hector)) {
toBeAdded = true;
if (verbose > 1)
edm::LogVerbatim("SimG4CoreGenerator3") << "GenParticle barcode = " << pitr->id() << " passed case 2"
Expand Down Expand Up @@ -477,7 +478,14 @@ void Generator3::particleAssignDaughters(G4PrimaryParticle *g4p, HepMC3::GenPart
LogDebug("SimG4CoreGenerator3::::particleAssignDaughters")
<< "Assigning a " << vpdec->pid() << " as daughter of a " << vp->pid() << " status=" << status;

if ((status == 2 || (status == 23 && std::abs(vp->pid()) == 1000015) || (status > 50 && status < 100)) &&
bool isInList;
if (fSlepton) {
std::vector<int> fParticleList = {1000011, 1000013, 1000015, 2000011, 2000013, 2000015};
isInList = (std::find(fParticleList.begin(), fParticleList.end(), std::abs(vp->pdg_id())) != fParticleList.end());
} else {
isInList = std::abs(vp->pdg_id()) == 1000015;
}
if ((status == 2 || (status == 23 && isInList) || (status > 50 && status < 100)) &&
vpdec->end_vertex() != nullptr) {
double x2 = vpdec->end_vertex()->position().x();
double y2 = vpdec->end_vertex()->position().y();
Expand Down