Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTDiJetAveEtaFilter
0002  *
0003  *
0004  *  \author Tomasz Fruboes
0005  *     based on HLTDiJetAveFilter
0006  */
0007 
0008 #include "HLTrigger/JetMET/interface/HLTDiJetAveEtaFilter.h"
0009 #include "DataFormats/Common/interface/Ref.h"
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/Math/interface/deltaPhi.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "FWCore/Utilities/interface/InputTag.h"
0019 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0020 
0021 //
0022 // constructors and destructor
0023 //
0024 template <typename T>
0025 HLTDiJetAveEtaFilter<T>::HLTDiJetAveEtaFilter(const edm::ParameterSet& iConfig)
0026     : HLTFilter(iConfig),
0027       inputJetTag_(iConfig.template getParameter<edm::InputTag>("inputJetTag")),
0028       minPtJet_(iConfig.template getParameter<double>("minPtJet")),
0029       minPtAve_(iConfig.template getParameter<double>("minPtAve")),
0030       //minPtJet3_   (iConfig.template getParameter<double> ("minPtJet3")),
0031       minDphi_(iConfig.template getParameter<double>("minDphi")),
0032       tagEtaMin_(iConfig.template getParameter<double>("minTagEta")),
0033       tagEtaMax_(iConfig.template getParameter<double>("maxTagEta")),
0034       probeEtaMin_(iConfig.template getParameter<double>("minProbeEta")),
0035       probeEtaMax_(iConfig.template getParameter<double>("maxProbeEta")),
0036       triggerType_(iConfig.template getParameter<int>("triggerType")) {
0037   m_theJetToken = consumes<std::vector<T>>(inputJetTag_);
0038   LogDebug("") << "HLTDiJetAveEtaFilter: Input/minPtAve/minDphi/triggerType : " << inputJetTag_.encode() << " "
0039                << minPtAve_
0040                << " "
0041                //<< minPtJet3_ << " "
0042                << minDphi_ << " " << triggerType_;
0043 }
0044 
0045 template <typename T>
0046 HLTDiJetAveEtaFilter<T>::~HLTDiJetAveEtaFilter() = default;
0047 
0048 template <typename T>
0049 void HLTDiJetAveEtaFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0050   edm::ParameterSetDescription desc;
0051   makeHLTFilterDescription(desc);
0052   desc.add<edm::InputTag>("inputJetTag", edm::InputTag("hltIterativeCone5CaloJets"));
0053   desc.add<double>("minPtAve", 100.0);
0054   desc.add<double>("minPtJet", 50.0);
0055   //desc.add<double>("minPtJet3",99999.0);
0056   desc.add<double>("minDphi", -1.0);
0057   desc.add<double>("minTagEta", -1.);
0058   desc.add<double>("maxTagEta", 1.4);
0059   desc.add<double>("minProbeEta", 2.7);
0060   desc.add<double>("maxProbeEta", 5.5);
0061   desc.add<int>("triggerType", trigger::TriggerJet);
0062   descriptions.add(defaultModuleLabel<HLTDiJetAveEtaFilter<T>>(), desc);
0063 }
0064 
0065 // ------------ method called to produce the data  ------------
0066 template <typename T>
0067 bool HLTDiJetAveEtaFilter<T>::hltFilter(edm::Event& iEvent,
0068                                         const edm::EventSetup& iSetup,
0069                                         trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0070   using namespace std;
0071   using namespace edm;
0072   using namespace reco;
0073   using namespace trigger;
0074 
0075   typedef vector<T> TCollection;
0076   typedef Ref<TCollection> TRef;
0077 
0078   // The filter object
0079   if (saveTags())
0080     filterproduct.addCollectionTag(inputJetTag_);
0081 
0082   // get hold of collection of objects
0083   Handle<TCollection> objects;
0084   iEvent.getByToken(m_theJetToken, objects);
0085 
0086   int n(0);
0087 
0088   if (objects->size() > 1) {  // events with two or more jets
0089     typename TCollection::const_iterator iProbe(objects->begin());
0090     typename TCollection::const_iterator iEnd(objects->end());
0091     for (; iProbe != iEnd; ++iProbe) {
0092       if (iProbe->pt() < minPtJet_)
0093         continue;
0094 
0095       // for easier trigger efficiency evaluation save all probe/tag
0096       // objects passing the minPT/eta criteria (outer loop)
0097       float eta = std::abs(iProbe->eta());
0098       bool isGood = false;  // probe or tag
0099       bool isProbe = false;
0100       if (eta > probeEtaMin_ && eta < probeEtaMax_) {
0101         isGood = true;
0102         isProbe = true;
0103       }
0104       if (eta > tagEtaMin_ && eta < tagEtaMax_) {
0105         isGood = true;
0106       }
0107       if (isGood) {
0108         filterproduct.addObject(triggerType_, TRef(objects, distance(objects->begin(), iProbe)));
0109       }
0110 
0111       if (!isProbe)
0112         continue;
0113 
0114       typename TCollection::const_iterator iTag(objects->begin());
0115       for (; iTag != iEnd; ++iTag) {
0116         if (iTag == iProbe)
0117           continue;
0118         if (iTag->pt() < minPtJet_)
0119           continue;
0120         float eta2 = std::abs(iTag->eta());
0121         if (eta2 < tagEtaMin_ || eta2 > tagEtaMax_)
0122           continue;
0123         double dphi = std::abs(deltaPhi(iProbe->phi(), iTag->phi()));
0124         if (dphi < minDphi_) {
0125           continue;
0126         }
0127 
0128         double ptAve = (iProbe->pt() + iTag->pt()) / 2;
0129         if (ptAve < minPtAve_) {
0130           continue;
0131         }
0132         ++n;
0133       }
0134     }
0135   }  // events with two or more jets
0136   // filter decision
0137   bool accept(n >= 1);
0138   return accept;
0139 }