CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
HFStripFilter Class Reference

#include <HFStripFilter.h>

Public Member Functions

 HFStripFilter (double stripThreshold, double maxThreshold, double timeMax, double maxStripTime, double wedgeCut, int seedHitIetaMax, int gap, int lstrips, int verboseLevel)
 
void runFilter (HFRecHitCollection &rec, const HcalChannelQuality *myqual) const
 
 ~HFStripFilter ()
 

Static Public Member Functions

static edm::ParameterSetDescription fillDescription ()
 
static std::unique_ptr< HFStripFilterparseParameterSet (const edm::ParameterSet &ps)
 

Private Attributes

int gap_
 
int lstrips_
 
double maxStripTime_
 
double maxThreshold_
 
int seedHitIetaMax_
 
double stripThreshold_
 
double timeMax_
 
int verboseLevel_
 
double wedgeCut_
 

Detailed Description

Definition at line 10 of file HFStripFilter.h.

Constructor & Destructor Documentation

HFStripFilter::HFStripFilter ( double  stripThreshold,
double  maxThreshold,
double  timeMax,
double  maxStripTime,
double  wedgeCut,
int  seedHitIetaMax,
int  gap,
int  lstrips,
int  verboseLevel 
)

Definition at line 10 of file HFStripFilter.cc.

References verboseLevel_.

20  gap_(gap),
23 {
24  // For the description of CMSSW message logging, see
25  // https://twiki.cern.ch/twiki/bin/view/CMSPublic/SWGuideMessageLogger
26  if (verboseLevel_ >= 20)
27  edm::LogInfo("HFStripFilter") << "constructor called";
28 }
double stripThreshold_
Definition: HFStripFilter.h:33
double maxThreshold_
Definition: HFStripFilter.h:34
double wedgeCut_
Definition: HFStripFilter.h:37
double maxStripTime_
Definition: HFStripFilter.h:36
HFStripFilter::~HFStripFilter ( )

Definition at line 30 of file HFStripFilter.cc.

References verboseLevel_.

31 {
32  if (verboseLevel_ >= 20)
33  edm::LogInfo("HFStripFilter") << "destructor called";
34 }

Member Function Documentation

edm::ParameterSetDescription HFStripFilter::fillDescription ( )
static

Definition at line 343 of file HFStripFilter.cc.

References edm::ParameterSetDescription::add(), and edm::ParameterSetDescription::addUntracked().

Referenced by HFPhase1Reconstructor::fillDescriptions().

