File indexing completed on 2024-04-06 12:18:20
0001
0002
0003
0004
0005
0006
0007
0008 #include "HLTEgammaEtFilter.h"
0009
0010 #include "DataFormats/Common/interface/Handle.h"
0011
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015
0016 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0017
0018
0019
0020
0021 HLTEgammaEtFilter::HLTEgammaEtFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0022 inputTag_ = iConfig.getParameter<edm::InputTag>("inputTag");
0023 etcutEB_ = iConfig.getParameter<double>("etcutEB");
0024 etcutEE_ = iConfig.getParameter<double>("etcutEE");
0025 minEtaCut_ = iConfig.getParameter<double>("minEtaCut");
0026 maxEtaCut_ = iConfig.getParameter<double>("maxEtaCut");
0027 ncandcut_ = iConfig.getParameter<int>("ncandcut");
0028 l1EGTag_ = iConfig.getParameter<edm::InputTag>("l1EGCand");
0029 inputToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(inputTag_);
0030 }
0031
0032 void HLTEgammaEtFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0033 edm::ParameterSetDescription desc;
0034 makeHLTFilterDescription(desc);
0035 desc.add<edm::InputTag>("inputTag", edm::InputTag("HLTEgammaL1MatchFilter"));
0036 desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
0037 desc.add<double>("etcutEB", 1.0);
0038 desc.add<double>("etcutEE", 1.0);
0039 desc.add<double>("minEtaCut", -9999.0);
0040 desc.add<double>("maxEtaCut", 9999.0);
0041 desc.add<int>("ncandcut", 1);
0042 descriptions.add("hltEgammaEtFilter", desc);
0043 }
0044
0045 HLTEgammaEtFilter::~HLTEgammaEtFilter() = default;
0046
0047
0048 bool HLTEgammaEtFilter::hltFilter(edm::Event& iEvent,
0049 const edm::EventSetup& iSetup,
0050 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0051 using namespace trigger;
0052
0053
0054 if (saveTags()) {
0055 filterproduct.addCollectionTag(l1EGTag_);
0056 ;
0057 }
0058
0059
0060 edm::Ref<reco::RecoEcalCandidateCollection> ref;
0061
0062
0063
0064 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0065 iEvent.getByToken(inputToken_, PrevFilterOutput);
0066
0067 std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
0068 PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
0069 if (recoecalcands.empty())
0070 PrevFilterOutput->getObjects(TriggerPhoton,
0071 recoecalcands);
0072
0073
0074 int n(0);
0075
0076 for (auto& recoecalcand : recoecalcands) {
0077 ref = recoecalcand;
0078
0079 if ((ref->eta() < minEtaCut_) or (ref->eta() > maxEtaCut_))
0080 continue;
0081
0082 if ((fabs(ref->eta()) < 1.479 && ref->et() >= etcutEB_) || (fabs(ref->eta()) >= 1.479 && ref->et() >= etcutEE_)) {
0083 n++;
0084
0085 filterproduct.addObject(TriggerCluster, ref);
0086 }
0087 }
0088
0089
0090 bool accept(n >= ncandcut_);
0091
0092 return accept;
0093 }
0094
0095
0096 #include "FWCore/Framework/interface/MakerMacros.h"
0097 DEFINE_FWK_MODULE(HLTEgammaEtFilter);