CMS 3D CMS Logo

HLTCSCAcceptBusyFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTCSCAcceptBusyFilter
4 // Class: HLTCSCAcceptBusyFilter
5 //
13 //
14 // Original Author: Ingo Bloch
15 // Created: Mon Mar 15 11:39:08 CDT 2010
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
25 //include "FWCore/Framework/interface/EDFilter.h"
26 
28 
31 
35 
38 
40 
41 #include <string>
42 
43 //
44 // class declaration
45 //
46 
48 
49 public:
51  ~HLTCSCAcceptBusyFilter() override;
52  bool hltFilter(edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterproduct) const override;
53  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
54 
55 private:
57 
58  // ----------member data ---------------------------
61  bool invert;
62  unsigned int maxRecHitsPerChamber;
63 
64 };
65 
66 //
67 // constants, enums and typedefs
68 //
69 
70 //
71 // static data member definitions
72 //
73 
74 //
75 // constructors and destructor
76 //
78 {
79  //now do what ever initialization is needed
80  cscrechitsTag = iConfig.getParameter<edm::InputTag>("cscrechitsTag");
81  invert = iConfig.getParameter<bool>("invert");
82  maxRecHitsPerChamber = iConfig.getParameter<unsigned int>("maxRecHitsPerChamber");
83  cscrechitsToken = consumes<CSCRecHit2DCollection>(cscrechitsTag);
84 }
85 
86 
88 {
89 
90  // do anything here that needs to be done at desctruction time
91  // (e.g. close files, deallocate resources etc.)
92 
93 }
94 
95 
96 void
100  desc.add<edm::InputTag>("cscrechitsTag",edm::InputTag("hltCsc2DRecHits"));
101  desc.add<bool>("invert",true);
102  desc.add<unsigned int>("maxRecHitsPerChamber",200);
103  descriptions.add("hltCSCAcceptBusyFilter",desc);
104 }
105 
106 //
107 // member functions
108 //
109 
110 // ------------ method called on each new Event ------------
112 
113  using namespace edm;
114 
115  // Get the RecHits collection :
117  iEvent.getByToken(cscrechitsToken,recHits);
118 
120  return (!invert);
121  } else {
122  return ( invert);
123  }
124 
125 }
126 
127 
128 // ------------ method to find chamber with nMax hits
130 
131  unsigned int maxNRecHitsPerChamber(0);
132 
133  const unsigned int nEndcaps(2);
134  const unsigned int nStations(4);
135  const unsigned int nRings(4);
136  const unsigned int nChambers(36);
137  unsigned int allRechits[nEndcaps][nStations][nRings][nChambers];
138  for(auto & allRechit : allRechits){
139  for(unsigned int iS = 0;iS<nStations;++iS){
140  for(unsigned int iR = 0;iR<nRings;++iR){
141  for(unsigned int iC = 0;iC<nChambers;++iC){
142  allRechit[iS][iR][iC] = 0;
143  }
144  }
145  }
146  }
147 
148  for(auto const & it : *recHits) {
149  ++allRechits
150  [it.cscDetId().endcap()-1]
151  [it.cscDetId().station()-1]
152  [it.cscDetId().ring()-1]
153  [it.cscDetId().chamber()-1];
154  }
155 
156  for(auto & allRechit : allRechits){
157  for(unsigned int iS = 0;iS<nStations;++iS){
158  for(unsigned int iR = 0;iR<nRings;++iR){
159  for(unsigned int iC = 0;iC<nChambers;++iC){
160  if(allRechit[iS][iR][iC] > maxNRecHitsPerChamber) {
161  maxNRecHitsPerChamber = allRechit[iS][iR][iC];
162  }
163  if(maxNRecHitsPerChamber > maxRecHitsPerChamber) {
164  return true;
165  }
166  }
167  }
168  }
169  }
170 
171  return false;
172 
173 }
174 
175 // declare this class as a framework plugin
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool AcceptManyHitsInChamber(unsigned int maxRecHitsPerChamber, const edm::Handle< CSCRecHit2DCollection > &recHits) const
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
static const int nEndcaps
Definition: GenABIO.cc:115
edm::EDGetTokenT< CSCRecHit2DCollection > cscrechitsToken
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
HLTCSCAcceptBusyFilter(const edm::ParameterSet &)