CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HBHENoiseFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HBHENoiseFilter
4 // Class: HBHENoiseFilter
5 //
14 //
15 // Original Author: John Paul Chou (Brown)
16 //
17 //
18 
19 #include <iostream>
20 
22 
25 
26 //
27 // constants, enums and typedefs
28 //
29 
30 //
31 // static data member definitions
32 //
33 
34 //
35 // constructors and destructor
36 //
37 
39 {
40  //now do what ever initialization is needed
41  noiselabel_ = iConfig.getParameter<edm::InputTag>("noiselabel");
42  minRatio_ = iConfig.getParameter<double>("minRatio");
43  maxRatio_ = iConfig.getParameter<double>("maxRatio");
44  minHPDHits_ = iConfig.getParameter<int>("minHPDHits");
45  minRBXHits_ = iConfig.getParameter<int>("minRBXHits");
46  minHPDNoOtherHits_ = iConfig.getParameter<int>("minHPDNoOtherHits");
47  minZeros_ = iConfig.getParameter<int>("minZeros");
48  minHighEHitTime_ = iConfig.getParameter<double>("minHighEHitTime");
49  maxHighEHitTime_ = iConfig.getParameter<double>("maxHighEHitTime");
50  maxRBXEMF_ = iConfig.getParameter<double>("maxRBXEMF");
51  minNumIsolatedNoiseChannels_ = iConfig.getParameter<int>("minNumIsolatedNoiseChannels");
52  minIsolatedNoiseSumE_ = iConfig.getParameter<double>("minIsolatedNoiseSumE");
53  minIsolatedNoiseSumEt_ = iConfig.getParameter<double>("minIsolatedNoiseSumEt");
54  useTS4TS5_ = iConfig.getParameter<bool>("useTS4TS5");
55 
56  IgnoreTS4TS5ifJetInLowBVRegion_ = iConfig.getParameter<bool>("IgnoreTS4TS5ifJetInLowBVRegion");
57  jetlabel_ = iConfig.getParameter<edm::InputTag>("jetlabel");
58  maxjetindex_ = iConfig.getParameter<int>("maxjetindex");
59  maxNHF_ = iConfig.getParameter<double>("maxNHF");
60 
61 }
62 
63 
65 {
66 
67 }
68 
69 
70 //
71 // member functions
72 //
73 
74 // ------------ method called on each new Event ------------
75 bool
77 {
78  using namespace edm;
79 
80  // get the Noise summary object
82  iEvent.getByLabel(noiselabel_, summary_h);
83  if(!summary_h.isValid()) {
84  throw edm::Exception(edm::errors::ProductNotFound) << " could not find HcalNoiseSummary.\n";
85  return true;
86  }
87  const HcalNoiseSummary summary = *summary_h;
88 
89  // if(summary.HasBadRBXTS4TS5() == true) std::cout << "TS4TS5 rejection!" << std::endl;
90  // else std::cout << "TS4TS5 passing!" << std::endl;
91 
92 
93  // One HBHE region has HPD with low BV, which increases hits in HPD, and
94  // often causes otherwise good events to be flagged as noise.
95 
96  bool goodJetFoundInLowBVRegion=false; // checks whether a jet is in a low BV region, where false noise flagging rate is higher.
98  {
100  iEvent.getByLabel(jetlabel_, pfjet_h);
101  if (pfjet_h.isValid())
102  {
103  // Loop over all jets up to (and including) maxjetindex.
104  // If jet found in low-BV region, set goodJetFoundInLowBVRegion = true
105  int jetindex=0;
106  for( reco::PFJetCollection::const_iterator jet = pfjet_h->begin(); jet != pfjet_h->end(); ++jet)
107  {
108  if (jetindex>maxjetindex_) break; // only look at first N jets (N specified by user)
109  // Check whether jet is in low-BV region (0<eta<1.4, -1.8<phi<-1.4)
110  if (jet->eta()>0 && jet->eta()<1.4 &&
111  jet->phi()>-1.8 && jet->phi()<-1.4)
112  {
113  if (maxNHF_<0 || (maxNHF_>0 && jet->neutralHadronEnergyFraction()<maxNHF_))
114  {
115  //std::cout <<"JET eta = "<<jet->eta()<<" phi = "<<jet->phi()<<" Energy = "<<jet->energy()<<" neutral had fraction = "<<jet->neutralHadronEnergy()/jet->energy() <<std::endl;
116  goodJetFoundInLowBVRegion=true;
117  break;
118  }
119  }
120  ++jetindex;
121  }
122  } //if (pfjet_h.isValid())
123  else // no valid jet collection found
124  {
125  // If no jet collection found, do we want to throw a fatal exception? Or just proceed normally, not treating the lowBV region as special?
126  //throw edm::Exception(edm::errors::ProductNotFound) << " could not find PFJetCollection with label "<<jetlabel_<<".\n";
127  }
128  } // if (IgnoreTS4TS5ifJetInLowBVRegion_==true)
129 
130  if(summary.minE2Over10TS()<minRatio_) return false;
131  if(summary.maxE2Over10TS()>maxRatio_) return false;
132  if(summary.maxHPDHits()>=minHPDHits_) return false;
133  if(summary.maxRBXHits()>=minRBXHits_) return false;
134  if(summary.maxHPDNoOtherHits()>=minHPDNoOtherHits_) return false;
135  if(summary.maxZeros()>=minZeros_) return false;
136  if(summary.min25GeVHitTime()<minHighEHitTime_) return false;
137  if(summary.max25GeVHitTime()>maxHighEHitTime_) return false;
138  if(summary.minRBXEMF()<maxRBXEMF_) return false;
139  if(summary.numIsolatedNoiseChannels()>=minNumIsolatedNoiseChannels_) return false;
140  if(summary.isolatedNoiseSumE()>=minIsolatedNoiseSumE_) return false;
141  if(summary.isolatedNoiseSumEt()>=minIsolatedNoiseSumEt_) return false;
142  // Only use TS4TS5 test if jet is not in low BV region
143  if(useTS4TS5_ == true && summary.HasBadRBXTS4TS5() == true && goodJetFoundInLowBVRegion==false) return false;
144 
145  return true;
146 }
147 
148 //define this as a plug-in
T getParameter(std::string const &) const
virtual bool filter(edm::Event &, const edm::EventSetup &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::InputTag jetlabel_
int iEvent
Definition: GenABIO.cc:243
double minIsolatedNoiseSumE_
bool IgnoreTS4TS5ifJetInLowBVRegion_
edm::InputTag noiselabel_
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
HBHENoiseFilter(const edm::ParameterSet &)
double minIsolatedNoiseSumEt_
int minNumIsolatedNoiseChannels_