CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HBHEIsolatedNoiseAlgos.cc
Go to the documentation of this file.
1 /*
2 Description: Isolation algorithms used to identify anomalous noise in the HB/HE.
3  These algorithms will be used to reflag HB/HE rechits as noise.
4 
5  There are 4 objects implemented here:
6  1) ObjectValidator
7  2) PhysicsTowerOrganizer
8  3) HBHEHitMap
9  4) HBHEHitMapOrganizer
10  See comments below for details.
11 
12 Original Author: John Paul Chou (Brown University)
13  Thursday, September 2, 2010
14 */
15 
17 
20 
27 
35 
37 //
38 // ObjectValidator
39 //
41 
43 {
44  HBThreshold_ = iConfig.getParameter<double>("HBThreshold");
45  HESThreshold_ = iConfig.getParameter<double>("HESThreshold");
46  HEDThreshold_ = iConfig.getParameter<double>("HEDThreshold");
47  EBThreshold_ = iConfig.getParameter<double>("EBThreshold");
48  EEThreshold_ = iConfig.getParameter<double>("EEThreshold");
49 
50  HcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("HcalAcceptSeverityLevel");
51  EcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("EcalAcceptSeverityLevel");
52  UseHcalRecoveredHits_ = iConfig.getParameter<bool>("UseHcalRecoveredHits");
53  UseEcalRecoveredHits_ = iConfig.getParameter<bool>("UseEcalRecoveredHits");
54 
55  MinValidTrackPt_ = iConfig.getParameter<double>("MinValidTrackPt");
56  MinValidTrackPtBarrel_ = iConfig.getParameter<double>("MinValidTrackPtBarrel");
57  MinValidTrackNHits_ = iConfig.getParameter<int>("MinValidTrackNHits");
58 
65 
66  return;
67 }
68 
70 {
71 }
72 
74 {
76 
77  // require the hit to pass a certain energy threshold
78  if(hit.id().subdet()==HcalBarrel && hit.energy()<HBThreshold_) return false;
79  else if(hit.id().subdet()==HcalEndcap && hit.id().ietaAbs()<=20 && hit.energy()<HESThreshold_) return false;
80  else if(hit.id().subdet()==HcalEndcap && hit.id().ietaAbs()>20 && hit.energy()<HEDThreshold_) return false;
81 
82  // determine if the hit is good, bad, or recovered
83  const DetId id = hit.detid();
84  const uint32_t recHitFlag = hit.flags();
85  const uint32_t dbStatusFlag = theHcalChStatus_->getValues(id)->getValue();
86  int severityLevel = theHcalSevLvlComputer_->getSeverityLevel(id, recHitFlag, dbStatusFlag);
87  bool isRecovered = theHcalSevLvlComputer_->recoveredRecHit(id, recHitFlag);
88 
89  if(severityLevel==0) return true;
90  if(isRecovered) return UseHcalRecoveredHits_;
91  if(severityLevel>static_cast<int>(HcalAcceptSeverityLevel_)) return false;
92  else return true;
93 }
94 
96 {
98 
99  // require the hit to pass a certain energy threshold
100  const DetId id = hit.detid();
101  if(id.subdetId()==EcalBarrel && hit.energy()<EBThreshold_) return false;
102  else if(id.subdetId()==EcalEndcap && hit.energy()<EEThreshold_) return false;
103 
104  // determine if the hit is good, bad, or recovered
105  int severityLevel = 999;
106  if (id.subdetId() == EcalBarrel && theEBRecHitCollection_!=0) severityLevel = theEcalSevLvlAlgo_->severityLevel(hit);//id, *theEBRecHitCollection_, *theEcalChStatus_, 5., EcalSeverityLevelAlgo::kSwissCross, 0.95, 2., 15., 0.999);
107  else if(id.subdetId() == EcalEndcap && theEERecHitCollection_!=0) severityLevel = theEcalSevLvlAlgo_->severityLevel(hit);//id, *theEERecHitCollection_, *theEcalChStatus_, 5., EcalSeverityLevelAlgo::kSwissCross, 0.95, 2., 15., 0.999);
108  else return false;
109 
110  if(severityLevel == EcalSeverityLevel::kGood) return true;
111  if(severityLevel == EcalSeverityLevel::kRecovered) return UseEcalRecoveredHits_;
112  if(severityLevel > static_cast<int>(EcalAcceptSeverityLevel_)) return false;
113  else return true;
114 }
115 
117 {
118  if(trk.pt()<MinValidTrackPt_) return false;
119  if(trk.pt()<MinValidTrackPtBarrel_ && std::fabs(trk.momentum().eta())<1.479) return false;
120  if(trk.numberOfValidHits()<MinValidTrackNHits_) return false;
121  return true;
122 }
123 
125 //
126 // PhysicsTowerOrganizer
127 //
129 
131  const edm::EventSetup& evSetup,
132  const edm::Handle<HBHERecHitCollection>& hbhehitcoll_h,
133  const edm::Handle<EcalRecHitCollection>& ebhitcoll_h,
134  const edm::Handle<EcalRecHitCollection>& eehitcoll_h,
135  const edm::Handle<std::vector<reco::TrackExtrapolation> >& trackextrapcoll_h,
136  const ObjectValidatorAbs& objectvalidator,
137  const CaloTowerConstituentsMap& ctcm)
138 {
139  // get some geometries
141  evSetup.get<CaloGeometryRecord>().get(pG);
142  const CaloGeometry* geo = pG.product();
145 
146  // do the HCAL hits
147  for(HBHERecHitCollection::const_iterator it=hbhehitcoll_h->begin(); it!=hbhehitcoll_h->end(); ++it) {
148  const HBHERecHit* hit=&(*it);
149 
150  // check that the hit is valid
151  if(!objectvalidator.validHit(*hit)) continue;
152 
153  // add the hit to the organizer
154  CaloTowerDetId tid = ctcm.towerOf(hit->id());
155  insert_(tid, hit);
156  }
157 
158  // do the EB hits
159  for(EcalRecHitCollection::const_iterator it=ebhitcoll_h->begin(); it!=ebhitcoll_h->end(); ++it) {
160  const EcalRecHit* hit=&(*it);
161 
162  if(!objectvalidator.validHit(*hit)) continue;
163  CaloTowerDetId tid = ctcm.towerOf(hit->id());
164  insert_(tid, hit);
165  }
166 
167  // do the EE hits
168  for(EcalRecHitCollection::const_iterator it=eehitcoll_h->begin(); it!=eehitcoll_h->end(); ++it) {
169  const EcalRecHit* hit=&(*it);
170 
171  if(!objectvalidator.validHit(*hit)) continue;
172  CaloTowerDetId tid = ctcm.towerOf(hit->id());
173  insert_(tid, hit);
174  }
175 
176  // do the tracks
177  for(std::vector<reco::TrackExtrapolation>::const_iterator it=trackextrapcoll_h->begin(); it!=trackextrapcoll_h->end(); ++it) {
178  const reco::TrackExtrapolation* extrap=&(*it);
179  const reco::Track* track = &(*(extrap->track()));
180 
181  // validate track
182  if(!objectvalidator.validTrack(*track)) continue;
183 
184  // get the point
185  if ( extrap->positions().size()==0 ) continue;
186  const GlobalPoint point(extrap->positions().front().x(),
187  extrap->positions().front().y(),
188  extrap->positions().front().z());
189 
190 
191  if(std::fabs(point.eta())<1.479) {
192  EBDetId cell = gEB->getClosestCell(point);
193  CaloTowerDetId tid = ctcm.towerOf(cell);
194  insert_(tid, track);
195  } else {
196  EEDetId cell = gEE->getClosestCell(point);
197  CaloTowerDetId tid = ctcm.towerOf(cell);
198  insert_(tid, track);
199  }
200  }
201 
202  return;
203 }
204 
206 {
207  // create dummy PhysicsTower
208  PhysicsTower dummy;
209 
210  // correct for the merging of the |ieta|=28-29 towers
211  if(id.ietaAbs()==29) dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
212  else dummy.id=id;
213 
214  // search on the dummy
215  std::set<PhysicsTower, towercmp>::iterator it=towers_.find(dummy);
216 
217  if(it==towers_.end()) return 0;
218 
219  // for whatever reason, I can't get a non-const out of the find method
220  PhysicsTower &twr = const_cast<PhysicsTower&>(*it);
221  return &twr;
222 }
223 
225 {
226  // create dummy PhysicsTower
227  PhysicsTower dummy;
228 
229  // correct for the merging of the |ieta|=28-29 towers
230  if(id.ietaAbs()==29) dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
231  else dummy.id=id;
232 
233  // search on the dummy
234  std::set<PhysicsTower, towercmp>::iterator it=towers_.find(dummy);
235 
236  if(it==towers_.end()) return 0;
237  return &(*it);
238 }
239 
240 const PhysicsTower* PhysicsTowerOrganizer::findTower(int ieta, int iphi) const
241 {
242  CaloTowerDetId tid(ieta, iphi);
243  return findTower(tid);
244 }
245 
247 {
248  CaloTowerDetId tid(ieta, iphi);
249  return findTower(tid);
250 }
251 
252 void PhysicsTowerOrganizer::findNeighbors(const CaloTowerDetId& tempid, std::set<const PhysicsTower*>& neighbors) const
253 {
254  // correct for the merging of the |ieta|=28-29 towers
255  CaloTowerDetId id(tempid);
256  if(tempid.ietaAbs()==29) id = CaloTowerDetId((tempid.ietaAbs()-1)*tempid.zside(), tempid.iphi());
257 
258  std::vector<CaloTowerDetId> ids;
259  // get the neighbor with higher iphi
260  if(id.ietaAbs()<=20) {
261  if(id.iphi()==72) ids.push_back(CaloTowerDetId(id.ieta(), 1));
262  else ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()+1));
263  } else {
264  if(id.iphi()==71) ids.push_back(CaloTowerDetId(id.ieta(), 1));
265  else ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()+2));
266  }
267 
268  // get the neighbor with the lower iphi
269  if(id.ietaAbs()<=20) {
270  if(id.iphi()==1) ids.push_back(CaloTowerDetId(id.ieta(), 72));
271  else ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()-1));
272  } else {
273  if(id.iphi()==1) ids.push_back(CaloTowerDetId(id.ieta(), 71));
274  else ids.push_back(CaloTowerDetId(id.ieta(), id.iphi()-2));
275  }
276 
277  // get the neighbor with the higher ietaAbs
278  if(id.ietaAbs()==20 && (id.iphi()%2)==0)
279  ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()-1));
280  else
281  ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()));
282 
283  // get the neighbor(s) with the lower ietaAbs
284  if(id.ietaAbs()==21) {
285  ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()));
286  ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()+1));
287  } else if(id.ietaAbs()==1) {
288  ids.push_back(CaloTowerDetId(-id.ieta(), id.iphi()));
289  } else {
290  ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()));
291  }
292 
293  // get the neighbor with higher ieta and higher iphi
294  if(id.ietaAbs()<=19 || (id.ietaAbs()==20 && (id.iphi()%2)==0)) {
295  if(id.iphi()==72) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 1));
296  else ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()+1));
297  } else if(id.ietaAbs()>=21) {
298  if(id.iphi()==71) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 1));
299  else ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()+2));
300  }
301 
302  // get the neighbor with higher ieta and lower iphi
303  if(id.ietaAbs()<=19) {
304  if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 72));
305  else ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()-1));
306  } else if(id.ietaAbs()>=21 || (id.ietaAbs()==20 && (id.iphi()%2)==1)) {
307  if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), 71));
308  else ids.push_back(CaloTowerDetId((id.ietaAbs()+1)*id.zside(), id.iphi()-2));
309  }
310 
311  // get the neighbor with lower ieta and higher iphi
312  if(id.ietaAbs()==1) {
313  if(id.iphi()==72) ids.push_back(CaloTowerDetId(-id.ieta(), 1));
314  else ids.push_back(CaloTowerDetId(-id.ieta(), id.iphi()+1));
315  } else if(id.ietaAbs()<=20) {
316  if(id.iphi()==72) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 1));
317  else ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()+1));
318  } else if(id.ietaAbs()>=21) {
319  if(id.iphi()==71) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 1));
320  else ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()+2));
321  }
322 
323  // get the neighbor with lower ieta and lower iphi
324  if(id.ietaAbs()==1) {
325  if(id.iphi()==1) ids.push_back(CaloTowerDetId(-id.ieta(), 72));
326  else ids.push_back(CaloTowerDetId(-id.ieta(), id.iphi()-1));
327  } else if(id.ietaAbs()<=20) {
328  if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 72));
329  else ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()-1));
330  } else if(id.ietaAbs()>=22) {
331  if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 71));
332  else ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()-2));
333  } else if(id.ietaAbs()==21) {
334  if(id.iphi()==1) ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), 72));
335  else ids.push_back(CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi()-1));
336  }
337 
338  // clear neighbors
339  neighbors.clear();
340 
341  // find the neighbors and add them to the eponymous set
342  for(std::vector<CaloTowerDetId>::const_iterator it=ids.begin(); it!=ids.end(); ++it) {
343  const PhysicsTower* twr=findTower(*it);
344  if(twr) neighbors.insert(twr);
345  }
346 
347  return;
348 }
349 
350 void PhysicsTowerOrganizer::findNeighbors(const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors) const
351 {
352  findNeighbors(twr->id, neighbors);
353  return;
354 }
355 
356 void PhysicsTowerOrganizer::findNeighbors(int ieta, int iphi, std::set<const PhysicsTower*>& neighbors) const
357 {
358  findNeighbors(CaloTowerDetId(ieta, iphi), neighbors);
359  return;
360 }
361 
363 {
364  PhysicsTower* twr=findTower(id);
365  if(twr==0) {
366  PhysicsTower dummy;
367  if(id.ietaAbs()==29)
368  dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
369  else
370  dummy.id = id;
371  dummy.hcalhits.insert(hit);
372  towers_.insert(dummy);
373  } else {
374  twr->hcalhits.insert(hit);
375  }
376  return;
377 }
378 
380 {
381  PhysicsTower* twr=findTower(id);
382  if(twr==0) {
383  PhysicsTower dummy;
384  if(id.ietaAbs()==29)
385  dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
386  else
387  dummy.id = id;
388  dummy.ecalhits.insert(hit);
389  towers_.insert(dummy);
390  } else {
391  twr->ecalhits.insert(hit);
392  }
393  return;
394 }
395 
397 {
398  PhysicsTower* twr=findTower(id);
399  if(twr==0) {
400  PhysicsTower dummy;
401  if(id.ietaAbs()==29)
402  dummy.id = CaloTowerDetId((id.ietaAbs()-1)*id.zside(), id.iphi());
403  else
404  dummy.id = id;
405  dummy.tracks.insert(track);
406  towers_.insert(dummy);
407  } else {
408  twr->tracks.insert(track);
409  }
410  return;
411 }
412 
413 
415 //
416 // HBHEHitMap
417 //
419 
421 {
423  nHits_=-999;
428 
429 }
430 
431 double HBHEHitMap::hitEnergy(void) const
432 {
433  if(hitEnergy_<-900) calcHits_();
434  return hitEnergy_;
435 }
436 
437 int HBHEHitMap::nHits(void) const
438 {
439  if(nHits_<-900) calcHits_();
440  return nHits_;
441 }
442 
444 {
445  if(hitEnergyTrkFid_<-900) calcHits_();
446  return hitEnergyTrkFid_;
447 }
448 
449 
451 {
453  return hcalEnergySameTowers_;
454 }
455 
457 {
459  return nHcalHitsSameTowers_;
460 }
461 
463 {
465  return ecalEnergySameTowers_;
466 }
467 
469 {
471  return nEcalHitsSameTowers_;
472 }
473 
475 {
477  return trackEnergySameTowers_;
478 }
479 
481 {
483  return nTracksSameTowers_;
484 }
485 
486 void HBHEHitMap::hcalHitsSameTowers(std::set<const HBHERecHit*>& v) const
487 {
488  v.clear();
489  for(hitmap_const_iterator it1=beginHits(); it1!=endHits(); ++it1) {
490  for(std::set<const HBHERecHit*>::const_iterator it2=it1->second->hcalhits.begin(); it2!=it1->second->hcalhits.end(); ++it2) {
491  const HBHERecHit* hit=(*it2);
492  // if the hit in the tower is already in the hitmap, don't include it
493  if(findHit(hit)==endHits()) v.insert(hit);
494  }
495  }
496  return;
497 }
498 
499 void HBHEHitMap::ecalHitsSameTowers(std::set<const EcalRecHit*>& v) const
500 {
501  v.clear();
502  for(hitmap_const_iterator it1=beginHits(); it1!=endHits(); ++it1) {
503  v.insert(it1->second->ecalhits.begin(), it1->second->ecalhits.end());
504  }
505  return;
506 }
507 
508 void HBHEHitMap::tracksSameTowers(std::set<const reco::Track*>& v) const
509 {
510  v.clear();
511  for(hitmap_const_iterator it1=beginHits(); it1!=endHits(); ++it1) {
512  v.insert(it1->second->tracks.begin(), it1->second->tracks.end());
513  }
514  return;
515 }
516 
518 {
521 }
522 
524 {
527 }
528 
530 {
533 }
534 
536 {
539 }
540 
542 {
545 }
546 
548 {
550  return nTracksNeighborTowers_;
551 }
552 
553 void HBHEHitMap::hcalHitsNeighborTowers(std::set<const HBHERecHit*>& v) const
554 {
555  v.clear();
556  for(neighbor_const_iterator it1=beginNeighbors(); it1!=endNeighbors(); ++it1) {
557  const PhysicsTower* twr=(*it1);
558  v.insert(twr->hcalhits.begin(), twr->hcalhits.end());
559  }
560  return;
561 }
562 
563 void HBHEHitMap::ecalHitsNeighborTowers(std::set<const EcalRecHit*>& v) const
564 {
565  v.clear();
566  for(neighbor_const_iterator it1=beginNeighbors(); it1!=endNeighbors(); ++it1) {
567  const PhysicsTower* twr=(*it1);
568  v.insert(twr->ecalhits.begin(), twr->ecalhits.end());
569  }
570 
571  return;
572 }
573 
574 void HBHEHitMap::tracksNeighborTowers(std::set<const reco::Track*>& v) const
575 {
576  v.clear();
577  for(neighbor_const_iterator it1=beginNeighbors(); it1!=endNeighbors(); ++it1) {
578  const PhysicsTower* twr=(*it1);
579  v.insert(twr->tracks.begin(), twr->tracks.end());
580  }
581  return;
582 }
583 
584 void HBHEHitMap::byTowers(std::vector<twrinfo>& v) const
585 {
586  assert(0);
587 }
588 
589 void HBHEHitMap::insert(const HBHERecHit* hit, const PhysicsTower* twr, std::set<const PhysicsTower*>& neighbors)
590 {
591  hits_[hit]=twr;
592  neighbors_.insert(neighbors.begin(), neighbors.end());
593 
594  // make sure none of the neighbors are also are part of the hitmap
595  for(hitmap_const_iterator it=beginHits(); it!=endHits(); ++it) {
596  const PhysicsTower* t=it->second;
598 
599  // if a hit is also a neighbor, remove the neighbor
600  if(find!=endNeighbors()) neighbors_.erase(find);
601  }
602  return;
603 }
604 
605 
606 void HBHEHitMap::calcHits_(void) const
607 {
608  hitEnergy_=0;
609  nHits_=0;
611  for(hitmap_const_iterator it=hits_.begin(); it!=hits_.end(); ++it) {
612  const HBHERecHit* hit=it->first;
613  if(hit->id().ietaAbs()<=26) hitEnergyTrkFid_+=hit->energy();
614  hitEnergy_+=hit->energy();
615  ++nHits_;
616  }
617  return;
618 }
619 
621 {
624  std::set<const HBHERecHit*> v;
626  for(std::set<const HBHERecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
627  const HBHERecHit* hit=(*it);
630  }
631  return;
632 }
633 
635 {
638  std::set<const EcalRecHit*> v;
640  for(std::set<const EcalRecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
641  const EcalRecHit* hit=(*it);
644  }
645  return;
646 }
647 
649 {
652  std::set<const reco::Track*> v;
653  tracksSameTowers(v);
654  for(std::set<const reco::Track*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
655  const reco::Track* trk=(*it);
656  trackEnergySameTowers_+=trk->p();
658  }
659  return;
660 }
661 
663 {
666  std::set<const HBHERecHit*> v;
668  for(std::set<const HBHERecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
669  const HBHERecHit* hit=(*it);
672  }
673  return;
674 }
675 
677 {
680  std::set<const EcalRecHit*> v;
682  for(std::set<const EcalRecHit*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
683  const EcalRecHit* hit=(*it);
686  }
687  return;
688 }
689 
691 {
694  std::set<const reco::Track*> v;
696  for(std::set<const reco::Track*>::const_iterator it=v.begin(); it!=v.end(); ++it) {
697  const reco::Track* trk=(*it);
700  }
701  return;
702 }
703 
705 //
706 // HBHEHitMapOrganizer
707 //
709 
710 
712  const ObjectValidatorAbs& objvalidator,
713  const PhysicsTowerOrganizer& pto)
714 {
715  // loop over the hits
716  for(HBHERecHitCollection::const_iterator it=hbhehitcoll_h->begin(); it!=hbhehitcoll_h->end(); ++it) {
717  const HBHERecHit *hit=&(*it);
718  if(!objvalidator.validHit(*hit)) continue;
719 
720  // get the Physics Tower and the neighbors
721  const PhysicsTower* tower=pto.findTower(hit->id().ieta(), hit->id().iphi());
722 
723  std::set<const PhysicsTower*> neighbors;
724  pto.findNeighbors(hit->id().ieta(), hit->id().iphi(), neighbors);
725 
726  // organize the RBXs
727  int rbxidnum = HcalHPDRBXMap::indexRBX(hit->id());
728  rbxs_[rbxidnum].insert(hit, tower, neighbors);
729 
730  // organize the HPDs
731  int hpdidnum = HcalHPDRBXMap::indexHPD(hit->id());
732  hpds_[hpdidnum].insert(hit, tower, neighbors);
733 
734 
735  // organize the dihits
736  std::vector<const HBHERecHit*> hpdneighbors;
737  getHPDNeighbors(hit, hpdneighbors, pto);
738 
739  if(hpdneighbors.size()==1) {
740  std::vector<const HBHERecHit*> hpdneighborsneighbors;
741  getHPDNeighbors(hpdneighbors[0], hpdneighborsneighbors, pto);
742 
743  if(hpdneighborsneighbors.size()==1 && hpdneighborsneighbors[0]==hit && hit->energy()>hpdneighbors[0]->energy()) {
744  // we've found two hits who are neighbors in the same HPD, but who have no other
745  // neighbors (in the same HPD) in common. In order not to double-count, we
746  // require that the first hit has more energy
747 
748  const PhysicsTower* tower2=pto.findTower(hpdneighbors[0]->id().ieta(), hpdneighbors[0]->id().iphi());
749  std::set<const PhysicsTower*> neighbors2;
750  pto.findNeighbors(hpdneighbors[0]->id().ieta(), hpdneighbors[0]->id().iphi(), neighbors2);
751 
752  HBHEHitMap dihit;
753  dihit.insert(hit, tower, neighbors);
754  dihit.insert(hpdneighbors[0], tower2, neighbors2);
755  dihits_.push_back(dihit);
756  }
757  } else if(hpdneighbors.size()==0) {
758 
759  // organize the monohits
760  HBHEHitMap monohit;
761  monohit.insert(hit, tower, neighbors);
762  monohits_.push_back(monohit);
763  }
764 
765  } // finished looping over HBHERecHits
766  return;
767 }
768 
769 void HBHEHitMapOrganizer::getRBXs(std::vector<HBHEHitMap>& v, double energy) const
770 {
771  for(std::map<int, HBHEHitMap>::const_iterator it=rbxs_.begin(); it!=rbxs_.end(); ++it) {
772  const HBHEHitMap &map=it->second;
773  if(map.hitEnergy()>energy) v.push_back(map);
774  }
775  return;
776 }
777 
778 void HBHEHitMapOrganizer::getHPDs(std::vector<HBHEHitMap>& v, double energy) const
779 {
780  for(std::map<int, HBHEHitMap>::const_iterator it=hpds_.begin(); it!=hpds_.end(); ++it) {
781  const HBHEHitMap &map=it->second;
782  if(map.hitEnergy()>energy) v.push_back(map);
783  }
784  return;
785 }
786 
787 void HBHEHitMapOrganizer::getDiHits(std::vector<HBHEHitMap>& v, double energy) const
788 {
789  for(std::vector<HBHEHitMap>::const_iterator it=dihits_.begin(); it!=dihits_.end(); ++it) {
790  if(it->hitEnergy()>energy) v.push_back(*it);
791  }
792  return;
793 }
794 
795 void HBHEHitMapOrganizer::getMonoHits(std::vector<HBHEHitMap>& v, double energy) const
796 {
797  for(std::vector<HBHEHitMap>::const_iterator it=monohits_.begin(); it!=monohits_.end(); ++it) {
798  if(it->hitEnergy()>energy) v.push_back(*it);
799  }
800  return;
801 }
802 
803 
804 
805 void HBHEHitMapOrganizer::getHPDNeighbors(const HBHERecHit* hit, std::vector<const HBHERecHit*>& neighbors, const PhysicsTowerOrganizer& pto)
806 {
807  std::set<const PhysicsTower*> temp;
808  pto.findNeighbors(hit->id().ieta(), hit->id().iphi(), temp);
809 
810  // make sure to include the same tower that the hit is in
811  temp.insert(pto.findTower(hit->id().ieta(), hit->id().iphi()));
812 
813  // loop over the rechits in the temp neighbors
814  for(std::set<const PhysicsTower*>::const_iterator it1=temp.begin(); it1!=temp.end(); ++it1) {
815  for(std::set<const HBHERecHit*>::const_iterator it2=(*it1)->hcalhits.begin(); it2!=(*it1)->hcalhits.end(); ++it2) {
816  const HBHERecHit* hit2(*it2);
817  if(hit!=hit2 && HcalHPDRBXMap::indexHPD(hit->id())==HcalHPDRBXMap::indexHPD(hit2->id())) {
818  neighbors.push_back(hit2);
819  }
820  }
821  }
822  return;
823 }
824 
825 
826 
void getHPDNeighbors(const HBHERecHit *hit, std::vector< const HBHERecHit * > &neighbors, const PhysicsTowerOrganizer &pto)
double p() const
momentum vector magnitude
Definition: TrackBase.h:568
double hcalEnergySameTowers(void) const
std::set< const HBHERecHit * > hcalhits
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
void tracksNeighborTowers(std::set< const reco::Track * > &v) const
int nHcalHitsSameTowers(void) const
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:30
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id) const
Evaluate status from id use channelStatus from DB.
std::map< const HBHERecHit *, const PhysicsTower * >::const_iterator hitmap_const_iterator
int ietaAbs() const
get the absolute value of the tower ieta
std::vector< HBHEHitMap > dihits_
hitmap_const_iterator findHit(const HBHERecHit *hit) const
const DetId & detid() const
Definition: CaloRecHit.h:20
static int indexRBX(const HcalDetId &)
void getHPDs(std::vector< HBHEHitMap > &v, double energy) const
HcalDetId id() const
get the id
Definition: HBHERecHit.h:23
assert(m_qm.get())
neighbor_const_iterator findNeighbor(const PhysicsTower *twr) const
const DetId & detid() const
Definition: EcalRecHit.h:71
std::vector< HBHERecHit >::const_iterator const_iterator
virtual bool validTrack(const reco::Track &) const =0
int nEcalHitsNeighborTowers(void) const
const Item * getValues(DetId fId, bool throwOnFail=true) const
bool validTrack(const reco::Track &) const
neighbor_const_iterator beginNeighbors(void) const
int zside(DetId const &)
double hcalEnergyNeighborTowers_
void getRBXs(std::vector< HBHEHitMap > &v, double energy) const
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:628
ObjectValidator(const edm::ParameterSet &)
void ecalHitsNeighborTowers(std::set< const EcalRecHit * > &v) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
hitmap_const_iterator endHits(void) const
std::map< const HBHERecHit *, const PhysicsTower * > hits_
void getMonoHits(std::vector< HBHEHitMap > &v, double energy) const
void ecalHitsSameTowers(std::set< const EcalRecHit * > &v) const
double trackEnergyNeighborTowers_
PhysicsTowerOrganizer(const edm::Event &iEvent, const edm::EventSetup &evSetup, const edm::Handle< HBHERecHitCollection > &hbhehitcoll_h, const edm::Handle< EcalRecHitCollection > &ebhitcoll_h, const edm::Handle< EcalRecHitCollection > &eehitcoll_h, const edm::Handle< std::vector< reco::TrackExtrapolation > > &trackextrapcoll_h, const ObjectValidatorAbs &objectvalidator, const CaloTowerConstituentsMap &ctcm)
void calcTracksSameTowers_(void) const
double ecalEnergyNeighborTowers(void) const
std::set< const PhysicsTower * > neighbors_
void insert(const HBHERecHit *hit, const PhysicsTower *twr, std::set< const PhysicsTower * > &neighbors)
const PhysicsTower * findTower(const CaloTowerDetId &id) const
void tracksSameTowers(std::set< const reco::Track * > &v) const
int iEvent
Definition: GenABIO.cc:230
int nHcalHitsNeighborTowers(void) const
CaloTowerDetId towerOf(const DetId &id) const
Get the tower id for this det id (or null if not known)
void calcHcalSameTowers_(void) const
std::map< int, HBHEHitMap > rbxs_
float energy() const
Definition: CaloRecHit.h:17
double pt() const
track transverse momentum
Definition: TrackBase.h:574
const EcalChannelStatus * theEcalChStatus_
bool recoveredRecHit(const DetId &myid, const uint32_t &myflag) const
const EcalRecHitCollection * theEERecHitCollection_
uint32_t flags() const
Definition: CaloRecHit.h:21
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
double hcalEnergyNeighborTowers(void) const
int iphi() const
get the tower iphi
double trackEnergySameTowers(void) const
void insert_(CaloTowerDetId &id, const HBHERecHit *hit)
float energy() const
Definition: EcalRecHit.h:68
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:773
const EcalRecHitCollection * theEBRecHitCollection_
std::set< const PhysicsTower * >::const_iterator neighbor_const_iterator
std::map< int, HBHEHitMap > hpds_
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:34
virtual DetId getClosestCell(const GlobalPoint &r) const
std::set< PhysicsTower, towercmp > towers_
int iphi() const
get the cell iphi
Definition: HcalDetId.h:38
Definition: DetId.h:18
double hitEnergy(void) const
void calcEcalSameTowers_(void) const
DetId id() const
get the id
Definition: EcalRecHit.h:76
virtual bool validHit(const HBHERecHit &) const =0
double ecalEnergyNeighborTowers_
void calcTracksNeighborTowers_(void) const
void calcHcalNeighborTowers_(void) const
const T & get() const
Definition: EventSetup.h:55
HBHEHitMapOrganizer(const edm::Handle< HBHERecHitCollection > &hbhehitcoll_h, const ObjectValidatorAbs &objvalidator, const PhysicsTowerOrganizer &pto)
T const * product() const
Definition: ESHandle.h:86
int zside() const
get the z-side of the tower (1/-1)
const EcalSeverityLevelAlgo * theEcalSevLvlAlgo_
int getSeverityLevel(const DetId &myid, const uint32_t &myflag, const uint32_t &mystatus) const
std::set< const EcalRecHit * > ecalhits
reco::TrackRef const & track() const
int nHits(void) const
double trackEnergyNeighborTowers(void) const
double hitEnergyTrackFiducial(void) const
void calcHits_(void) const
std::set< const reco::Track * > tracks
void calcEcalNeighborTowers_(void) const
std::vector< HBHEHitMap > monohits_
int nTracksSameTowers(void) const
void getDiHits(std::vector< HBHEHitMap > &v, double energy) const
std::vector< Point > const & positions() const
double ecalEnergySameTowers(void) const
const HcalChannelQuality * theHcalChStatus_
int nEcalHitsSameTowers(void) const
static int indexHPD(const HcalDetId &)
int nTracksNeighborTowers(void) const
uint32_t getValue() const
void hcalHitsSameTowers(std::set< const HBHERecHit * > &v) const
bool validHit(const HBHERecHit &) const
const HcalSeverityLevelComputer * theHcalSevLvlComputer_
void byTowers(std::vector< twrinfo > &v) const
neighbor_const_iterator endNeighbors(void) const
void findNeighbors(const CaloTowerDetId &id, std::set< const PhysicsTower * > &neighbors) const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
hitmap_const_iterator beginHits(void) const
void hcalHitsNeighborTowers(std::set< const HBHERecHit * > &v) const