CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTHFAsymmetryFilter.cc
Go to the documentation of this file.
5 
7 {
8  HFHits_ = iConfig.getParameter<edm::InputTag>("HFHitCollection");
9  eCut_HF_ = iConfig.getParameter<double>("ECut_HF");
10  os_asym_ = iConfig.getParameter<double>("OS_Asym_max");
11  ss_asym_ = iConfig.getParameter<double>("SS_Asym_min");
12 }
13 
14 
16 {
17 }
18 
19 
20 // ------------ method called to produce the data ------------
21 bool
23 {
25 
26  iEvent.getByLabel(HFHits_,HFRecHitsH);
27 
28  double asym_hf_1 = -1;
29  double asym_hf_2 = -1;
30 
31  int n_hf[2] = {0,0};
32  double e_hf[2] = {0.,0.};
33 
34  //Select interesting HFRecHits
35  for (HFRecHitCollection::const_iterator it=HFRecHitsH->begin(); it!=HFRecHitsH->end(); it++) {
36  if (it->energy()>eCut_HF_)
37  {
38  if (it->id().zside() == -1)
39  {
40  ++n_hf[0];
41  e_hf[0] += it->energy();
42  }
43  else
44  {
45  ++n_hf[1];
46  e_hf[1] += it->energy();
47  }
48  }
49  }
50 
51 
52  for (int i=0;i<2;++i)
53  {
54  if (n_hf[i])
55  e_hf[i] /= n_hf[i];
56  }
57 
58  if (e_hf[0]+e_hf[1] != 0)
59  {
60  asym_hf_1 = e_hf[0]/(e_hf[0]+e_hf[1]);
61  asym_hf_2 = e_hf[1]/(e_hf[0]+e_hf[1]);
62  }
63  else
64  {
65  return false;
66  }
67 
68  bool pkam_1 = (asym_hf_1 >= ss_asym_ || asym_hf_1 <= os_asym_);
69  bool pkam_2 = (asym_hf_2 >= ss_asym_ || asym_hf_2 <= os_asym_);
70 
71  if (pkam_1 || pkam_2)
72  return true;
73 
74  return false;
75 }
76 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< T >::const_iterator const_iterator
int iEvent
Definition: GenABIO.cc:243
virtual bool filter(edm::Event &, const edm::EventSetup &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
HLTHFAsymmetryFilter(const edm::ParameterSet &)