File indexing completed on 2024-04-06 12:18:32
0001
0002
0003
0004
0005
0006
0007
0008 #include "HLTrigger/JetMET/interface/HLTDiJetAveFilter.h"
0009 #include "DataFormats/Common/interface/Ref.h"
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/Math/interface/deltaPhi.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "FWCore/Utilities/interface/InputTag.h"
0019 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0020
0021
0022
0023
0024 template <typename T>
0025 HLTDiJetAveFilter<T>::HLTDiJetAveFilter(const edm::ParameterSet& iConfig)
0026 : HLTFilter(iConfig),
0027 inputJetTag_(iConfig.template getParameter<edm::InputTag>("inputJetTag")),
0028 minPtAve_(iConfig.template getParameter<double>("minPtAve")),
0029 minPtJet3_(iConfig.template getParameter<double>("minPtJet3")),
0030 minDphi_(iConfig.template getParameter<double>("minDphi")),
0031 triggerType_(iConfig.template getParameter<int>("triggerType")) {
0032 m_theJetToken = consumes<std::vector<T>>(inputJetTag_);
0033 LogDebug("") << "HLTDiJetAveFilter: Input/minPtAve/minPtJet3/minDphi/triggerType : " << inputJetTag_.encode() << " "
0034 << minPtAve_ << " " << minPtJet3_ << " " << minDphi_ << " " << triggerType_;
0035 }
0036
0037 template <typename T>
0038 HLTDiJetAveFilter<T>::~HLTDiJetAveFilter() = default;
0039
0040 template <typename T>
0041 void HLTDiJetAveFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0042 edm::ParameterSetDescription desc;
0043 makeHLTFilterDescription(desc);
0044 desc.add<edm::InputTag>("inputJetTag", edm::InputTag("hltIterativeCone5CaloJets"));
0045 desc.add<double>("minPtAve", 100.0);
0046 desc.add<double>("minPtJet3", 99999.0);
0047 desc.add<double>("minDphi", -1.0);
0048 desc.add<int>("triggerType", trigger::TriggerJet);
0049 descriptions.add(defaultModuleLabel<HLTDiJetAveFilter<T>>(), desc);
0050 }
0051
0052
0053 template <typename T>
0054 bool HLTDiJetAveFilter<T>::hltFilter(edm::Event& iEvent,
0055 const edm::EventSetup& iSetup,
0056 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0057 using namespace std;
0058 using namespace edm;
0059 using namespace reco;
0060 using namespace trigger;
0061
0062 typedef vector<T> TCollection;
0063 typedef Ref<TCollection> TRef;
0064
0065
0066 if (saveTags())
0067 filterproduct.addCollectionTag(inputJetTag_);
0068
0069
0070 Handle<TCollection> objects;
0071 iEvent.getByToken(m_theJetToken, objects);
0072
0073
0074 int n(0);
0075
0076 if (objects->size() > 1) {
0077
0078
0079 double ptjet1 = 0., ptjet2 = 0., ptjet3 = 0.;
0080 double phijet1 = 0., phijet2 = 0;
0081 int countjets = 0;
0082
0083 int nmax = 1;
0084 if (objects->size() > 2)
0085 nmax = 2;
0086
0087 TRef JetRef1, JetRef2;
0088
0089 typename TCollection::const_iterator i(objects->begin());
0090 for (; i <= (objects->begin() + nmax); i++) {
0091 if (countjets == 0) {
0092 ptjet1 = i->pt();
0093 phijet1 = i->phi();
0094 JetRef1 = TRef(objects, distance(objects->begin(), i));
0095 }
0096 if (countjets == 1) {
0097 ptjet2 = i->pt();
0098 phijet2 = i->phi();
0099 JetRef2 = TRef(objects, distance(objects->begin(), i));
0100 }
0101 if (countjets == 2) {
0102 ptjet3 = i->pt();
0103 }
0104 ++countjets;
0105 }
0106
0107 double PtAve = (ptjet1 + ptjet2) / 2.;
0108 double Dphi = std::abs(deltaPhi(phijet1, phijet2));
0109
0110 if (PtAve > minPtAve_ && ptjet3 < minPtJet3_ && Dphi > minDphi_) {
0111 filterproduct.addObject(triggerType_, JetRef1);
0112 filterproduct.addObject(triggerType_, JetRef2);
0113 ++n;
0114 }
0115
0116 }
0117
0118
0119 bool accept(n >= 1);
0120
0121 return accept;
0122 }