Skip to content

Commit 2ad1972

Browse files
committed
Adding optional deltaR cuts to HLTPMMassFilter
1 parent 344d809 commit 2ad1972

File tree

2 files changed

+9
-1
lines changed

2 files changed

+9
-1
lines changed

HLTrigger/Egamma/plugins/HLTPMMassFilter.cc

Lines changed: 7 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -20,6 +20,8 @@ HLTPMMassFilter::HLTPMMassFilter(const edm::ParameterSet& iConfig) : HLTFilter(i
2020

2121
lowerMassCut_ = iConfig.getParameter<double>("lowerMassCut");
2222
upperMassCut_ = iConfig.getParameter<double>("upperMassCut");
23+
lowerdRCut_ = iConfig.getParameter<double>("lowerdRCut");
24+
upperdRCut_ = iConfig.getParameter<double>("upperdRCut");
2325
nZcandcut_ = iConfig.getParameter<int>("nZcandcut");
2426
reqOppCharge_ = iConfig.getUntrackedParameter<bool>("reqOppCharge", false);
2527
isElectron1_ = iConfig.getUntrackedParameter<bool>("isElectron1", true);
@@ -38,6 +40,8 @@ void HLTPMMassFilter::fillDescriptions(edm::ConfigurationDescriptions& descripti
3840
desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOfflineBeamSpot"));
3941
desc.add<double>("lowerMassCut", 8.0);
4042
desc.add<double>("upperMassCut", 11.0);
43+
desc.add<double>("lowerdRCut", -1.0);
44+
desc.add<double>("upperdRCut", 9999.0);
4145
desc.add<int>("nZcandcut", 1);
4246
desc.addUntracked<bool>("reqOppCharge", true);
4347
desc.addUntracked<bool>("isElectron1", false);
@@ -157,7 +161,9 @@ bool HLTPMMassFilter::isGoodPair(TLorentzVector const& v1, TLorentzVector const&
157161
return false;
158162

159163
auto const mass = (v1 + v2).M();
160-
return (mass >= lowerMassCut_ and mass <= upperMassCut_);
164+
auto const dr = v1.DeltaR(v2);
165+
166+
return (mass >= lowerMassCut_ and mass <= upperMassCut_ and dr >= lowerdRCut_ and dr <= upperdRCut_);
161167
}
162168

163169
TLorentzVector HLTPMMassFilter::approxMomAtVtx(const MagneticField& magField,

HLTrigger/Egamma/plugins/HLTPMMassFilter.h

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -52,6 +52,8 @@ class HLTPMMassFilter : public HLTFilter {
5252
edm::EDGetTokenT<reco::BeamSpot> beamSpotToken_;
5353
double lowerMassCut_;
5454
double upperMassCut_;
55+
double lowerdRCut_;
56+
double upperdRCut_;
5557
int nZcandcut_; // number of Z candidates required
5658
bool reqOppCharge_;
5759

0 commit comments

Comments
 (0)