File indexing completed on 2024-04-06 12:18:21
0001 #include "HLTElectronPFMTFilter.h"
0002 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0003
0004 template <typename T>
0005 HLTElectronPFMTFilter<T>::HLTElectronPFMTFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0006
0007 inputMetTag_ = iConfig.getParameter<edm::InputTag>("inputMetTag");
0008 minMht_ = iConfig.getParameter<double>("minMht");
0009
0010
0011 inputEleTag_ = iConfig.getParameter<edm::InputTag>("inputEleTag");
0012 lowerMTCut_ = iConfig.getParameter<double>("lowerMTCut");
0013 upperMTCut_ = iConfig.getParameter<double>("upperMTCut");
0014 minN_ = iConfig.getParameter<int>("minN");
0015 l1EGTag_ = iConfig.getParameter<edm::InputTag>("l1EGCand");
0016
0017 inputMetToken_ = consumes<reco::METCollection>(inputMetTag_);
0018 inputEleToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(inputEleTag_);
0019 }
0020
0021 template <typename T>
0022 HLTElectronPFMTFilter<T>::~HLTElectronPFMTFilter() = default;
0023
0024 template <typename T>
0025 void HLTElectronPFMTFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0026 edm::ParameterSetDescription desc;
0027 makeHLTFilterDescription(desc);
0028 desc.add<edm::InputTag>("inputMetTag", edm::InputTag("hltPFMHT"));
0029 desc.add<edm::InputTag>("inputEleTag", edm::InputTag("hltEle25CaloIdVTTrkIdTCaloIsoTTrkIsoTTrackIsolFilter"));
0030 desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
0031 desc.add<int>("minN", 0);
0032 desc.add<double>("minMht", 0.0);
0033 desc.add<double>("lowerMTCut", 0.0);
0034 desc.add<double>("upperMTCut", 9999.0);
0035 descriptions.add(defaultModuleLabel<HLTElectronPFMTFilter<T>>(), desc);
0036 }
0037
0038 template <typename T>
0039 bool HLTElectronPFMTFilter<T>::hltFilter(edm::Event& iEvent,
0040 const edm::EventSetup& iSetup,
0041 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0042 using namespace std;
0043 using namespace edm;
0044 using namespace reco;
0045 using namespace trigger;
0046
0047 if (saveTags()) {
0048 filterproduct.addCollectionTag(inputMetTag_);
0049 filterproduct.addCollectionTag(l1EGTag_);
0050 }
0051
0052
0053 edm::Handle<reco::METCollection> pfMHT;
0054 iEvent.getByToken(inputMetToken_, pfMHT);
0055
0056
0057 if (!pfMHT.isValid()) {
0058 edm::LogError("HLTElectronPFMTFilter") << "missing input Met collection!";
0059 }
0060
0061 const METCollection* metcol = pfMHT.product();
0062 const MET* met;
0063 met = &(metcol->front());
0064
0065 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0066 iEvent.getByToken(inputEleToken_, PrevFilterOutput);
0067
0068 int nW = 0;
0069
0070 vector<Ref<vector<T>>> refEleCollection;
0071
0072 PrevFilterOutput->getObjects(TriggerElectron, refEleCollection);
0073 int trigger_type = trigger::TriggerElectron;
0074 if (refEleCollection.empty()) {
0075 PrevFilterOutput->getObjects(TriggerCluster, refEleCollection);
0076 trigger_type = trigger::TriggerCluster;
0077 if (refEleCollection.empty()) {
0078 PrevFilterOutput->getObjects(TriggerPhoton, refEleCollection);
0079 trigger_type = trigger::TriggerPhoton;
0080 }
0081 }
0082
0083 TLorentzVector pMET(met->px(), met->py(), 0.0, sqrt(met->px() * met->px() + met->py() * met->py()));
0084
0085 for (unsigned int i = 0; i < refEleCollection.size(); i++) {
0086 TLorentzVector pThisEle(
0087 refEleCollection.at(i)->px(), refEleCollection.at(i)->py(), 0.0, refEleCollection.at(i)->et());
0088 TLorentzVector pTot = pMET + pThisEle;
0089 double mass = pTot.M();
0090
0091 if (mass >= lowerMTCut_ && mass <= upperMTCut_ && pMET.E() >= minMht_) {
0092 nW++;
0093 filterproduct.addObject(trigger_type, refEleCollection.at(i));
0094 }
0095 }
0096
0097
0098 const bool accept(nW >= minN_);
0099 return accept;
0100 }
0101
0102
0103 #include "FWCore/Framework/interface/MakerMacros.h"
0104
0105 using HLTGsfElectronPFMTFilter = HLTElectronPFMTFilter<reco::Electron>;
0106 DEFINE_FWK_MODULE(HLTGsfElectronPFMTFilter);
0107
0108 using HLTEcalCandidatePFMTFilter = HLTElectronPFMTFilter<reco::RecoEcalCandidate>;
0109 DEFINE_FWK_MODULE(HLTEcalCandidatePFMTFilter);