Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTJetCollectionsVBFFilter
0002  *
0003  *
0004  *  \author Leonard Apanasevich
0005  *
0006  */
0007 
0008 #include "HLTrigger/JetMET/interface/HLTJetCollectionsVBFFilter.h"
0009 #include "DataFormats/Common/interface/Handle.h"
0010 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "DataFormats/Math/interface/deltaPhi.h"
0015 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0017 #include "FWCore/Utilities/interface/InputTag.h"
0018 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0019 
0020 //
0021 // constructors and destructor
0022 //
0023 template <typename T>
0024 HLTJetCollectionsVBFFilter<T>::HLTJetCollectionsVBFFilter(const edm::ParameterSet& iConfig)
0025     : HLTFilter(iConfig),
0026       inputTag_(iConfig.getParameter<edm::InputTag>("inputTag")),
0027       originalTag_(iConfig.getParameter<edm::InputTag>("originalTag")),
0028       softJetPt_(iConfig.getParameter<double>("SoftJetPt")),
0029       hardJetPt_(iConfig.getParameter<double>("HardJetPt")),
0030       minDeltaEta_(iConfig.getParameter<double>("MinDeltaEta")),
0031       thirdJetPt_(iConfig.getParameter<double>("ThirdJetPt")),
0032       maxAbsJetEta_(iConfig.getParameter<double>("MaxAbsJetEta")),
0033       maxAbsThirdJetEta_(iConfig.getParameter<double>("MaxAbsThirdJetEta")),
0034       minNJets_(iConfig.getParameter<unsigned int>("MinNJets")),
0035       triggerType_(iConfig.getParameter<int>("TriggerType")) {
0036   typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
0037       TCollectionVector;
0038   m_theJetToken = consumes<TCollectionVector>(inputTag_);
0039 }
0040 
0041 template <typename T>
0042 HLTJetCollectionsVBFFilter<T>::~HLTJetCollectionsVBFFilter() = default;
0043 
0044 template <typename T>
0045 void HLTJetCollectionsVBFFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0046   edm::ParameterSetDescription desc;
0047   makeHLTFilterDescription(desc);
0048   desc.add<edm::InputTag>("inputTag", edm::InputTag("hltIterativeCone5CaloJets"));
0049   desc.add<edm::InputTag>("originalTag", edm::InputTag("hltIterativeCone5CaloJets"));
0050   desc.add<double>("SoftJetPt", 25.0);
0051   desc.add<double>("HardJetPt", 35.0);
0052   desc.add<double>("MinDeltaEta", 3.0);
0053   desc.add<double>("ThirdJetPt", 20.0);
0054   desc.add<double>("MaxAbsJetEta", 9999.);
0055   desc.add<double>("MaxAbsThirdJetEta", 2.6);
0056   desc.add<unsigned int>("MinNJets", 2);
0057   desc.add<int>("TriggerType", trigger::TriggerJet);
0058   descriptions.add(defaultModuleLabel<HLTJetCollectionsVBFFilter<T>>(), desc);
0059 }
0060 
0061 // ------------ method called to produce the data  ------------
0062 template <typename T>
0063 bool HLTJetCollectionsVBFFilter<T>::hltFilter(edm::Event& iEvent,
0064                                               const edm::EventSetup& iSetup,
0065                                               trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0066   using namespace std;
0067   using namespace edm;
0068   using namespace reco;
0069   using namespace trigger;
0070 
0071   typedef vector<T> TCollection;
0072   typedef Ref<TCollection> TRef;
0073   typedef edm::RefVector<TCollection> TRefVector;
0074   typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
0075       TCollectionVector;
0076 
0077   // The filter object
0078   if (saveTags())
0079     filterproduct.addCollectionTag(originalTag_);
0080 
0081   Handle<TCollectionVector> theJetCollectionsHandle;
0082   iEvent.getByToken(m_theJetToken, theJetCollectionsHandle);
0083   const TCollectionVector& theJetCollections = *theJetCollectionsHandle;
0084   // filter decision
0085   bool accept(false);
0086   std::vector<TRef> goodJetRefs;
0087 
0088   for (unsigned int collection = 0; collection < theJetCollections.size(); ++collection) {
0089     const TRefVector& refVector = theJetCollections[collection];
0090     if (refVector.size() < minNJets_)
0091       continue;
0092 
0093     // VBF decision
0094     bool thereAreVBFJets(false);
0095     // 3rd Jet check decision
0096     bool goodThirdJet(false);
0097     if (minNJets_ < 3)
0098       goodThirdJet = true;
0099 
0100     //empty the good jets collection
0101     goodJetRefs.clear();
0102 
0103     TRef refOne;
0104     TRef refTwo;
0105     typename TRefVector::const_iterator jetOne(refVector.begin());
0106     int firstJetIndex = 100, secondJetIndex = 100, thirdJetIndex = 100;
0107 
0108     // Cycle to look for VBF jets
0109     for (; jetOne != refVector.end(); jetOne++) {
0110       TRef jetOneRef(*jetOne);
0111 
0112       if (thereAreVBFJets)
0113         break;
0114       if (jetOneRef->pt() < hardJetPt_)
0115         break;
0116       if (std::abs(jetOneRef->eta()) > maxAbsJetEta_)
0117         continue;
0118 
0119       typename TRefVector::const_iterator jetTwo = jetOne + 1;
0120       secondJetIndex = firstJetIndex;
0121       for (; jetTwo != refVector.end(); jetTwo++) {
0122         TRef jetTwoRef(*jetTwo);
0123 
0124         if (jetTwoRef->pt() < softJetPt_)
0125           break;
0126         if (std::abs(jetTwoRef->eta()) > maxAbsJetEta_)
0127           continue;
0128 
0129         if (std::abs(jetTwoRef->eta() - jetOneRef->eta()) < minDeltaEta_)
0130           continue;
0131 
0132         thereAreVBFJets = true;
0133         refOne = *jetOne;
0134         goodJetRefs.push_back(refOne);
0135         refTwo = *jetTwo;
0136         goodJetRefs.push_back(refTwo);
0137 
0138         firstJetIndex = (int)(jetOne - refVector.begin());
0139         secondJetIndex = (int)(jetTwo - refVector.begin());
0140 
0141         break;
0142       }
0143     }  // Close looop on VBF
0144 
0145     // Look for a third jet, if you've found the previous 2
0146     if (minNJets_ > 2 && thereAreVBFJets) {
0147       TRef refThree;
0148       typename TRefVector::const_iterator jetThree(refVector.begin());
0149       for (; jetThree != refVector.end(); jetThree++) {
0150         thirdJetIndex = (int)(jetThree - refVector.begin());
0151 
0152         TRef jetThreeRef(*jetThree);
0153 
0154         if (thirdJetIndex == firstJetIndex || thirdJetIndex == secondJetIndex)
0155           continue;
0156 
0157         if (jetThreeRef->pt() >= thirdJetPt_ && std::abs(jetThreeRef->eta()) <= maxAbsThirdJetEta_) {
0158           goodThirdJet = true;
0159           refThree = *jetThree;
0160           goodJetRefs.push_back(refThree);
0161           break;
0162         }
0163       }
0164     }
0165 
0166     if (thereAreVBFJets && goodThirdJet) {
0167       accept = true;
0168       break;
0169     }
0170   }
0171 
0172   //fill the filter object
0173   for (unsigned int refIndex = 0; refIndex < goodJetRefs.size(); ++refIndex) {
0174     filterproduct.addObject(triggerType_, goodJetRefs.at(refIndex));
0175   }
0176 
0177   return accept;
0178 }