Skip to content

Commit 0c082e9

Browse files
authored
Merge pull request #47165 from 24LopezR/FixDisplacedSUSY
Fix in HepMC - G4 interface for displaced SUSY simulation
2 parents 7b5e37d + a952fa9 commit 0c082e9

File tree

3 files changed

+12
-2
lines changed

3 files changed

+12
-2
lines changed

SimG4Core/Application/python/g4SimHits_cfi.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -285,6 +285,7 @@
285285
MinPhiCut = cms.double(-3.14159265359), ## (radians)
286286
MaxPhiCut = cms.double(3.14159265359), ## according to CMS conventions
287287
ApplyLumiMonitorCuts = cms.bool(False), ## primary for lumi monitors
288+
IsSmuon = cms.bool(False),
288289
Verbosity = cms.untracked.int32(0),
289290
PDGselection = cms.PSet(
290291
PDGfilterSel = cms.bool(False), ## filter out unwanted particles

SimG4Core/Generators/interface/Generator.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -45,6 +45,7 @@ class Generator {
4545
bool fEtaCuts;
4646
bool fPhiCuts;
4747
bool fFiductialCuts;
48+
bool fSmuon;
4849
double theMinPhiCut;
4950
double theMaxPhiCut;
5051
double theMinEtaCut;

SimG4Core/Generators/src/Generator.cc

Lines changed: 10 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -25,6 +25,7 @@ Generator::Generator(const ParameterSet &p)
2525
fPtransCut(p.getParameter<bool>("ApplyPtransCut")),
2626
fEtaCuts(p.getParameter<bool>("ApplyEtaCuts")),
2727
fPhiCuts(p.getParameter<bool>("ApplyPhiCuts")),
28+
fSmuon(p.getParameter<bool>("IsSmuon")),
2829
theMinPhiCut(p.getParameter<double>("MinPhiCut")), // in radians (CMS standard)
2930
theMaxPhiCut(p.getParameter<double>("MaxPhiCut")),
3031
theMinEtaCut(p.getParameter<double>("MinEtaCut")),
@@ -365,7 +366,7 @@ void Generator::HepMC2G4(const HepMC::GenEvent *evt_orig, G4Event *g4evt) {
365366
// Decay chain outside the fiducial cylinder defined by theRDecLenCut
366367
// are used for Geant4 tracking with predefined decay channel
367368
// In the case of decay in vacuum particle is not tracked by Geant4
368-
} else if (2 == status && x2 * x2 + y2 * y2 >= theDecRCut2 && std::abs(z2) < Z_hector) {
369+
} else if (2 == status && x2 * x2 + y2 * y2 >= theDecRCut2 && (fSmuon || std::abs(z2) < Z_hector)) {
369370
toBeAdded = true;
370371
if (verbose > 1)
371372
edm::LogVerbatim("SimG4CoreGenerator") << "GenParticle barcode = " << (*pitr)->barcode() << " passed case 2"
@@ -476,7 +477,14 @@ void Generator::particleAssignDaughters(G4PrimaryParticle *g4p, HepMC::GenPartic
476477
LogDebug("SimG4CoreGenerator::::particleAssignDaughters")
477478
<< "Assigning a " << (*vpdec)->pdg_id() << " as daughter of a " << vp->pdg_id() << " status=" << status;
478479

479-
if ((status == 2 || (status == 23 && std::abs(vp->pdg_id()) == 1000015) || (status > 50 && status < 100)) &&
480+
bool isInList;
481+
if (fSmuon) {
482+
std::vector<int> fParticleList = {1000011, 1000013, 1000015, 2000011, 2000013, 2000015};
483+
isInList = (std::find(fParticleList.begin(), fParticleList.end(), std::abs(vp->pdg_id())) != fParticleList.end());
484+
} else {
485+
isInList = std::abs(vp->pdg_id()) == 1000015;
486+
}
487+
if ((status == 2 || (status == 23 && isInList) || (status > 50 && status < 100)) &&
480488
(*vpdec)->end_vertex() != nullptr) {
481489
double x2 = (*vpdec)->end_vertex()->position().x();
482490
double y2 = (*vpdec)->end_vertex()->position().y();

0 commit comments

Comments
 (0)