Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:05:54

0001 /** \class HLTJetTagWithMatching
0002  *
0003  *  This class is an HLTFilter (a spcialized EDFilter) implementing
0004  *  tagged multi-jet trigger for b and tau.
0005  *  It should be run after the normal multi-jet trigger.
0006  *  It is like HLTJetTag, but it loop on jets collection instead of jetTag collection.
0007  *
0008  */
0009 
0010 #include <vector>
0011 #include <string>
0012 
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "DataFormats/Common/interface/Handle.h"
0015 #include "DataFormats/Common/interface/RefToBase.h"
0016 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0017 #include "DataFormats/BTauReco/interface/JetTag.h"
0018 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0019 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0020 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0022 #include "DataFormats/Math/interface/deltaR.h"
0023 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0024 
0025 #include "HLTJetTagWithMatching.h"
0026 
0027 //
0028 // constructors and destructor
0029 //
0030 
0031 template <typename T>
0032 HLTJetTagWithMatching<T>::HLTJetTagWithMatching(const edm::ParameterSet& config)
0033     : HLTFilter(config),
0034       m_Jets(config.getParameter<edm::InputTag>("Jets")),
0035       m_JetTags(config.getParameter<edm::InputTag>("JetTags")),
0036       m_MinTag(config.getParameter<double>("MinTag")),
0037       m_MaxTag(config.getParameter<double>("MaxTag")),
0038       m_MinJets(config.getParameter<int>("MinJets")),
0039       m_TriggerType(config.getParameter<int>("TriggerType")),
0040       m_deltaR(config.getParameter<double>("deltaR")) {
0041   m_JetsToken = consumes<std::vector<T>>(m_Jets), m_JetTagsToken = consumes<reco::JetTagCollection>(m_JetTags),
0042 
0043   edm::LogInfo("") << " (HLTJetTagWithMatching) trigger cuts: " << std::endl
0044                    << "\ttype of        jets used: " << m_Jets.encode() << std::endl
0045                    << "\ttype of tagged jets used: " << m_JetTags.encode() << std::endl
0046                    << "\tmin/max tag value: [" << m_MinTag << ".." << m_MaxTag << "]" << std::endl
0047                    << "\tmin no. tagged jets: " << m_MinJets << "\tTriggerType: " << m_TriggerType << std::endl;
0048 }
0049 
0050 template <typename T>
0051 HLTJetTagWithMatching<T>::~HLTJetTagWithMatching() = default;
0052 
0053 template <typename T>
0054 float HLTJetTagWithMatching<T>::findCSV(const typename std::vector<T>::const_iterator& jet,
0055                                         const reco::JetTagCollection& jetTags,
0056                                         float minDr) {
0057   float tmpCSV = -20;
0058   for (auto jetb = jetTags.begin(); (jetb != jetTags.end()); ++jetb) {
0059     float tmpDr = reco::deltaR(*jet, *(jetb->first));
0060 
0061     if (tmpDr < minDr) {
0062       minDr = tmpDr;
0063       tmpCSV = jetb->second;
0064     }
0065   }
0066   return tmpCSV;
0067 }
0068 
0069 template <typename T>
0070 void HLTJetTagWithMatching<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0071   edm::ParameterSetDescription desc;
0072   makeHLTFilterDescription(desc);
0073   desc.add<edm::InputTag>("Jets", edm::InputTag("hltJetCollection"));
0074   desc.add<edm::InputTag>("JetTags", edm::InputTag("HLTJetTagWithMatchingCollection"));
0075   desc.add<double>("MinTag", 2.0);
0076   desc.add<double>("MaxTag", 999999.0);
0077   desc.add<int>("MinJets", 1);
0078   desc.add<int>("TriggerType", 0);
0079   desc.add<double>("deltaR", 0.1);
0080   descriptions.add(defaultModuleLabel<HLTJetTagWithMatching<T>>(), desc);
0081 }
0082 
0083 //
0084 // member functions
0085 //
0086 
0087 // ------------ method called to produce the data  ------------
0088 template <typename T>
0089 bool HLTJetTagWithMatching<T>::hltFilter(edm::Event& event,
0090                                          const edm::EventSetup& setup,
0091                                          trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0092   using namespace std;
0093   using namespace edm;
0094   using namespace reco;
0095 
0096   typedef vector<T> TCollection;
0097   typedef Ref<TCollection> TRef;
0098 
0099   edm::Handle<TCollection> h_Jets;
0100   event.getByToken(m_JetsToken, h_Jets);
0101   if (saveTags())
0102     filterproduct.addCollectionTag(m_Jets);
0103 
0104   edm::Handle<JetTagCollection> h_JetTags;
0105   event.getByToken(m_JetTagsToken, h_JetTags);
0106 
0107   // check if the product this one depends on is available
0108   auto const& handle = h_JetTags;
0109   auto const& dependent = handle->keyProduct();
0110   if (not dependent.isNull() and not dependent.hasCache()) {
0111     // only an empty AssociationVector can have a invalid dependent collection
0112     edm::StableProvenance const& dependent_provenance = event.getStableProvenance(dependent.id());
0113     if (dependent_provenance.branchDescription().dropped())
0114       // FIXME the error message should be made prettier
0115       throw edm::Exception(edm::errors::ProductNotFound)
0116           << "Product " << handle.provenance()->branchName() << " requires product "
0117           << dependent_provenance.branchName() << ", which has been dropped";
0118   }
0119 
0120   TRef jetRef;
0121 
0122   // Look at all jets in decreasing order of Et.
0123   int nJet = 0;
0124   int nTag = 0;
0125   for (typename TCollection::const_iterator jet = h_Jets->begin(); jet != h_Jets->end(); ++jet) {
0126     jetRef = TRef(h_Jets, nJet);
0127     LogTrace("") << "Jet " << nJet << " : Et = " << jet->et()
0128                  << " , tag value = " << findCSV(jet, *h_JetTags, m_deltaR);
0129     ++nJet;
0130     // Check if jet is tagged.
0131     if ((m_MinTag <= findCSV(jet, *h_JetTags, m_deltaR)) and (findCSV(jet, *h_JetTags, m_deltaR) <= m_MaxTag)) {
0132       ++nTag;
0133 
0134       // Store a reference to the jets which passed tagging cuts
0135       filterproduct.addObject(m_TriggerType, jetRef);
0136     }
0137   }
0138 
0139   // filter decision
0140   bool accept = (nTag >= m_MinJets);
0141 
0142   edm::LogInfo("") << " trigger accept ? = " << accept << " nTag/nJet = " << nTag << "/" << nJet << std::endl;
0143 
0144   return accept;
0145 }