CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTPixelActivityHFSumEnergyFilter.cc
Go to the documentation of this file.
2 
6 
7 //
8 // constructors and destructor
9 //
10 
12  inputTag_ (config.getParameter<edm::InputTag>("inputTag")),
13  HFHits_ (config.getParameter<edm::InputTag>("HFHitCollection")),
14  eCut_HF_ (config.getParameter<double>("eCut_HF")),
15  eMin_HF_ (config.getParameter<double>("eMin_HF")),
16  offset_ (config.getParameter<double>("offset")),
17  slope_ (config.getParameter<double>("slope"))
18 {
19  inputToken_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(inputTag_);
20  HFHitsToken_ = consumes<HFRecHitCollection>(HFHits_);
21 
22  LogDebug("") << "Using the " << inputTag_ << " input collection";
23 }
24 
26 {
27 }
28 
29 void
32  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltSiPixelClusters"));
33  desc.add<edm::InputTag>("HFHitCollection",edm::InputTag("hltHfreco"));
34  desc.add<double>("eCut_HF",0);
35  desc.add<double>("eMin_HF",10000.);
36  desc.add<double>("offset",-1000.);
37  desc.add<double>("slope",0.5);
38  descriptions.add("hltPixelActivityHFSumEnergyFilter",desc);
39 }
40 
41 //
42 // member functions
43 //
44 
45 // ------------ method called to produce the data ------------
47 {
48 
49  // get hold of products from Event
51  iEvent.getByToken(inputToken_, clusterColl);
52 
53  unsigned int clusterSize = clusterColl->dataSize();
54  LogDebug("") << "Number of clusters: " << clusterSize;
55 
57  iEvent.getByToken(HFHitsToken_,HFRecHitsH);
58 
59  double sumE = 0.;
60 
61  for (HFRecHitCollection::const_iterator it=HFRecHitsH->begin(); it!=HFRecHitsH->end(); it++) {
62  if (it->energy()>eCut_HF_) {
63  sumE += it->energy();
64  }
65  }
66 
67  bool accept = false;
68 
69  double thres = offset_ + slope_ * sumE;
70  double diff = clusterSize - thres; //diff = clustersize - (correlation line + offset)
71  if(sumE>eMin_HF_ && diff < 0.) accept = true;
72 
73  // return with final filter decision
74  return accept;
75 }
#define LogDebug(id)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< HFRecHit >::const_iterator const_iterator
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
edm::EDGetTokenT< HFRecHitCollection > HFHitsToken_
int iEvent
Definition: GenABIO.cc:230
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > inputToken_
HLTPixelActivityHFSumEnergyFilter(const edm::ParameterSet &)
virtual bool filter(edm::Event &iEvent, const edm::EventSetup &iSetup) override