344 {
346 
347  desc.add<double>("stripThreshold", 40.0);
348  desc.add<double>("maxThreshold", 100.0);
349  desc.add<double>("timeMax", 6.0);
350  desc.add<double>("maxStripTime", 10.0);
351  desc.add<double>("wedgeCut", 0.05);
352  desc.add<int>("seedHitIetaMax", 35);
353  desc.add<int>("gap", 2);
354  desc.add<int>("lstrips", 2);
355  desc.addUntracked<int>("verboseLevel", 0);
356 
357  return desc;
358 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::unique_ptr< HFStripFilter > HFStripFilter::parseParameterSet ( const edm::ParameterSet ps)
static

Definition at line 327 of file HFStripFilter.cc.

References edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

Referenced by HFPhase1Reconstructor::HFPhase1Reconstructor().

329 {
330  return std::make_unique<HFStripFilter>(
331  ps.getParameter<double>("stripThreshold"),
332  ps.getParameter<double>("maxThreshold"),
333  ps.getParameter<double>("timeMax"),
334  ps.getParameter<double>("maxStripTime"),
335  ps.getParameter<double>("wedgeCut"),
336  ps.getParameter<int>("seedHitIetaMax"),
337  ps.getParameter<int>("gap"),
338  ps.getParameter<int>("lstrips"),
339  ps.getUntrackedParameter<int>("verboseLevel")
340  );
341 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void HFStripFilter::runFilter ( HFRecHitCollection rec,
const HcalChannelQuality myqual 
) const

Definition at line 36 of file HFStripFilter.cc.

References funct::abs(), HcalTopology::decIPhi(), HcalDetId::depth(), CaloRecHit::energy(), gap_, HcalForward, HcalPhase1FlagLabels::HFAnomalousHit, triggerObjects_cff::id, HFRecHit::id(), HcalDetId::ieta(), HcalTopology::incIPhi(), createfilelist::int, HcalDetId::iphi(), LogDebug, lstrips_, maxStripTime_, maxThreshold_, seedHitIetaMax_, CaloRecHit::setEnergy(), stripThreshold_, CaloRecHit::time(), timeMax_, HcalCondObjectContainerBase::topo(), mitigatedMETSequence_cff::U, verboseLevel_, and wedgeCut_.

38 {
39  if (verboseLevel_ >= 20)
40  edm::LogInfo("HFStripFilter") << "runFilter called";
41 
42  std::vector<HFRecHit> d1strip;
43  std::vector<HFRecHit> d2strip;
44 
45  HFRecHit d1max;
46  HFRecHit d2max;
47 
48  d1max.setEnergy(-10);
49  // find d1 and d2 seed hits with max energy and time < timeMax_
50  for (auto& it : rec) {
51  if (it.time() > timeMax_ || it.time() < 0 || it.energy() < stripThreshold_) continue;
52  // find HF hit with maximum signal in depth = 1
53  if (it.id().depth() == 1) {
54  if (it.energy() > d1max.energy() && std::abs(it.id().ieta()) < seedHitIetaMax_) {
55  d1max = it;
56  }
57  }
58  }
59 
60  bool signStripIeta = false;
61  int stripIphiMax = 0;
62  int stripIetaMax = 0;
63  int stripDepthMax = 0;
64 
65  if (d1max.energy() > 0) {
66  d1strip.push_back(d1max);
67  signStripIeta = std::signbit(d1max.id().ieta());
68  stripIphiMax = d1max.id().iphi();
69  stripIetaMax = d1max.id().ieta();
70  stripDepthMax = d1max.id().depth();
71  }
72 
73  d2max.setEnergy(-10);
74  for (auto& it : rec) {
75  if (it.time() > timeMax_ || it.time() < 0 || it.energy() < stripThreshold_) continue;
76  // find HFhit with maximum signal in depth = 2
77  if (it.id().depth() == 2) {
78  if (it.energy() > d2max.energy() && std::abs(it.id().ieta()) < seedHitIetaMax_) {
79  if (d1max.energy() > 0) {
80  bool signIeta = std::signbit(it.id().ieta());
81  if (it.id().iphi() == stripIphiMax && signIeta == signStripIeta) {
82  d2max = it;
83  }
84  } else {
85  d2max = it;
86  }
87  }
88  }
89  }
90 
91  if (d2max.energy() > 0) d2strip.push_back(d2max);
92 
93  // check if at least one of possible seed hits has enough energy
94  if (d1max.energy() < maxThreshold_ && d2max.energy() < maxThreshold_) return;
95 
96  if (d1max.energy() <= 0 && d2max.energy() > 0) {
97  signStripIeta = std::signbit(d2max.id().ieta());
98  stripIphiMax = d2max.id().iphi();
99  stripIetaMax = d2max.id().ieta();
100  stripDepthMax = d2max.id().depth();
101  }
102 
103  std::stringstream ss;
104  if (verboseLevel_ >= 30) {
105  if (d1max.energy() > 0) {
106  ss << " MaxHit in Depth 1: ieta = " << d1max.id().ieta() << " iphi = " << stripIphiMax
107  << " energy = " << d1max.energy() << " time = " << d1max.time() << std::endl;
108  }
109  if (d2max.energy() > 0) {
110  ss << " MaxHit in Depth 2: ieta = " << d2max.id().ieta() << " iphi = " << stripIphiMax
111  << " energy = " << d2max.energy() << " time = " << d2max.time() << std::endl;
112  }
113  ss << " stripThreshold_ = " << stripThreshold_ << std::endl;
114  }
115 
116  // prepare the strips: all hits along given ieta in one wedge (d1strip and d2strip)
117  //---------------------------------------------------------------------------------
118  for (auto& it : rec) {
119  if (it.energy() < stripThreshold_) continue;
120  bool signIeta = std::signbit(it.id().ieta());
121 
122  if (verboseLevel_ >= 30) {
123  ss << " HF hit: ieta = " << it.id().ieta() << "\t iphi = " << it.id().iphi()
124  << "\t depth = " << it.id().depth() << "\t time = " << it.time() << "\t energy = "
125  << it.energy() << "\t flags = " << it.flags() << std::endl;
126  }
127 
128  // collect hits with the same iphi but different ieta into strips
129  if (it.id().iphi() == stripIphiMax && signIeta == signStripIeta && it.time() < maxStripTime_) {
130  if (it.id().depth() == 1) {
131  // check if hit = it is already in d1strip
132  bool pass = false;
133  if (d1strip.empty()) {
134  if (std::abs(it.id().iphi() - stripIetaMax) <= gap_) {
135  d1strip.push_back(it);
136  pass = true;
137  }
138  } else {
139  for (auto& it1 : d1strip) {
140  if (it.id().ieta() == it1.id().ieta() && it.energy() == it1.energy()) {
141  pass = true;
142  break;
143  }
144  }
145  }
146  if (pass) continue;
147  // add hit = it to d1strip if distance to the closest hit in d1strip <= gap_
148  for (auto& it1 : d1strip) {
149  // check distance along Ieta to the closest hit
150  if (std::abs(it1.id().ieta() - it.id().ieta()) <= gap_) {
151  d1strip.push_back(it);
152  break;
153  }
154  }
155  }
156  else if (it.id().depth() == 2) {
157  // check if hit = it is already in d2strip
158  bool pass= false;
159  if (d2strip.empty()) {
160  if (std::abs(it.id().ieta() - stripIetaMax) <= gap_) {
161  d2strip.push_back(it);
162  pass = true;
163  }
164  } else {
165  for (auto& it1 : d2strip) {
166  if (it.id().ieta() == it1.id().ieta() && it.energy() == it1.energy()) {
167  pass = true;
168  break;
169  }
170  }
171  }
172  if (pass) continue;
173 
174  // add hit = it to d2strip if distance to the closest hit in d1strip <= gap_
175  for (auto& it1 : d2strip) {
176  // check distance along Ieta to the closest hit
177  if (std::abs(it1.id().ieta() - it.id().ieta()) <= gap_) {
178  d2strip.push_back(it);
179  break;
180  }
181  }
182  }
183  }
184  }
185 
186  if (verboseLevel_ >= 30) {
187  ss << " Lstrip1 = " << (int)d1strip.size() << " (iphi = " << stripIphiMax
188  << ") Lstrip2 = " << (int)d2strip.size() << std::endl << " Strip1: ";
189  for (auto& it1 : d1strip) {
190  ss << it1.energy() << " (" << it1.id().ieta() << ") ";
191  }
192  ss << std::endl << " Strip2: ";
193  for (auto& it1 : d2strip) {
194  ss << it1.energy() << " (" << it1.id().ieta() << ") ";
195  }
196  ss << std::endl;
197  }
198 
199  // check if at least one of strips in depth1 or depth2 is not less than lstrips_
200  if ((int)d1strip.size() < lstrips_ && (int)d2strip.size() < lstrips_) {
201  if (verboseLevel_ >= 30) {
202  LogDebug("HFSFilter") << ss.str();
203  }
204  return;
205  }
206 
207  // define range of strips in ieta
208  int ietaMin1 = 1000; // for d1strip
209  int ietaMax1 = -1000;
210  for (auto& it1 : d1strip) {
211  if (it1.id().ieta() < ietaMin1) ietaMin1 = it1.id().ieta();
212  if (it1.id().ieta() > ietaMax1) ietaMax1 = it1.id().ieta();
213  }
214  int ietaMin2 = 1000; // for d2strip
215  int ietaMax2 = -1000;
216  for (auto& it1 : d2strip) {
217  if (it1.id().ieta() < ietaMin2) ietaMin2 = it1.id().ieta();
218  if (it1.id().ieta() > ietaMax2) ietaMax2 = it1.id().ieta();
219  }
220 
221  // define ietamin and ietamax - common area for d1strip and d2strip
222  int ietaMin = ietaMin1;
223  int ietaMax = ietaMax1;
224 
225  if (ietaMin2 >= ietaMin1 && ietaMin2 <= ietaMax1) {
226  if (ietaMax2 > ietaMax1) ietaMax = ietaMax2;
227  } else if (ietaMin2 <= ietaMin1 && ietaMax2 >= ietaMin1) {
228  if (ietaMin2 < ietaMin1) ietaMin = ietaMin2;
229  if (ietaMax2 > ietaMax1) ietaMax = ietaMax2;
230  }
231 
232  // calculate the total energy in strips
233  double eStrip = 0;
234  for (auto& it1 : d1strip) {
235  if (it1.id().ieta() < ietaMin || it1.id().ieta() > ietaMax) continue;
236  eStrip += it1.energy();
237  }
238  for (auto& it1 : d2strip) {
239  if (it1.id().ieta() < ietaMin || it1.id().ieta() > ietaMax) continue;
240  eStrip += it1.energy();
241  }
242 
243  if (verboseLevel_ >= 30) {
244  ss << " ietaMin1 = " << ietaMin1 << " ietaMax1 = " << ietaMax1 << std::endl
245  << " ietaMin2 = " << ietaMin2 << " ietaMax2 = " << ietaMax2 << std::endl
246  << " Common strip: ietaMin = " << ietaMin << " ietaMax = " << ietaMax << std::endl;
247  }
248 
249  // energies in the neighboring wedges
250  double energyIphi1 = 0;
251  double energyIphi2 = 0;
252  int iphi1 = 0, iphi2 = 0;
253 
254  // get information about the neighboring wedges from HcalTopology
255  HcalDetId neighbour;
256 
257  for (auto& it : rec) {
258  if (it.energy() < stripThreshold_) continue;
259  if (it.id().ieta() < ietaMin || it.id().ieta() > ietaMax) continue;
260  HcalDetId id(HcalForward, it.id().ieta(), stripIphiMax, stripDepthMax);
261  bool neigh = myqual->topo()->decIPhi(id, neighbour);
262  iphi1 = (neigh) ? neighbour.iphi() : 0;
263  neigh = myqual->topo()->incIPhi(id, neighbour);
264  iphi2 = (neigh) ? neighbour.iphi() : 0;
265  if (it.id().iphi() == iphi1) energyIphi1 += it.energy(); // Energy iphi1
266  else if (it.id().iphi() == iphi2) energyIphi2 += it.energy(); // Energy iphi2
267  }
268 
269  double ratio1 = eStrip > 0 ? energyIphi1/eStrip : 0;
270  double ratio2 = eStrip > 0 ? energyIphi2/eStrip : 0;
271 
272  if (verboseLevel_ >= 30) {
273  ss << " iphi = " << stripIphiMax << " iphi1 = " << iphi1 << " iphi2 = "
274  << iphi2 << std::endl
275  << " Estrip = " << eStrip << " EnergyIphi1 = " << energyIphi1
276  << " Ratio = " << ratio1 << std::endl
277  << " " << " EnergyIphi2 = " << energyIphi2
278  << " Ratio = " << ratio2 << std::endl;
279  }
280 
281  // check if our wedge does not have substantial leak into adjacent wedges
282  if (ratio1 < wedgeCut_ && ratio2 < wedgeCut_) { // noise event with strips (d1 and/or d2)
283 
284  if (verboseLevel_ >= 30) {
285  ss << " stripPass = false" << std::endl
286  << " mark hits in strips now " << std::endl;
287  }
288 
289  // Figure out which rechits need to be tagged (d1strip)
290  for (auto& it1 : d1strip) {
291  if (it1.id().ieta() < ietaMin || it1.id().ieta() > ietaMax) continue;
292  HFRecHitCollection::iterator hit = rec.find(it1.id());
293  if (hit != rec.end()) {
294  // tag a rechit with the anomalous hit flag
295  if (verboseLevel_ >= 30) {
296  ss << " d1strip: marked hit = " << (*hit) << std::endl;
297  }
298  hit->setFlagField(1U, HcalPhase1FlagLabels::HFAnomalousHit);
299  }
300  }
301 
302  // Figure out which rechits need to be tagged (dstrip)
303  for (auto& it1 : d2strip) {
304  if (it1.id().ieta() < ietaMin || it1.id().ieta() > ietaMax) continue;
305  HFRecHitCollection::iterator hit = rec.find(it1.id());
306  if (hit != rec.end()) {
307  // tag a rechit with the anomalous hit flag
308  if (verboseLevel_ >= 30) {
309  ss << " d2strip: marked hit = " << (*hit) << std::endl;
310  }
311  hit->setFlagField(1U, HcalPhase1FlagLabels::HFAnomalousHit);
312  }
313  }
314  if (verboseLevel_ >= 30) {
315  LogDebug("HFSFilter") << ss.str();
316  }
317  }
318  else {
319  if (verboseLevel_ >= 30) {
320  ss << " stripPass = true" << std::endl;
321  LogDebug("HFSFilter") << ss.str();
322  }
323  }
324 }
#define LogDebug(id)
constexpr float energy() const
Definition: CaloRecHit.h:31
bool decIPhi(const HcalDetId &id, HcalDetId &neighbor) const
int depth() const
get the tower depth
Definition: HcalDetId.h:166
double stripThreshold_
Definition: HFStripFilter.h:33
bool incIPhi(const HcalDetId &id, HcalDetId &neighbor) const
constexpr float time() const
Definition: CaloRecHit.h:33
int ieta() const
get the cell ieta
Definition: HcalDetId.h:159
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double maxThreshold_
Definition: HFStripFilter.h:34
std::vector< T >::iterator iterator
int iphi() const
get the cell iphi
Definition: HcalDetId.h:161
double wedgeCut_
Definition: HFStripFilter.h:37
HcalDetId id() const
Definition: HFRecHit.h:31
double maxStripTime_
Definition: HFStripFilter.h:36
constexpr void setEnergy(float energy)
Definition: CaloRecHit.h:32
const HcalTopology * topo() const

Member Data Documentation

int HFStripFilter::gap_
private

Definition at line 39 of file HFStripFilter.h.

Referenced by runFilter().

int HFStripFilter::lstrips_
private

Definition at line 40 of file HFStripFilter.h.

Referenced by runFilter().

double HFStripFilter::maxStripTime_
private

Definition at line 36 of file HFStripFilter.h.

Referenced by runFilter().

double HFStripFilter::maxThreshold_
private

Definition at line 34 of file HFStripFilter.h.

Referenced by runFilter().

int HFStripFilter::seedHitIetaMax_
private

Definition at line 38 of file HFStripFilter.h.

Referenced by runFilter().

double HFStripFilter::stripThreshold_
private

Definition at line 33 of file HFStripFilter.h.

Referenced by runFilter().

double HFStripFilter::timeMax_
private

Definition at line 35 of file HFStripFilter.h.

Referenced by runFilter().

int HFStripFilter::verboseLevel_
private

Definition at line 41 of file HFStripFilter.h.

Referenced by HFStripFilter(), runFilter(), and ~HFStripFilter().

double HFStripFilter::wedgeCut_
private

Definition at line 37 of file HFStripFilter.h.

Referenced by runFilter().