Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTFatJetMassFilter
0002 *
0003 *
0004 *  \author Maurizio Pierini
0005 *
0006 */
0007 
0008 #include <vector>
0009 
0010 #include "HLTrigger/JetMET/interface/HLTFatJetMassFilter.h"
0011 #include "DataFormats/Common/interface/Ref.h"
0012 #include "DataFormats/Common/interface/Handle.h"
0013 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0014 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0015 #include "FWCore/Framework/interface/ESHandle.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0021 #include "FWCore/Utilities/interface/InputTag.h"
0022 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0023 
0024 #include "DataFormats/Math/interface/deltaR.h"
0025 
0026 //
0027 // constructors and destructor
0028 //
0029 template <typename jetType>
0030 HLTFatJetMassFilter<jetType>::HLTFatJetMassFilter(const edm::ParameterSet& iConfig)
0031     : HLTFilter(iConfig),
0032       inputJetTag_(iConfig.template getParameter<edm::InputTag>("inputJetTag")),
0033       minMass_(iConfig.template getParameter<double>("minMass")),
0034       fatJetDeltaR_(iConfig.template getParameter<double>("fatJetDeltaR")),
0035       maxDeltaEta_(iConfig.template getParameter<double>("maxDeltaEta")),
0036       maxJetEta_(iConfig.template getParameter<double>("maxJetEta")),
0037       minJetPt_(iConfig.template getParameter<double>("minJetPt")),
0038       triggerType_(iConfig.template getParameter<int>("triggerType")) {
0039   m_theJetToken = consumes<std::vector<jetType>>(inputJetTag_);
0040   LogDebug("") << "HLTFatJetMassFilter: Input/minMass/fatJetDeltaR/maxDeltaEta/maxJetEta/minJetPt/triggerType : "
0041                << inputJetTag_.encode() << " " << minMass_ << " " << fatJetDeltaR_ << " " << maxDeltaEta_ << " "
0042                << maxJetEta_ << " " << minJetPt_ << " " << triggerType_;
0043 }
0044 
0045 template <typename jetType>
0046 HLTFatJetMassFilter<jetType>::~HLTFatJetMassFilter() = default;
0047 
0048 template <typename jetType>
0049 void HLTFatJetMassFilter<jetType>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0050   edm::ParameterSetDescription desc;
0051   makeHLTFilterDescription(desc);
0052   desc.add<edm::InputTag>("inputJetTag", edm::InputTag("hltCollection"));
0053   desc.add<double>("minMass", 0.0);
0054   desc.add<double>("fatJetDeltaR", 1.1);
0055   desc.add<double>("maxDeltaEta", 10.0);
0056   desc.add<double>("maxJetEta", 3.0);
0057   desc.add<double>("minJetPt", 30.0);
0058   desc.add<int>("triggerType", trigger::TriggerJet);
0059   descriptions.add(defaultModuleLabel<HLTFatJetMassFilter<jetType>>(), desc);
0060 }
0061 
0062 // ------------ method called to produce the data  ------------
0063 template <typename jetType>
0064 bool HLTFatJetMassFilter<jetType>::hltFilter(edm::Event& iEvent,
0065                                              const edm::EventSetup& iSetup,
0066                                              trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0067   using namespace std;
0068   using namespace edm;
0069   using namespace reco;
0070 
0071   typedef vector<jetType> JetCollection;
0072   typedef Ref<JetCollection> JetRef;
0073 
0074   // The filter object
0075   if (saveTags())
0076     filterproduct.addCollectionTag(inputJetTag_);
0077 
0078   // All jets
0079   Handle<JetCollection> objects;
0080   iEvent.getByToken(m_theJetToken, objects);
0081 
0082   // Selected jets
0083   CaloJetCollection recojets;
0084   typename JetCollection::const_iterator i(objects->begin());
0085   for (; i != objects->end(); i++) {
0086     if (std::abs(i->eta()) < maxJetEta_ && i->pt() >= minJetPt_) {
0087       reco::CaloJet jet(i->p4(), i->vertex(), reco::CaloJet::Specific());
0088       recojets.push_back(jet);
0089     }
0090   }
0091 
0092   // events with at least two jets
0093   if (recojets.size() < 2)
0094     return false;
0095 
0096   math::PtEtaPhiMLorentzVector j1(0.1, 0., 0., 0.);
0097   math::PtEtaPhiMLorentzVector j2(0.1, 0., 0., 0.);
0098   double jetPt1 = 0.;
0099   double jetPt2 = 0.;
0100 
0101   // look for the two highest-pT jet
0102   CaloJetCollection::const_iterator recojet(recojets.begin());
0103   for (; recojet != recojets.end(); recojet++) {
0104     if (recojet->pt() > jetPt1) {
0105       // downgrade the 1st jet to 2nd jet
0106       j2 = j1;
0107       jetPt2 = j2.pt();
0108       // promote this jet to 1st jet
0109       j1 = recojet->p4();
0110       jetPt1 = recojet->pt();
0111     } else if (recojet->pt() > jetPt2) {
0112       // promote this jet to 2nd jet
0113       j2 = recojet->p4();
0114       jetPt2 = recojet->pt();
0115     }
0116   }
0117 
0118   // apply DeltaEta cut
0119   double DeltaEta = std::abs(j1.eta() - j2.eta());
0120   if (DeltaEta > maxDeltaEta_)
0121     return false;
0122 
0123   math::PtEtaPhiMLorentzVector fj1;
0124   math::PtEtaPhiMLorentzVector fj2;
0125 
0126   // apply radiation recovery
0127   for (recojet = recojets.begin(); recojet != recojets.end(); recojet++) {
0128     double DeltaR1sq = reco::deltaR2(*recojet, j1);
0129     double DeltaR2sq = reco::deltaR2(*recojet, j2);
0130     if (DeltaR1sq < DeltaR2sq && DeltaR1sq < fatJetDeltaR_ * fatJetDeltaR_) {
0131       fj1 += recojet->p4();
0132     } else if (DeltaR2sq < fatJetDeltaR_ * fatJetDeltaR_) {
0133       fj2 += recojet->p4();
0134     }
0135   }
0136 
0137   // Apply mass cut
0138   fj1 += fj2;
0139   if (fj1.mass() < minMass_)
0140     return false;
0141 
0142   return true;
0143 }