Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \ class HLTAcoFilter
0002  *
0003  *
0004  *  \author Dominique J. Mangeol
0005  *
0006 * acoplanar dphi(jet1,jet2),dphi(jet2,met),dphi(jet1,met)
0007 
0008 */
0009 
0010 #include <string>
0011 #include "HLTrigger/JetMET/interface/HLTAcoFilter.h"
0012 
0013 #include "DataFormats/Common/interface/Handle.h"
0014 
0015 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0016 
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 
0019 #include "FWCore/Framework/interface/ESHandle.h"
0020 #include "FWCore/Framework/interface/EventSetup.h"
0021 
0022 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0023 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0024 #include "FWCore/Utilities/interface/InputTag.h"
0025 
0026 //
0027 // constructors and destructor
0028 //
0029 HLTAcoFilter::HLTAcoFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0030   inputJetTag_ = iConfig.getParameter<edm::InputTag>("inputJetTag");
0031   inputMETTag_ = iConfig.getParameter<edm::InputTag>("inputMETTag");
0032   minDPhi_ = iConfig.getParameter<double>("minDeltaPhi");
0033   maxDPhi_ = iConfig.getParameter<double>("maxDeltaPhi");
0034   minEtjet1_ = iConfig.getParameter<double>("minEtJet1");
0035   minEtjet2_ = iConfig.getParameter<double>("minEtJet2");
0036   AcoString_ = iConfig.getParameter<std::string>("Acoplanar");
0037 
0038   m_theJetToken = consumes<reco::CaloJetCollection>(inputJetTag_);
0039   m_theMETToken = consumes<trigger::TriggerFilterObjectWithRefs>(inputMETTag_);
0040 }
0041 
0042 HLTAcoFilter::~HLTAcoFilter() = default;
0043 
0044 void HLTAcoFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0045   edm::ParameterSetDescription desc;
0046   makeHLTFilterDescription(desc);
0047   desc.add<edm::InputTag>("inputJetTag", edm::InputTag("IterativeCone5CaloJets"));
0048   desc.add<edm::InputTag>("inputMETTag", edm::InputTag("MET"));
0049   desc.add<double>("minDeltaPhi", 0.0);
0050   desc.add<double>("maxDeltaPhi", 2.0);
0051   desc.add<double>("minEtJet1", 20.0);
0052   desc.add<double>("minEtJet2", 20.0);
0053   desc.add<std::string>("Acoplanar", "Jet1Jet2");
0054   descriptions.add("hltAcoFilter", desc);
0055 }
0056 
0057 // ------------ method called to produce the data  ------------
0058 bool HLTAcoFilter::hltFilter(edm::Event& iEvent,
0059                              const edm::EventSetup& iSetup,
0060                              trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0061   using namespace std;
0062   using namespace edm;
0063   using namespace reco;
0064   using namespace trigger;
0065 
0066   // The filter object
0067   if (saveTags()) {
0068     filterproduct.addCollectionTag(inputJetTag_);
0069     filterproduct.addCollectionTag(inputMETTag_);
0070   }
0071 
0072   Handle<CaloJetCollection> recocalojets;
0073   iEvent.getByToken(m_theJetToken, recocalojets);
0074   Handle<trigger::TriggerFilterObjectWithRefs> metcal;
0075   iEvent.getByToken(m_theMETToken, metcal);
0076 
0077   // look at all candidates,  check cuts and add to filter object
0078   int n(0);
0079   int JetNum = recocalojets->size();
0080 
0081   // events with two or more jets
0082 
0083   double etjet1 = 0.;
0084   double etjet2 = 0.;
0085   double phijet1 = 0.;
0086   double phijet2 = 0.;
0087   //double etmiss=0.;
0088   double phimiss = 0.;
0089 
0090   VRcalomet vrefMET;
0091   metcal->getObjects(TriggerMET, vrefMET);
0092   CaloMETRef metRef = vrefMET.at(0);
0093   //etmiss  = vrefMET.at(0)->et();
0094   phimiss = vrefMET.at(0)->phi();
0095 
0096   CaloJetRef ref1, ref2;
0097 
0098   if (JetNum > 0) {
0099     auto recocalojet = recocalojets->begin();
0100 
0101     etjet1 = recocalojet->et();
0102     phijet1 = recocalojet->phi();
0103     ref1 = CaloJetRef(recocalojets, distance(recocalojets->begin(), recocalojet));
0104 
0105     if (JetNum > 1) {
0106       recocalojet++;
0107       etjet2 = recocalojet->et();
0108       phijet2 = recocalojet->phi();
0109       ref2 = CaloJetRef(recocalojets, distance(recocalojets->begin(), recocalojet));
0110     }
0111     double Dphi = -1.;
0112     int JetSel = 0;
0113 
0114     if (AcoString_ == "Jet2Met") {
0115       Dphi = std::abs(phimiss - phijet2);
0116       if (JetNum >= 2 && etjet1 > minEtjet1_ && etjet2 > minEtjet2_) {
0117         JetSel = 1;
0118       }
0119     }
0120     if (AcoString_ == "Jet1Jet2") {
0121       Dphi = std::abs(phijet1 - phijet2);
0122       if (JetNum >= 2 && etjet1 > minEtjet1_ && etjet2 > minEtjet2_) {
0123         JetSel = 1;
0124       }
0125     }
0126     if (AcoString_ == "Jet1Met") {
0127       Dphi = std::abs(phimiss - phijet1);
0128       if (JetNum >= 1 && etjet1 > minEtjet1_) {
0129         JetSel = 1;
0130       }
0131     }
0132 
0133     if (Dphi > M_PI) {
0134       Dphi = 2.0 * M_PI - Dphi;
0135     }
0136     if (JetSel > 0 && Dphi >= minDPhi_ && Dphi <= maxDPhi_) {
0137       if (AcoString_ == "Jet2Met" || AcoString_ == "Jet1Met") {
0138         filterproduct.addObject(TriggerMET, metRef);
0139       }
0140       if (AcoString_ == "Jet1Met" || AcoString_ == "Jet1Jet2") {
0141         filterproduct.addObject(TriggerJet, ref1);
0142       }
0143       if (AcoString_ == "Jet2Met" || AcoString_ == "Jet1Jet2") {
0144         filterproduct.addObject(TriggerJet, ref2);
0145       }
0146       n++;
0147     }
0148 
0149   }  // at least one jet
0150 
0151   // filter decision
0152   bool accept(n >= 1);
0153 
0154   return accept;
0155 }