CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTPixelActivityFilter.cc
Go to the documentation of this file.
2 
6 
8 
9 //
10 // class declaration
11 //
12 
14 public:
17  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
18 
19 private:
20  virtual bool hltFilter(edm::Event&, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterproduct) const override;
21 
22  edm::InputTag inputTag_; // input tag identifying product containing pixel clusters
23  unsigned int min_clusters_; // minimum number of clusters
24  unsigned int max_clusters_; // maximum number of clusters
26 
27 };
28 
35 
36 //
37 // constructors and destructor
38 //
39 
41  inputTag_ (config.getParameter<edm::InputTag>("inputTag")),
42  min_clusters_ (config.getParameter<unsigned int>("minClusters")),
43  max_clusters_ (config.getParameter<unsigned int>("maxClusters"))
44 {
45  inputToken_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(inputTag_);
46  LogDebug("") << "Using the " << inputTag_ << " input collection";
47  LogDebug("") << "Requesting at least " << min_clusters_ << " clusters";
48  if(max_clusters_ > 0)
49  LogDebug("") << "...but no more than " << max_clusters_ << " clusters";
50 }
51 
53 {
54 }
55 
56 void
60  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltSiPixelClusters"));
61  desc.add<unsigned int>("minClusters",3);
62  desc.add<unsigned int>("maxClusters",0);
63  descriptions.add("hltPixelActivityFilter",desc);
64 }
65 
66 //
67 // member functions
68 //
69 
70 // ------------ method called to produce the data ------------
72 {
73  // All HLT filters must create and fill an HLT filter object,
74  // recording any reconstructed physics objects satisfying (or not)
75  // this HLT filter, and place it in the Event.
76 
77  // The filter object
78  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
79 
80  // get hold of products from Event
82  event.getByToken(inputToken_, clusterColl);
83 
84  unsigned int clusterSize = clusterColl->dataSize();
85  LogDebug("") << "Number of clusters accepted: " << clusterSize;
86  bool accept = (clusterSize >= min_clusters_);
87  if(max_clusters_ > 0)
88  accept &= (clusterSize <= max_clusters_);
89 
90  // return with final filter decision
91  return accept;
92 }
93 
94 // define as a framework module
#define LogDebug(id)
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > inputToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
HLTPixelActivityFilter(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool saveTags() const
Definition: HLTFilter.h:45