File indexing completed on 2024-09-07 04:36:36
0001
0002
0003
0004
0005
0006
0007 #include "HLTEgammaL1MatchFilterRegional.h"
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0010 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
0011 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/Framework/interface/EventSetup.h"
0016
0017 #define TWOPI 6.283185308
0018
0019
0020
0021 HLTEgammaL1MatchFilterRegional::HLTEgammaL1MatchFilterRegional(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0022 candIsolatedTag_ = iConfig.getParameter<edm::InputTag>("candIsolatedTag");
0023 l1IsolatedTag_ = iConfig.getParameter<edm::InputTag>("l1IsolatedTag");
0024 candNonIsolatedTag_ = iConfig.getParameter<edm::InputTag>("candNonIsolatedTag");
0025 l1NonIsolatedTag_ = iConfig.getParameter<edm::InputTag>("l1NonIsolatedTag");
0026 L1SeedFilterTag_ = iConfig.getParameter<edm::InputTag>("L1SeedFilterTag");
0027 l1CenJetsTag_ = iConfig.getParameter<edm::InputTag>("l1CenJetsTag");
0028 ncandcut_ = iConfig.getParameter<int>("ncandcut");
0029 doIsolated_ = iConfig.getParameter<bool>("doIsolated");
0030 region_eta_size_ = iConfig.getParameter<double>("region_eta_size");
0031 region_eta_size_ecap_ = iConfig.getParameter<double>("region_eta_size_ecap");
0032 region_phi_size_ = iConfig.getParameter<double>("region_phi_size");
0033 barrel_end_ = iConfig.getParameter<double>("barrel_end");
0034 endcap_end_ = iConfig.getParameter<double>("endcap_end");
0035
0036 candIsolatedToken_ = consumes<reco::RecoEcalCandidateCollection>(candIsolatedTag_);
0037 if (!doIsolated_)
0038 candNonIsolatedToken_ = consumes<reco::RecoEcalCandidateCollection>(candNonIsolatedTag_);
0039 L1SeedFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(L1SeedFilterTag_);
0040 }
0041
0042 HLTEgammaL1MatchFilterRegional::~HLTEgammaL1MatchFilterRegional() = default;
0043
0044 void HLTEgammaL1MatchFilterRegional::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0045 edm::ParameterSetDescription desc;
0046 makeHLTFilterDescription(desc);
0047 desc.add<edm::InputTag>("candIsolatedTag", edm::InputTag("hltRecoIsolatedEcalCandidate"));
0048 desc.add<edm::InputTag>("l1IsolatedTag", edm::InputTag("l1extraParticles", "Isolated"));
0049 desc.add<edm::InputTag>("candNonIsolatedTag", edm::InputTag("hltRecoNonIsolatedEcalCandidate"));
0050 desc.add<edm::InputTag>("l1NonIsolatedTag", edm::InputTag("l1extraParticles", "NonIsolated"));
0051 desc.add<edm::InputTag>("L1SeedFilterTag", edm::InputTag("theL1SeedFilter"));
0052 desc.add<edm::InputTag>("l1CenJetsTag", edm::InputTag("hltL1extraParticles", "Central"));
0053 desc.add<int>("ncandcut", 1);
0054 desc.add<bool>("doIsolated", true);
0055 desc.add<double>("region_eta_size", 0.522);
0056 desc.add<double>("region_eta_size_ecap", 1.0);
0057 desc.add<double>("region_phi_size", 1.044);
0058 desc.add<double>("barrel_end", 1.4791);
0059 desc.add<double>("endcap_end", 2.65);
0060 descriptions.add("hltEgammaL1MatchFilterRegional", desc);
0061 }
0062
0063
0064
0065
0066
0067 bool HLTEgammaL1MatchFilterRegional::hltFilter(edm::Event& iEvent,
0068 const edm::EventSetup& iSetup,
0069 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0070
0071 using namespace trigger;
0072 using namespace l1extra;
0073
0074 if (saveTags()) {
0075 filterproduct.addCollectionTag(l1IsolatedTag_);
0076 if (not doIsolated_)
0077 filterproduct.addCollectionTag(l1NonIsolatedTag_);
0078 filterproduct.addCollectionTag(l1CenJetsTag_);
0079 }
0080
0081 edm::Ref<reco::RecoEcalCandidateCollection> ref;
0082
0083
0084 int n(0);
0085
0086
0087 edm::Handle<reco::RecoEcalCandidateCollection> recoIsolecalcands;
0088 iEvent.getByToken(candIsolatedToken_, recoIsolecalcands);
0089
0090 edm::Handle<trigger::TriggerFilterObjectWithRefs> L1SeedOutput;
0091 iEvent.getByToken(L1SeedFilterToken_, L1SeedOutput);
0092
0093 std::vector<l1extra::L1EmParticleRef> l1EGIso;
0094 L1SeedOutput->getObjects(TriggerL1IsoEG, l1EGIso);
0095
0096 std::vector<l1extra::L1EmParticleRef> l1EGNonIso;
0097 L1SeedOutput->getObjects(TriggerL1NoIsoEG, l1EGNonIso);
0098
0099 std::vector<l1extra::L1JetParticleRef> l1Jets;
0100 L1SeedOutput->getObjects(TriggerL1CenJet, l1Jets);
0101
0102 for (auto recoecalcand = recoIsolecalcands->begin(); recoecalcand != recoIsolecalcands->end(); recoecalcand++) {
0103 if (fabs(recoecalcand->eta()) < endcap_end_) {
0104
0105
0106 bool matchedSCIso = matchedToL1Cand(l1EGIso, recoecalcand->eta(), recoecalcand->phi());
0107
0108
0109
0110
0111
0112 bool matchedSCNonIso = false;
0113 if (!doIsolated_) {
0114 matchedSCNonIso = matchedToL1Cand(l1EGNonIso, recoecalcand->eta(), recoecalcand->phi());
0115 }
0116
0117 bool matchedSCJet = matchedToL1Cand(l1Jets, recoecalcand->eta(), recoecalcand->phi());
0118
0119
0120 if (matchedSCIso || matchedSCNonIso || matchedSCJet) {
0121 n++;
0122
0123 ref = edm::Ref<reco::RecoEcalCandidateCollection>(recoIsolecalcands,
0124 distance(recoIsolecalcands->begin(), recoecalcand));
0125 filterproduct.addObject(TriggerCluster, ref);
0126 }
0127
0128 }
0129
0130 }
0131
0132
0133
0134
0135
0136 if (!doIsolated_ && !candNonIsolatedTag_.label().empty()) {
0137 edm::Handle<reco::RecoEcalCandidateCollection> recoNonIsolecalcands;
0138 iEvent.getByToken(candNonIsolatedToken_, recoNonIsolecalcands);
0139
0140 for (auto recoecalcand = recoNonIsolecalcands->begin(); recoecalcand != recoNonIsolecalcands->end();
0141 recoecalcand++) {
0142 if (fabs(recoecalcand->eta()) < endcap_end_) {
0143 bool matchedSCNonIso = matchedToL1Cand(l1EGNonIso, recoecalcand->eta(), recoecalcand->phi());
0144
0145 bool matchedSCJet = matchedToL1Cand(l1Jets, recoecalcand->eta(), recoecalcand->phi());
0146
0147 if (matchedSCNonIso || matchedSCJet) {
0148 n++;
0149 ref = edm::Ref<reco::RecoEcalCandidateCollection>(recoNonIsolecalcands,
0150 distance(recoNonIsolecalcands->begin(), recoecalcand));
0151 filterproduct.addObject(TriggerCluster, ref);
0152 }
0153
0154 }
0155
0156 }
0157 }
0158
0159
0160 bool accept(n >= ncandcut_);
0161
0162 return accept;
0163 }
0164
0165 bool HLTEgammaL1MatchFilterRegional::matchedToL1Cand(const std::vector<l1extra::L1EmParticleRef>& l1Cands,
0166 const float scEta,
0167 const float scPhi) const {
0168 for (auto const& l1Cand : l1Cands) {
0169
0170 double etaBinLow = 0.;
0171 double etaBinHigh = 0.;
0172 if (fabs(scEta) < barrel_end_) {
0173 etaBinLow = l1Cand->eta() - region_eta_size_ / 2.;
0174 etaBinHigh = etaBinLow + region_eta_size_;
0175 } else {
0176 etaBinLow = l1Cand->eta() - region_eta_size_ecap_ / 2.;
0177 etaBinHigh = etaBinLow + region_eta_size_ecap_;
0178 }
0179
0180 float deltaphi = fabs(scPhi - l1Cand->phi());
0181 if (deltaphi > TWOPI)
0182 deltaphi -= TWOPI;
0183 if (deltaphi > TWOPI / 2.)
0184 deltaphi = TWOPI - deltaphi;
0185
0186 if (scEta < etaBinHigh && scEta > etaBinLow && deltaphi < region_phi_size_ / 2.) {
0187 return true;
0188 }
0189 }
0190 return false;
0191 }
0192
0193 bool HLTEgammaL1MatchFilterRegional::matchedToL1Cand(const std::vector<l1extra::L1JetParticleRef>& l1Cands,
0194 const float scEta,
0195 const float scPhi) const {
0196 for (auto const& l1Cand : l1Cands) {
0197
0198 double etaBinLow = 0.;
0199 double etaBinHigh = 0.;
0200 if (fabs(scEta) < barrel_end_) {
0201 etaBinLow = l1Cand->eta() - region_eta_size_ / 2.;
0202 etaBinHigh = etaBinLow + region_eta_size_;
0203 } else {
0204 etaBinLow = l1Cand->eta() - region_eta_size_ecap_ / 2.;
0205 etaBinHigh = etaBinLow + region_eta_size_ecap_;
0206 }
0207
0208 float deltaphi = fabs(scPhi - l1Cand->phi());
0209 if (deltaphi > TWOPI)
0210 deltaphi -= TWOPI;
0211 if (deltaphi > TWOPI / 2.)
0212 deltaphi = TWOPI - deltaphi;
0213
0214 if (scEta < etaBinHigh && scEta > etaBinLow && deltaphi < region_phi_size_ / 2.) {
0215 return true;
0216 }
0217 }
0218 return false;
0219 }
0220
0221
0222 #include "FWCore/Framework/interface/MakerMacros.h"
0223 DEFINE_FWK_MODULE(HLTEgammaL1MatchFilterRegional);