File indexing completed on 2024-04-06 12:18:39
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024
0025 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0026
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0033
0034 #include "FWCore/Framework/interface/EventSetup.h"
0035 #include "FWCore/Framework/interface/ESHandle.h"
0036
0037 #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
0038
0039 #include <string>
0040
0041
0042
0043
0044
0045 class HLTCSCAcceptBusyFilter : public HLTFilter {
0046 public:
0047 explicit HLTCSCAcceptBusyFilter(const edm::ParameterSet&);
0048 ~HLTCSCAcceptBusyFilter() override;
0049 bool hltFilter(edm::Event&,
0050 const edm::EventSetup&,
0051 trigger::TriggerFilterObjectWithRefs& filterproduct) const override;
0052 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0053
0054 private:
0055 bool AcceptManyHitsInChamber(unsigned int maxRecHitsPerChamber,
0056 const edm::Handle<CSCRecHit2DCollection>& recHits) const;
0057
0058
0059 edm::EDGetTokenT<CSCRecHit2DCollection> cscrechitsToken;
0060 edm::InputTag cscrechitsTag;
0061 bool invert;
0062 unsigned int maxRecHitsPerChamber;
0063 };
0064
0065
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076 HLTCSCAcceptBusyFilter::HLTCSCAcceptBusyFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0077
0078 cscrechitsTag = iConfig.getParameter<edm::InputTag>("cscrechitsTag");
0079 invert = iConfig.getParameter<bool>("invert");
0080 maxRecHitsPerChamber = iConfig.getParameter<unsigned int>("maxRecHitsPerChamber");
0081 cscrechitsToken = consumes<CSCRecHit2DCollection>(cscrechitsTag);
0082 }
0083
0084 HLTCSCAcceptBusyFilter::~HLTCSCAcceptBusyFilter() {
0085
0086
0087 }
0088
0089 void HLTCSCAcceptBusyFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0090 edm::ParameterSetDescription desc;
0091 makeHLTFilterDescription(desc);
0092 desc.add<edm::InputTag>("cscrechitsTag", edm::InputTag("hltCsc2DRecHits"));
0093 desc.add<bool>("invert", true);
0094 desc.add<unsigned int>("maxRecHitsPerChamber", 200);
0095 descriptions.add("hltCSCAcceptBusyFilter", desc);
0096 }
0097
0098
0099
0100
0101
0102
0103 bool HLTCSCAcceptBusyFilter::hltFilter(edm::Event& iEvent,
0104 const edm::EventSetup& iSetup,
0105 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0106 using namespace edm;
0107
0108
0109 Handle<CSCRecHit2DCollection> recHits;
0110 iEvent.getByToken(cscrechitsToken, recHits);
0111
0112 if (AcceptManyHitsInChamber(maxRecHitsPerChamber, recHits)) {
0113 return (!invert);
0114 } else {
0115 return (invert);
0116 }
0117 }
0118
0119
0120 bool HLTCSCAcceptBusyFilter::AcceptManyHitsInChamber(unsigned int maxRecHitsPerChamber,
0121 const edm::Handle<CSCRecHit2DCollection>& recHits) const {
0122 unsigned int maxNRecHitsPerChamber(0);
0123
0124 const unsigned int nEndcaps(2);
0125 const unsigned int nStations(4);
0126 const unsigned int nRings(4);
0127 const unsigned int nChambers(36);
0128 unsigned int allRechits[nEndcaps][nStations][nRings][nChambers];
0129 for (auto& allRechit : allRechits) {
0130 for (unsigned int iS = 0; iS < nStations; ++iS) {
0131 for (unsigned int iR = 0; iR < nRings; ++iR) {
0132 for (unsigned int iC = 0; iC < nChambers; ++iC) {
0133 allRechit[iS][iR][iC] = 0;
0134 }
0135 }
0136 }
0137 }
0138
0139 for (auto const& it : *recHits) {
0140 ++allRechits[it.cscDetId().endcap() - 1][it.cscDetId().station() - 1][it.cscDetId().ring() - 1]
0141 [it.cscDetId().chamber() - 1];
0142 }
0143
0144 for (auto& allRechit : allRechits) {
0145 for (unsigned int iS = 0; iS < nStations; ++iS) {
0146 for (unsigned int iR = 0; iR < nRings; ++iR) {
0147 for (unsigned int iC = 0; iC < nChambers; ++iC) {
0148 if (allRechit[iS][iR][iC] > maxNRecHitsPerChamber) {
0149 maxNRecHitsPerChamber = allRechit[iS][iR][iC];
0150 }
0151 if (maxNRecHitsPerChamber > maxRecHitsPerChamber) {
0152 return true;
0153 }
0154 }
0155 }
0156 }
0157 }
0158
0159 return false;
0160 }
0161
0162
0163 DEFINE_FWK_MODULE(HLTCSCAcceptBusyFilter);