Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:36:37

0001 /** \class HLTJetVBFFilter
0002  *
0003  *
0004  *  \author Monica Vazquez Acosta (CERN)
0005  *  \modifier Phst Srimanobhas (srimanob@mail.cern.ch)
0006  *
0007  */
0008 
0009 #include "FWCore/Framework/interface/ESHandle.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "DataFormats/Common/interface/Ref.h"
0017 #include "DataFormats/Common/interface/Handle.h"
0018 #include "HLTrigger/JetMET/interface/HLTJetVBFFilter.h"
0019 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0020 
0021 //
0022 // constructors and destructor
0023 //
0024 template <typename T>
0025 HLTJetVBFFilter<T>::HLTJetVBFFilter(const edm::ParameterSet& iConfig)
0026     : HLTFilter(iConfig),
0027       inputTag_(iConfig.template getParameter<edm::InputTag>("inputTag")),
0028       minPtLow_(iConfig.template getParameter<double>("minPtLow")),
0029       minPtHigh_(iConfig.template getParameter<double>("minPtHigh")),
0030       etaOpposite_(iConfig.template getParameter<bool>("etaOpposite")),
0031       minDeltaEta_(iConfig.template getParameter<double>("minDeltaEta")),
0032       minInvMass_(iConfig.template getParameter<double>("minInvMass")),
0033       maxEta_(iConfig.template getParameter<double>("maxEta")),
0034       leadingJetOnly_(iConfig.template getParameter<bool>("leadingJetOnly")),
0035       triggerType_(iConfig.template getParameter<int>("triggerType")) {
0036   m_theObjectToken = consumes<std::vector<T>>(inputTag_);
0037   LogDebug("") << "HLTJetVBFFilter: Input/minPtLow_/minPtHigh_/triggerType : " << inputTag_.encode() << " " << minPtLow_
0038                << " " << minPtHigh_ << " " << triggerType_;
0039 }
0040 
0041 template <typename T>
0042 HLTJetVBFFilter<T>::~HLTJetVBFFilter() = default;
0043 
0044 template <typename T>
0045 void HLTJetVBFFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0046   edm::ParameterSetDescription desc;
0047   makeHLTFilterDescription(desc);
0048   desc.add<edm::InputTag>("inputTag", edm::InputTag("hltAntiKT5ConvPFJets"));
0049   desc.add<double>("minPtLow", 40.);
0050   desc.add<double>("minPtHigh", 40.);
0051   desc.add<bool>("etaOpposite", false);
0052   desc.add<double>("minDeltaEta", 4.0);
0053   desc.add<double>("minInvMass", 1000.);
0054   desc.add<double>("maxEta", 5.0);
0055   desc.add<bool>("leadingJetOnly", false);
0056   desc.add<int>("triggerType", trigger::TriggerJet);
0057   descriptions.add(defaultModuleLabel<HLTJetVBFFilter<T>>(), desc);
0058 }
0059 
0060 //
0061 // ------------ method called to produce the data  ------------
0062 //
0063 template <typename T>
0064 bool HLTJetVBFFilter<T>::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   using namespace trigger;
0071 
0072   typedef vector<T> TCollection;
0073   typedef Ref<TCollection> TRef;
0074 
0075   // The filter object
0076   if (saveTags())
0077     filterproduct.addCollectionTag(inputTag_);
0078 
0079   // get hold of collection of objects
0080   Handle<TCollection> objects;
0081   iEvent.getByToken(m_theObjectToken, objects);
0082 
0083   // look at all candidates, check cuts and add to filter object
0084   int n(0);
0085 
0086   // events with two or more jets
0087   if (objects->size() > 1) {
0088     double ejet1 = 0.;
0089     double pxjet1 = 0.;
0090     double pyjet1 = 0.;
0091     double pzjet1 = 0.;
0092     double ptjet1 = 0.;
0093     double etajet1 = 0.;
0094 
0095     double ejet2 = 0.;
0096     double pxjet2 = 0.;
0097     double pyjet2 = 0.;
0098     double pzjet2 = 0.;
0099     double ptjet2 = 0.;
0100     double etajet2 = 0.;
0101 
0102     // loop on all jets
0103     int countJet1(0);
0104     int countJet2(0);
0105     typename TCollection::const_iterator jet1(objects->begin());
0106     for (; jet1 != objects->end(); jet1++) {
0107       countJet1++;
0108       if (leadingJetOnly_ == true && countJet1 > 2)
0109         break;
0110       //
0111       if (jet1->pt() < minPtHigh_)
0112         break;  //No need to go to the next jet (lower PT)
0113       if (std::abs(jet1->eta()) > maxEta_)
0114         continue;
0115       //
0116       countJet2 = countJet1 - 1;
0117       typename TCollection::const_iterator jet2(jet1 + 1);
0118       for (; jet2 != objects->end(); jet2++) {
0119         countJet2++;
0120         if (leadingJetOnly_ == true && countJet2 > 2)
0121           break;
0122         //
0123         if (jet2->pt() < minPtLow_)
0124           break;  //No need to go to the next jet (lower PT)
0125         if (std::abs(jet2->eta()) > maxEta_)
0126           continue;
0127         //
0128         ejet1 = jet1->energy();
0129         pxjet1 = jet1->px();
0130         pyjet1 = jet1->py();
0131         pzjet1 = jet1->pz();
0132         ptjet1 = jet1->pt();
0133         etajet1 = jet1->eta();
0134 
0135         ejet2 = jet2->energy();
0136         pxjet2 = jet2->px();
0137         pyjet2 = jet2->py();
0138         pzjet2 = jet2->pz();
0139         ptjet2 = jet2->pt();
0140         etajet2 = jet2->eta();
0141         //
0142         float deltaetajet = etajet1 - etajet2;
0143         float invmassjet = sqrt((ejet1 + ejet2) * (ejet1 + ejet2) - (pxjet1 + pxjet2) * (pxjet1 + pxjet2) -
0144                                 (pyjet1 + pyjet2) * (pyjet1 + pyjet2) - (pzjet1 + pzjet2) * (pzjet1 + pzjet2));
0145 
0146         // VBF cuts
0147         if ((ptjet1 > minPtHigh_) && (ptjet2 > minPtLow_) &&
0148             ((etaOpposite_ == true && etajet1 * etajet2 < 0) || (etaOpposite_ == false)) &&
0149             (std::abs(deltaetajet) > minDeltaEta_) && (std::abs(invmassjet) > minInvMass_)) {
0150           ++n;
0151           TRef ref1 = TRef(objects, distance(objects->begin(), jet1));
0152           TRef ref2 = TRef(objects, distance(objects->begin(), jet2));
0153           filterproduct.addObject(triggerType_, ref1);
0154           filterproduct.addObject(triggerType_, ref2);
0155         }  // VBF cuts
0156         //if(n>=1) break; //Store all possible pairs
0157       }
0158       //if(n>=1) break; //Store all possible pairs
0159     }  // loop on all jets
0160   }  // events with two or more jets
0161 
0162   // filter decision
0163   bool accept(n >= 1);
0164 
0165   return accept;
0166 }