File indexing completed on 2024-04-06 12:18:34
0001
0002
0003
0004
0005
0006
0007
0008 #include "HLTrigger/JetMET/interface/HLTRapGapFilter.h"
0009
0010 #include "DataFormats/Common/interface/Handle.h"
0011
0012 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0013
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015
0016 #include "FWCore/Framework/interface/ESHandle.h"
0017 #include "FWCore/Framework/interface/EventSetup.h"
0018
0019 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0021
0022
0023
0024
0025 HLTRapGapFilter::HLTRapGapFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0026 inputTag_ = iConfig.getParameter<edm::InputTag>("inputTag");
0027 absEtaMin_ = iConfig.getParameter<double>("minEta");
0028 absEtaMax_ = iConfig.getParameter<double>("maxEta");
0029 caloThresh_ = iConfig.getParameter<double>("caloThresh");
0030 m_theJetToken = consumes<reco::CaloJetCollection>(inputTag_);
0031 }
0032
0033 HLTRapGapFilter::~HLTRapGapFilter() = default;
0034
0035 void HLTRapGapFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0036 edm::ParameterSetDescription desc;
0037 makeHLTFilterDescription(desc);
0038 desc.add<edm::InputTag>("inputJetTag", edm::InputTag("iterativeCone5CaloJets"));
0039 desc.add<double>("minEta", 3.0);
0040 desc.add<double>("maxEta", 5.0);
0041 desc.add<double>("caloThresh", 20.);
0042 descriptions.add("hltRapGapFilter", desc);
0043 }
0044
0045
0046 bool HLTRapGapFilter::hltFilter(edm::Event& iEvent,
0047 const edm::EventSetup& iSetup,
0048 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0049 using namespace reco;
0050 using namespace trigger;
0051
0052
0053 if (saveTags())
0054 filterproduct.addCollectionTag(inputTag_);
0055
0056 edm::Handle<CaloJetCollection> recocalojets;
0057 iEvent.getByToken(m_theJetToken, recocalojets);
0058
0059
0060 int n(0);
0061
0062
0063
0064 if (recocalojets->size() > 1) {
0065
0066
0067 double etjet = 0.;
0068 double etajet = 0.;
0069 double sumets = 0.;
0070
0071 for (auto const& recocalojet : *recocalojets) {
0072 etjet = recocalojet.energy();
0073 etajet = recocalojet.eta();
0074
0075 if (std::abs(etajet) > absEtaMin_ && std::abs(etajet) < absEtaMax_) {
0076 sumets += etjet;
0077
0078
0079 }
0080 }
0081
0082
0083 if (sumets <= caloThresh_) {
0084
0085 for (auto recocalojet = recocalojets->begin(); recocalojet != (recocalojets->end()); recocalojet++) {
0086 CaloJetRef ref(CaloJetRef(recocalojets, distance(recocalojets->begin(), recocalojet)));
0087 filterproduct.addObject(TriggerJet, ref);
0088 n++;
0089 }
0090 }
0091
0092 }
0093
0094
0095 bool accept(n > 0);
0096
0097 return accept;
0098 }