Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:34

0001 /** \class HLTPhi2METFilter
0002  *
0003  *
0004  *  \author Dominique J. Mangeol
0005  *
0006  */
0007 
0008 #include "HLTrigger/JetMET/interface/HLTPhi2METFilter.h"
0009 
0010 #include "DataFormats/Common/interface/Handle.h"
0011 
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "DataFormats/METReco/interface/CaloMET.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 // constructors and destructor
0024 //
0025 HLTPhi2METFilter::HLTPhi2METFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0026   inputJetTag_ = iConfig.getParameter<edm::InputTag>("inputJetTag");
0027   inputMETTag_ = iConfig.getParameter<edm::InputTag>("inputMETTag");
0028   minDPhi_ = iConfig.getParameter<double>("minDeltaPhi");
0029   maxDPhi_ = iConfig.getParameter<double>("maxDeltaPhi");
0030   minEtjet1_ = iConfig.getParameter<double>("minEtJet1");
0031   minEtjet2_ = iConfig.getParameter<double>("minEtJet2");
0032   m_theJetToken = consumes<reco::CaloJetCollection>(inputJetTag_);
0033   m_theMETToken = consumes<trigger::TriggerFilterObjectWithRefs>(inputMETTag_);
0034 }
0035 
0036 HLTPhi2METFilter::~HLTPhi2METFilter() = default;
0037 
0038 void HLTPhi2METFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0039   edm::ParameterSetDescription desc;
0040   makeHLTFilterDescription(desc);
0041   desc.add<edm::InputTag>("inputJetTag", edm::InputTag("iterativeCone5CaloJets"));
0042   desc.add<edm::InputTag>("inputMETTag", edm::InputTag("hlt1MET60"));
0043   desc.add<double>("minDeltaPhi", 0.377);
0044   desc.add<double>("minEtJet2", 60.);
0045   descriptions.add("hltPhi2METFilter", desc);
0046 }
0047 
0048 // ------------ method called to produce the data  ------------
0049 bool HLTPhi2METFilter::hltFilter(edm::Event& iEvent,
0050                                  const edm::EventSetup& iSetup,
0051                                  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0052   using namespace std;
0053   using namespace edm;
0054   using namespace reco;
0055   using namespace trigger;
0056 
0057   // The filter object
0058   if (saveTags()) {
0059     filterproduct.addCollectionTag(inputJetTag_);
0060     filterproduct.addCollectionTag(inputMETTag_);
0061   }
0062 
0063   Handle<CaloJetCollection> recocalojets;
0064   iEvent.getByToken(m_theJetToken, recocalojets);
0065   Handle<trigger::TriggerFilterObjectWithRefs> metcal;
0066   iEvent.getByToken(m_theMETToken, metcal);
0067 
0068   // look at all candidates,  check cuts and add to filter object
0069   int n(0);
0070 
0071   VRcalomet vrefMET;
0072   metcal->getObjects(TriggerMET, vrefMET);
0073   CaloMETRef metRef = vrefMET.at(0);
0074   CaloJetRef ref1, ref2;
0075 
0076   if (recocalojets->size() > 1) {
0077     // events with two or more jets
0078 
0079     double etjet1 = 0.;
0080     double etjet2 = 0.;
0081     double phijet2 = 0.;
0082     // double etmiss  = vrefMET.at(0)->et();
0083     double phimiss = vrefMET.at(0)->phi();
0084     int countjets = 0;
0085 
0086     for (auto recocalojet = recocalojets->begin(); recocalojet <= (recocalojets->begin() + 1); recocalojet++) {
0087       if (countjets == 0) {
0088         etjet1 = recocalojet->et();
0089         ref1 = CaloJetRef(recocalojets, distance(recocalojets->begin(), recocalojet));
0090       }
0091       if (countjets == 1) {
0092         etjet2 = recocalojet->et();
0093         phijet2 = recocalojet->phi();
0094         ref2 = CaloJetRef(recocalojets, distance(recocalojets->begin(), recocalojet));
0095       }
0096       countjets++;
0097     }
0098     double Dphi = std::abs(phimiss - phijet2);
0099     if (Dphi > M_PI)
0100       Dphi = 2.0 * M_PI - Dphi;
0101     if (etjet1 > minEtjet1_ && etjet2 > minEtjet2_ && Dphi >= minDPhi_ && Dphi <= maxDPhi_) {
0102       filterproduct.addObject(TriggerMET, metRef);
0103       filterproduct.addObject(TriggerJet, ref1);
0104       filterproduct.addObject(TriggerJet, ref2);
0105       n++;
0106     }
0107 
0108   }  // events with two or more jets
0109 
0110   // filter decision
0111   bool accept(n >= 1);
0112 
0113   return accept;
0114 }