File indexing completed on 2024-04-06 12:18:41
0001 #include "HLTPixelActivityHFSumEnergyFilter.h"
0002
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006
0007
0008
0009
0010
0011 HLTPixelActivityHFSumEnergyFilter::HLTPixelActivityHFSumEnergyFilter(const edm::ParameterSet& config)
0012 : inputTag_(config.getParameter<edm::InputTag>("inputTag")),
0013 HFHits_(config.getParameter<edm::InputTag>("HFHitCollection")),
0014 eCut_HF_(config.getParameter<double>("eCut_HF")),
0015 eMin_HF_(config.getParameter<double>("eMin_HF")),
0016 offset_(config.getParameter<double>("offset")),
0017 slope_(config.getParameter<double>("slope")) {
0018 inputToken_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(inputTag_);
0019 HFHitsToken_ = consumes<HFRecHitCollection>(HFHits_);
0020
0021 LogDebug("") << "Using the " << inputTag_ << " input collection";
0022 }
0023
0024 HLTPixelActivityHFSumEnergyFilter::~HLTPixelActivityHFSumEnergyFilter() = default;
0025
0026 void HLTPixelActivityHFSumEnergyFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0027 edm::ParameterSetDescription desc;
0028 desc.add<edm::InputTag>("inputTag", edm::InputTag("hltSiPixelClusters"));
0029 desc.add<edm::InputTag>("HFHitCollection", edm::InputTag("hltHfreco"));
0030 desc.add<double>("eCut_HF", 0);
0031 desc.add<double>("eMin_HF", 10000.);
0032 desc.add<double>("offset", -1000.);
0033 desc.add<double>("slope", 0.5);
0034 descriptions.add("hltPixelActivityHFSumEnergyFilter", desc);
0035 }
0036
0037
0038
0039
0040
0041
0042 bool HLTPixelActivityHFSumEnergyFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0043
0044 edm::Handle<edmNew::DetSetVector<SiPixelCluster> > clusterColl;
0045 iEvent.getByToken(inputToken_, clusterColl);
0046
0047 unsigned int clusterSize = clusterColl->dataSize();
0048 LogDebug("") << "Number of clusters: " << clusterSize;
0049
0050 edm::Handle<HFRecHitCollection> HFRecHitsH;
0051 iEvent.getByToken(HFHitsToken_, HFRecHitsH);
0052
0053 double sumE = 0.;
0054
0055 for (auto const& it : *HFRecHitsH) {
0056 if (it.energy() > eCut_HF_) {
0057 sumE += it.energy();
0058 }
0059 }
0060
0061 bool accept = false;
0062
0063 double thres = offset_ + slope_ * sumE;
0064 double diff = clusterSize - thres;
0065 if (sumE > eMin_HF_ && diff < 0.)
0066 accept = true;
0067
0068
0069 return accept;
0070 }
0071
0072
0073 #include "FWCore/Framework/interface/MakerMacros.h"
0074 DEFINE_FWK_MODULE(HLTPixelActivityHFSumEnergyFilter);