Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 #include "HLTJetPairDzMatchFilter.h"
0003 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "DataFormats/Common/interface/RefToBase.h"
0007 #include "DataFormats/JetReco/interface/CaloJet.h"
0008 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0009 #include "DataFormats/JetReco/interface/PFJet.h"
0010 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0012 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0013 #include "DataFormats/Math/interface/deltaPhi.h"
0014 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0015 
0016 template <typename T>
0017 HLTJetPairDzMatchFilter<T>::HLTJetPairDzMatchFilter(const edm::ParameterSet& conf) : HLTFilter(conf) {
0018   m_jetTag = conf.getParameter<edm::InputTag>("JetSrc");
0019   m_jetToken = consumes<std::vector<T>>(m_jetTag);
0020   m_jetMinPt = conf.getParameter<double>("JetMinPt");
0021   m_jetMaxEta = conf.getParameter<double>("JetMaxEta");
0022   m_jetMinDR = conf.getParameter<double>("JetMinDR");
0023   m_jetMaxDZ = conf.getParameter<double>("JetMaxDZ");
0024   m_triggerType = conf.getParameter<int>("TriggerType");
0025 
0026   // set the minimum DR between jets, so that one never has a chance
0027   // to create a pair out of the same Jet replicated with two
0028   // different vertices
0029   if (m_jetMinDR < 0.1)
0030     m_jetMinDR = 0.1;
0031 }
0032 
0033 template <typename T>
0034 void HLTJetPairDzMatchFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0035   edm::ParameterSetDescription desc;
0036   makeHLTFilterDescription(desc);
0037   desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltMatchL2Tau30ToPixelTrk5"));
0038   desc.add<double>("JetMinPt", 25.0);
0039   desc.add<double>("JetMaxEta", 2.4);
0040   desc.add<double>("JetMinDR", 0.2);
0041   desc.add<double>("JetMaxDZ", 0.2);
0042   desc.add<int>("TriggerType", 84);
0043   descriptions.add(defaultModuleLabel<HLTJetPairDzMatchFilter<T>>(), desc);
0044 }
0045 
0046 template <typename T>
0047 HLTJetPairDzMatchFilter<T>::~HLTJetPairDzMatchFilter() = default;
0048 
0049 template <typename T>
0050 bool HLTJetPairDzMatchFilter<T>::hltFilter(edm::Event& ev,
0051                                            const edm::EventSetup& es,
0052                                            trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0053   using namespace std;
0054   using namespace edm;
0055   using namespace reco;
0056 
0057   typedef vector<T> TCollection;
0058   typedef Ref<TCollection> TRef;
0059 
0060   // The resuilting filter object to store in the Event
0061 
0062   if (saveTags())
0063     filterproduct.addCollectionTag(m_jetTag);
0064 
0065   // Ref to Candidate object to be recorded in the filter object
0066   TRef ref;
0067 
0068   // *** Pick up L2 tau jets which have been equipped with some meaningful vertices before that ***
0069 
0070   edm::Handle<TCollection> jetsHandle;
0071   ev.getByToken(m_jetToken, jetsHandle);
0072   const TCollection& jets = *jetsHandle;
0073   const size_t n_jets = jets.size();
0074 
0075   // *** Combine jets into pairs and check the dz matching ***
0076 
0077   size_t npairs = 0;
0078   if (n_jets > 1)
0079     for (size_t j1 = 0; j1 < n_jets; ++j1) {
0080       if (jets[j1].pt() < m_jetMinPt || std::abs(jets[j1].eta()) > m_jetMaxEta)
0081         continue;
0082 
0083       float mindz = 99.f;
0084       for (size_t j2 = j1 + 1; j2 < n_jets; ++j2) {
0085         if (jets[j2].pt() < m_jetMinPt || std::abs(jets[j2].eta()) > m_jetMaxEta)
0086           continue;
0087 
0088         float deta = jets[j1].eta() - jets[j2].eta();
0089         float dphi = reco::deltaPhi(jets[j1].phi(), jets[j2].phi());
0090         float dr2 = dphi * dphi + deta * deta;
0091         float dz = jets[j1].vz() - jets[j2].vz();
0092 
0093         // skip pairs of jets that are close
0094         if (dr2 < m_jetMinDR * m_jetMinDR) {
0095           continue;
0096         }
0097 
0098         if (std::abs(dz) < std::abs(mindz))
0099           mindz = dz;
0100 
0101         // do not form a pair if dz is too large
0102         if (std::abs(dz) > m_jetMaxDZ) {
0103           continue;
0104         }
0105 
0106         // add references to both jets
0107         ref = TRef(jetsHandle, j1);
0108         filterproduct.addObject(m_triggerType, ref);
0109 
0110         ref = TRef(jetsHandle, j2);
0111         filterproduct.addObject(m_triggerType, ref);
0112 
0113         ++npairs;
0114       }
0115     }
0116 
0117   return npairs > 0;
0118 }