Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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