Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTForwardBackwardJetsFilter
0002  *
0003  *
0004  *
0005  */
0006 
0007 #include "HLTrigger/JetMET/interface/HLTForwardBackwardJetsFilter.h"
0008 
0009 #include "DataFormats/Common/interface/Ref.h"
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0012 
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 #include "FWCore/Utilities/interface/InputTag.h"
0020 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0021 
0022 //
0023 // constructors and destructor
0024 //
0025 template <typename T>
0026 HLTForwardBackwardJetsFilter<T>::HLTForwardBackwardJetsFilter(const edm::ParameterSet& iConfig)
0027     : HLTFilter(iConfig),
0028       inputTag_(iConfig.template getParameter<edm::InputTag>("inputTag")),
0029       minPt_(iConfig.template getParameter<double>("minPt")),
0030       minEta_(iConfig.template getParameter<double>("minEta")),
0031       maxEta_(iConfig.template getParameter<double>("maxEta")),
0032       nNeg_(iConfig.template getParameter<unsigned int>("nNeg")),
0033       nPos_(iConfig.template getParameter<unsigned int>("nPos")),
0034       nTot_(iConfig.template getParameter<unsigned int>("nTot")),
0035       triggerType_(iConfig.template getParameter<int>("triggerType")) {
0036   m_theJetToken = consumes<std::vector<T>>(inputTag_);
0037   LogDebug("") << "HLTForwardBackwardJetsFilter: Input/minPt/minEta/maxEta/triggerType : " << inputTag_.encode() << " "
0038                << minPt_ << " " << minEta_ << " " << maxEta_ << " " << nNeg_ << " " << nPos_ << " " << nTot_ << " "
0039                << triggerType_;
0040 }
0041 
0042 template <typename T>
0043 HLTForwardBackwardJetsFilter<T>::~HLTForwardBackwardJetsFilter() = default;
0044 
0045 template <typename T>
0046 void HLTForwardBackwardJetsFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0047   edm::ParameterSetDescription desc;
0048   makeHLTFilterDescription(desc);
0049   desc.add<edm::InputTag>("inputTag", edm::InputTag("hltIterativeCone5CaloJetsRegional"));
0050   desc.add<double>("minPt", 15.0);
0051   desc.add<double>("minEta", 3.0);
0052   desc.add<double>("maxEta", 5.1);
0053   desc.add<unsigned int>("nNeg", 1);
0054   desc.add<unsigned int>("nPos", 1);
0055   desc.add<unsigned int>("nTot", 0);
0056   desc.add<int>("triggerType", trigger::TriggerJet);
0057   descriptions.add(defaultModuleLabel<HLTForwardBackwardJetsFilter<T>>(), desc);
0058 }
0059 
0060 // ------------ method called to produce the data  ------------
0061 template <typename T>
0062 bool HLTForwardBackwardJetsFilter<T>::hltFilter(edm::Event& iEvent,
0063                                                 const edm::EventSetup& iSetup,
0064                                                 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0065   using namespace std;
0066   using namespace edm;
0067   using namespace reco;
0068   using namespace trigger;
0069 
0070   typedef vector<T> TCollection;
0071   typedef Ref<TCollection> TRef;
0072 
0073   // The filter object
0074   if (saveTags())
0075     filterproduct.addCollectionTag(inputTag_);
0076 
0077   // get hold of collection of objects
0078   Handle<TCollection> objects;
0079   iEvent.getByToken(m_theJetToken, objects);
0080 
0081   // look at all candidates,  check cuts and add to filter object
0082   unsigned int nPosJets(0);
0083   unsigned int nNegJets(0);
0084 
0085   typename TCollection::const_iterator jet;
0086   // look for jets satifying pt and eta cuts; first on the plus side, then the minus side
0087 
0088   for (jet = objects->begin(); jet != objects->end(); jet++) {
0089     double ptjet = jet->pt();
0090     double etajet = jet->eta();
0091     if (ptjet >= minPt_) {
0092       if ((minEta_ <= etajet) && (etajet <= maxEta_)) {
0093         nPosJets++;
0094         TRef ref = TRef(objects, distance(objects->begin(), jet));
0095         filterproduct.addObject(triggerType_, ref);
0096       }
0097       if ((-maxEta_ <= etajet) && (etajet <= -minEta_)) {
0098         nNegJets++;
0099         TRef ref = TRef(objects, distance(objects->begin(), jet));
0100         filterproduct.addObject(triggerType_, ref);
0101       }
0102     }
0103   }
0104 
0105   // filter decision
0106   const bool accept((nNegJets >= nNeg_) && (nPosJets >= nPos_) && ((nNegJets + nPosJets) >= nTot_));
0107 
0108   return accept;
0109 }