Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    HLTHcalSimpleRecHitFilter
0004 // Class:      HLTHcalSimpleRecHitFilter
0005 //
0006 /**\class HLTHcalSimpleRecHitFilter HLTHcalSimpleRecHitFilter.cc Work/HLTHcalSimpleRecHitFilter/src/HLTHcalSimpleRecHitFilter.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Bryan DAHMES
0015 //         Created:  Wed Sep 19 16:21:29 CEST 2007
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0025 
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0032 
0033 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0034 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0035 
0036 //
0037 // class declaration
0038 //
0039 
0040 class HLTHcalSimpleRecHitFilter : public HLTFilter {
0041 public:
0042   explicit HLTHcalSimpleRecHitFilter(const edm::ParameterSet&);
0043   ~HLTHcalSimpleRecHitFilter() override;
0044   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0045 
0046 private:
0047   bool hltFilter(edm::Event&,
0048                  const edm::EventSetup&,
0049                  trigger::TriggerFilterObjectWithRefs& filterproduct) const override;
0050 
0051   // ----------member data ---------------------------
0052   edm::EDGetTokenT<HFRecHitCollection> HcalRecHitsToken_;
0053   edm::InputTag HcalRecHitCollection_;
0054   double threshold_;
0055   int minNHitsNeg_;
0056   int minNHitsPos_;
0057   bool doCoincidence_;
0058   std::vector<unsigned int> maskedList_;
0059 };
0060 
0061 //
0062 // constructors and destructor
0063 //
0064 HLTHcalSimpleRecHitFilter::HLTHcalSimpleRecHitFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0065   //now do what ever initialization is needed
0066   threshold_ = iConfig.getParameter<double>("threshold");
0067   minNHitsNeg_ = iConfig.getParameter<int>("minNHitsNeg");
0068   minNHitsPos_ = iConfig.getParameter<int>("minNHitsPos");
0069   doCoincidence_ = iConfig.getParameter<bool>("doCoincidence");
0070   if (iConfig.existsAs<std::vector<unsigned int> >("maskedChannels"))
0071     maskedList_ =
0072         iConfig.getParameter<std::vector<unsigned int> >("maskedChannels");  //this is using the raw DetId index
0073   else
0074     //worry about possible user menus with the old interface
0075     if (iConfig.existsAs<std::vector<int> >("maskedChannels")) {
0076       std::vector<int> tVec = iConfig.getParameter<std::vector<int> >("maskedChannels");
0077       if (!tVec.empty()) {
0078         edm::LogError("cfg error") << "masked list of channels missing from HLT menu. Migration from vector of ints to "
0079                                       "vector of uints needed for this release";
0080         cms::Exception("Invalid/obsolete masked list of channels");
0081       }
0082     }
0083   HcalRecHitCollection_ = iConfig.getParameter<edm::InputTag>("HFRecHitCollection");
0084   HcalRecHitsToken_ = consumes<HFRecHitCollection>(HcalRecHitCollection_);
0085 }
0086 
0087 HLTHcalSimpleRecHitFilter::~HLTHcalSimpleRecHitFilter() {
0088   // do anything here that needs to be done at desctruction time
0089   // (e.g. close files, deallocate resources etc.)
0090 }
0091 
0092 void HLTHcalSimpleRecHitFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0093   edm::ParameterSetDescription desc;
0094   makeHLTFilterDescription(desc);
0095   desc.add<edm::InputTag>("HFRecHitCollection", edm::InputTag("hltHfreco"));
0096   desc.add<double>("threshold", 3.0);
0097   desc.add<int>("minNHitsNeg", 1);
0098   desc.add<int>("minNHitsPos", 1);
0099   desc.add<bool>("doCoincidence", true);
0100   std::vector<unsigned int> temp;
0101   desc.add<std::vector<unsigned int> >("maskedChannels", temp)->setComment(" # now by raw detid, not hashed id");
0102   descriptions.add("hltHcalSimpleRecHitFilter", desc);
0103 }
0104 
0105 //
0106 // member functions
0107 //
0108 
0109 // ------------ method called on each new Event  ------------
0110 bool HLTHcalSimpleRecHitFilter::hltFilter(edm::Event& iEvent,
0111                                           const edm::EventSetup& iSetup,
0112                                           trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0113   // using namespace edm;
0114 
0115   // getting very basic uncalRH
0116   edm::Handle<HFRecHitCollection> crudeHits;
0117   try {
0118     iEvent.getByToken(HcalRecHitsToken_, crudeHits);
0119   } catch (std::exception& ex) {
0120     edm::LogWarning("HLTHcalSimpleRecHitFilter") << HcalRecHitCollection_ << " not available";
0121   }
0122 
0123   bool accept = false;
0124 
0125   int nHitsNeg = 0, nHitsPos = 0;
0126   for (auto hit : *crudeHits) {
0127     // masking noisy channels
0128     if (std::find(maskedList_.begin(), maskedList_.end(), hit.id().rawId()) != maskedList_.end())
0129       continue;
0130 
0131     // only count tower above threshold
0132     if (hit.energy() < threshold_)
0133       continue;
0134 
0135     // count
0136     if (hit.id().zside() < 0)
0137       ++nHitsNeg;
0138     else
0139       ++nHitsPos;
0140   }
0141 
0142   // Logic
0143   if (!doCoincidence_)
0144     accept = (nHitsNeg >= minNHitsNeg_) || (nHitsPos >= minNHitsPos_);
0145   else
0146     accept = (nHitsNeg >= minNHitsNeg_) && (nHitsPos >= minNHitsPos_);
0147   //  edm::LogInfo("HcalFilter")  << "at evet: " << iEvent.id().event()
0148   //    << " and run: " << iEvent.id().run()
0149   //    << " Total HF hits: " << crudeHits->size() << " Above Threshold - nNeg: " << nHitsNeg << " nPos " << nHitsPos
0150   //    << " doCoinc: " << doCoincidence_ << " accept: " << accept << std::endl;
0151 
0152   // result
0153   return accept;
0154 }
0155 
0156 // declare this class as a framework plugin
0157 DEFINE_FWK_MODULE(HLTHcalSimpleRecHitFilter);