Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-05-06 02:46:28

0001 /** \class HLTJetHbbFilter
0002  *
0003  * See header file for documentation
0004  *
0005  *  \author Ann Wang
0006  *
0007  */
0008 
0009 #include <vector>
0010 #include <string>
0011 
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "DataFormats/Math/interface/deltaPhi.h"
0018 #include "DataFormats/Math/interface/deltaR.h"
0019 #include "DataFormats/Common/interface/RefToBase.h"
0020 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0021 #include "DataFormats/JetReco/interface/CaloJet.h"
0022 #include "DataFormats/JetReco/interface/PFJet.h"
0023 #include "HLTrigger/JetMET/plugins/HLTJetHbbFilter.h"
0024 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0025 
0026 using namespace std;
0027 using namespace reco;
0028 using namespace trigger;
0029 //
0030 // constructors and destructor//
0031 //
0032 template <typename T>
0033 HLTJetHbbFilter<T>::HLTJetHbbFilter(const edm::ParameterSet& iConfig)
0034     : HLTFilter(iConfig),
0035       inputJets_(iConfig.getParameter<edm::InputTag>("inputJets")),
0036       inputJetTags_(iConfig.getParameter<edm::InputTag>("inputJetTags")),
0037       minmbb_(iConfig.getParameter<double>("minMbb")),
0038       maxmbb_(iConfig.getParameter<double>("maxMbb")),
0039       minptb1_(iConfig.getParameter<double>("minPtb1")),
0040       minptb2_(iConfig.getParameter<double>("minPtb2")),
0041       maxetab_(iConfig.getParameter<double>("maxEtab")),
0042       minptbb_(iConfig.getParameter<double>("minPtbb")),
0043       maxptbb_(iConfig.getParameter<double>("maxPtbb")),
0044       mintag1_(iConfig.getParameter<double>("minTag1")),
0045       mintag2_(iConfig.getParameter<double>("minTag2")),
0046       maxtag_(iConfig.getParameter<double>("maxTag")),
0047       triggerType_(iConfig.getParameter<int>("triggerType")) {
0048   m_theJetsToken = consumes<std::vector<T>>(inputJets_);
0049   m_theJetTagsToken = consumes<reco::JetTagCollection>(inputJetTags_);
0050 
0051   //put a dummy METCollection into the event, holding values for csv tag 1 and tag 2 values
0052   produces<reco::METCollection>();
0053 }
0054 
0055 template <typename T>
0056 HLTJetHbbFilter<T>::~HLTJetHbbFilter() {}
0057 
0058 template <typename T>
0059 void HLTJetHbbFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0060   edm::ParameterSetDescription desc;
0061   makeHLTFilterDescription(desc);
0062   desc.add<edm::InputTag>("inputJets", edm::InputTag("hltJetCollection"));
0063   desc.add<edm::InputTag>("inputJetTags", edm::InputTag(""));
0064   desc.add<double>("minMbb", 70);
0065   desc.add<double>("maxMbb", 200);
0066   desc.add<double>("minPtb1", -1);
0067   desc.add<double>("minPtb2", -1);
0068   desc.add<double>("maxEtab", 99999.0);
0069   desc.add<double>("minPtbb", -1);
0070   desc.add<double>("maxPtbb", -1);
0071   desc.add<double>("minTag1", 0.5);
0072   desc.add<double>("minTag2", 0.2);
0073   desc.add<double>("maxTag", 99999.0);
0074   desc.add<int>("triggerType", trigger::TriggerJet);
0075   descriptions.add(defaultModuleLabel<HLTJetHbbFilter<T>>(), desc);
0076 }
0077 
0078 template <typename T>
0079 float HLTJetHbbFilter<T>::findCSV(const typename std::vector<T>::const_iterator& jet,
0080                                   const reco::JetTagCollection& jetTags) {
0081   float minDr = 0.1;  //matching jet tag with jet
0082   float tmpCSV = -20;
0083   for (auto jetb = jetTags.begin(); (jetb != jetTags.end()); ++jetb) {
0084     float tmpDr = reco::deltaR(*jet, *(jetb->first));
0085     if (tmpDr < minDr) {
0086       minDr = tmpDr;
0087       tmpCSV = jetb->second;
0088     }
0089   }
0090   return tmpCSV;
0091 }
0092 //
0093 // member functions
0094 //
0095 
0096 // ------------ method called to produce the data  ------------
0097 template <typename T>
0098 bool HLTJetHbbFilter<T>::hltFilter(edm::Event& event,
0099                                    const edm::EventSetup& setup,
0100                                    trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0101   using namespace std;
0102   using namespace edm;
0103   using namespace reco;
0104   using namespace trigger;
0105 
0106   typedef vector<T> TCollection;
0107   typedef Ref<TCollection> TRef;
0108 
0109   bool accept(false);
0110   //const unsigned int nMax(15);
0111 
0112   Handle<TCollection> jets;
0113   event.getByToken(m_theJetsToken, jets);
0114   Handle<JetTagCollection> jetTags;
0115 
0116   event.getByToken(m_theJetTagsToken, jetTags);
0117 
0118   double tag1 = -99.;
0119   double tag2 = -99.;
0120 
0121   if (jetTags->size() < 2)
0122     return false;
0123 
0124   double ejet1 = -99.;
0125   double pxjet1 = -99.;
0126   double pyjet1 = -99.;
0127   double pzjet1 = -99.;
0128   double ptjet1 = -99.;
0129   double etajet1 = -99.;
0130 
0131   double ejet2 = -99.;
0132   double pxjet2 = -99.;
0133   double pyjet2 = -99.;
0134   double pzjet2 = -99.;
0135   double ptjet2 = -99.;
0136   double etajet2 = -99.;
0137 
0138   //looping through sets of jets
0139   for (auto jet1 = jets->begin(); (jet1 != jets->end()); ++jet1) {
0140     tag1 = findCSV(jet1, *jetTags);
0141     for (auto jet2 = (jet1 + 1); (jet2 != jets->end()); ++jet2) {
0142       tag2 = findCSV(jet2, *jetTags);
0143 
0144       ejet1 = jet1->energy();
0145       pxjet1 = jet1->px();
0146       pyjet1 = jet1->py();
0147       pzjet1 = jet1->pz();
0148       ptjet1 = jet1->pt();
0149       etajet1 = jet1->eta();
0150 
0151       ejet2 = jet2->energy();
0152       pxjet2 = jet2->px();
0153       pyjet2 = jet2->py();
0154       pzjet2 = jet2->pz();
0155       ptjet2 = jet2->pt();
0156       etajet2 = jet2->eta();
0157 
0158       if (((mintag1_ <= tag1) and (tag1 <= maxtag_)) &&
0159           ((mintag2_ <= tag2) and (tag2 <= maxtag_))) {                // if they're both b's
0160         if (fabs(etajet1) <= maxetab_ && fabs(etajet2) <= maxetab_) {  // if they satisfy the eta requirement
0161           if ((ptjet1 >= minptb1_ && ptjet2 >= minptb2_) ||
0162               (ptjet2 >= minptb1_ && ptjet1 >= minptb2_)) {  // if they satisfy the pt requirement
0163 
0164             double ptbb = sqrt((pxjet1 + pxjet2) * (pxjet1 + pxjet2) +
0165                                (pyjet1 + pyjet2) * (pyjet1 + pyjet2));  // pt of the two jets
0166 
0167             if (ptbb >= minptbb_ && (maxptbb_ < 0 || ptbb <= maxptbb_)) {  //if they satisfy the vector pt requirement
0168 
0169               double mbb = sqrt((ejet1 + ejet2) * (ejet1 + ejet2) - (pxjet1 + pxjet2) * (pxjet1 + pxjet2) -
0170                                 (pyjet1 + pyjet2) * (pyjet1 + pyjet2) -
0171                                 (pzjet1 + pzjet2) * (pzjet1 + pzjet2));  // mass of two jets
0172 
0173               if ((minmbb_ <= mbb) and (mbb <= maxmbb_)) {  // if they fit the mass requirement
0174                 accept = true;
0175 
0176                 TRef ref1 = TRef(jets, distance(jets->begin(), jet1));
0177                 TRef ref2 = TRef(jets, distance(jets->begin(), jet2));
0178 
0179                 if (saveTags())
0180                   filterproduct.addCollectionTag(inputJets_);
0181                 filterproduct.addObject(triggerType_, ref1);
0182                 filterproduct.addObject(triggerType_, ref2);
0183 
0184                 //create METCollection for storing csv tag1 and tag2 results
0185                 std::unique_ptr<reco::METCollection> csvObject(new reco::METCollection());
0186                 reco::MET::LorentzVector csvP4(tag1, tag2, 0, 0);
0187                 reco::MET::Point vtx(0, 0, 0);
0188                 reco::MET csvTags(csvP4, vtx);
0189                 csvObject->push_back(csvTags);
0190                 edm::RefProd<reco::METCollection> ref_before_put = event.getRefBeforePut<reco::METCollection>();
0191                 //put the METCollection into the event (necessary because of how addCollectionTag works...)
0192                 event.put(std::move(csvObject));
0193                 edm::Ref<reco::METCollection> csvRef(ref_before_put, 0);
0194                 if (saveTags())
0195                   filterproduct.addCollectionTag(edm::InputTag(*moduleLabel()));
0196                 filterproduct.addObject(trigger::TriggerMET, csvRef);  //give it the ID of a MET object
0197                 return accept;
0198               }
0199             }
0200           }
0201         }
0202       }
0203     }
0204   }
0205   return accept;
0206 }
0207 typedef HLTJetHbbFilter<CaloJet> HLTCaloJetHbbFilter;
0208 typedef HLTJetHbbFilter<PFJet> HLTPFJetHbbFilter;
0209 
0210 DEFINE_FWK_MODULE(HLTCaloJetHbbFilter);
0211 DEFINE_FWK_MODULE(HLTPFJetHbbFilter);