Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class HLTExclDiJetFilter
0002  *
0003  *
0004  *  \author Leonard Apanasevich
0005  *
0006  */
0007 
0008 #include <cmath>
0009 
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "DataFormats/Common/interface/Ref.h"
0018 #include "DataFormats/Common/interface/Handle.h"
0019 #include "HLTrigger/JetMET/interface/HLTExclDiJetFilter.h"
0020 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0021 
0022 //
0023 // constructors and destructor
0024 //
0025 template <typename T>
0026 HLTExclDiJetFilter<T>::HLTExclDiJetFilter(const edm::ParameterSet& iConfig)
0027     : HLTFilter(iConfig),
0028       inputJetTag_(iConfig.template getParameter<edm::InputTag>("inputJetTag")),
0029       caloTowerTag_(iConfig.template getParameter<edm::InputTag>("caloTowerTag")),
0030       minPtJet_(iConfig.template getParameter<double>("minPtJet")),
0031       minHFe_(iConfig.template getParameter<double>("minHFe")),
0032       HF_OR_(iConfig.template getParameter<bool>("HF_OR")),
0033       triggerType_(iConfig.template getParameter<int>("triggerType")) {
0034   m_theJetToken = consumes<std::vector<T>>(inputJetTag_);
0035   m_theCaloTowerCollectionToken = consumes<CaloTowerCollection>(caloTowerTag_);
0036   LogDebug("") << "HLTExclDiJetFilter: Input/minPtJet/minHFe/HF_OR/triggerType : " << inputJetTag_.encode() << " "
0037                << caloTowerTag_.encode() << " " << minPtJet_ << " " << minHFe_ << " " << HF_OR_ << " " << triggerType_;
0038 }
0039 
0040 template <typename T>
0041 HLTExclDiJetFilter<T>::~HLTExclDiJetFilter() = default;
0042 
0043 template <typename T>
0044 void HLTExclDiJetFilter<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0045   edm::ParameterSetDescription desc;
0046   makeHLTFilterDescription(desc);
0047   desc.add<edm::InputTag>("inputJetTag", edm::InputTag("hltMCJetCorJetIcone5HF07"));
0048   desc.add<edm::InputTag>("caloTowerTag", edm::InputTag("hltTowerMakerForAll"));
0049   desc.add<double>("minPtJet", 30.0);
0050   desc.add<double>("minHFe", 50.0);
0051   desc.add<bool>("HF_OR", false);
0052   desc.add<int>("triggerType", trigger::TriggerJet);
0053   descriptions.add(defaultModuleLabel<HLTExclDiJetFilter<T>>(), desc);
0054 }
0055 
0056 // ------------ method called to produce the data  ------------
0057 template <typename T>
0058 bool HLTExclDiJetFilter<T>::hltFilter(edm::Event& iEvent,
0059                                       const edm::EventSetup& iSetup,
0060                                       trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0061   using namespace std;
0062   using namespace edm;
0063   using namespace reco;
0064   using namespace trigger;
0065 
0066   typedef vector<T> TCollection;
0067   typedef Ref<TCollection> TRef;
0068 
0069   // The filter object
0070   if (saveTags())
0071     filterproduct.addCollectionTag(inputJetTag_);
0072 
0073   Handle<TCollection> recojets;  //recojets can be any jet collections
0074   iEvent.getByToken(m_theJetToken, recojets);
0075 
0076   // look at all candidates,  check cuts and add to filter object
0077   int n(0);
0078 
0079   double ptjet1 = 0., ptjet2 = 0.;
0080   double phijet1 = 0., phijet2 = 0.;
0081 
0082   if (recojets->size() > 1) {
0083     // events with two or more jets
0084 
0085     int countjets = 0;
0086 
0087     TRef JetRef1, JetRef2;
0088 
0089     typename TCollection::const_iterator recojet(recojets->begin());
0090     for (; recojet <= (recojets->begin() + 1); ++recojet) {
0091       //
0092       if (countjets == 0) {
0093         ptjet1 = recojet->pt();
0094         phijet1 = recojet->phi();
0095 
0096         JetRef1 = TRef(recojets, distance(recojets->begin(), recojet));
0097       }
0098       //
0099       if (countjets == 1) {
0100         ptjet2 = recojet->pt();
0101         phijet2 = recojet->phi();
0102 
0103         JetRef2 = TRef(recojets, distance(recojets->begin(), recojet));
0104       }
0105       //
0106       ++countjets;
0107     }
0108 
0109     if (ptjet1 > minPtJet_ && ptjet2 > minPtJet_) {
0110       double Dphi = std::abs(phijet1 - phijet2);
0111       if (Dphi > M_PI)
0112         Dphi = 2.0 * M_PI - Dphi;
0113       if (Dphi > 0.5 * M_PI) {
0114         filterproduct.addObject(triggerType_, JetRef1);
0115         filterproduct.addObject(triggerType_, JetRef2);
0116         ++n;
0117       }
0118     }  // pt(jet1,jet2) > minPtJet_
0119 
0120   }  // events with two or more jets
0121 
0122   // calotowers
0123   bool hf_accept = false;
0124 
0125   if (n > 0) {
0126     double ehfp(0.);
0127     double ehfm(0.);
0128 
0129     Handle<CaloTowerCollection> o;
0130     iEvent.getByToken(m_theCaloTowerCollectionToken, o);
0131     //     if( o.isValid()) {
0132     for (auto const& cc : *o) {
0133       if (std::abs(cc.ieta()) > 28 && cc.energy() < 4.0)
0134         continue;
0135       if (cc.ieta() > 28)
0136         ehfp += cc.energy();  // HF+ energy
0137       if (cc.ieta() < -28)
0138         ehfm += cc.energy();  // HF- energy
0139     }
0140     //    }
0141 
0142     bool hf_accept_and = (ehfp < minHFe_) && (ehfm < minHFe_);
0143     bool hf_accept_or = (ehfp < minHFe_) || (ehfm < minHFe_);
0144 
0145     hf_accept = HF_OR_ ? hf_accept_or : hf_accept_and;
0146 
0147   }  // n>0
0148 
0149   ////////////////////////////////////////////////////////
0150 
0151   // filter decision
0152   bool accept(n > 0 && hf_accept);
0153 
0154   return accept;
0155 }