CMS 3D CMS Logo

HGCalIsoCalculator.cc
Go to the documentation of this file.
1 /*
2  * HGCalIsoCalculator.cc
3  *
4  * Created on: 13 Oct 2017
5  * Author: jkiesele, ncsmith
6  */
7 
10 
12  : dr2_(0.15 * 0.15), mindr2_(0), rechittools_(nullptr), debug_(false), nlayers_(30) {
13  setNRings(5);
14 }
15 
17 
21  recHitsEE_ = hitsEE;
22  recHitsFH_ = hitsFH;
23  recHitsBH_ = hitsBH;
24 
25  if (!rechittools_)
26  throw cms::Exception("HGCalIsoCalculator::produceHGCalIso: rechittools not set");
27 
28  hitEtaPhiCache_.clear();
29  hitEtaPhiCache_.reserve(recHitsEE_->size() + recHitsFH_->size() + recHitsBH_->size());
30 
31  // Since HGCal is not projective and the rechits don't cache any
32  // eta,phi, we make our own here
33  auto makeEtaPhiPair = [this](const auto& hit) {
35  float eta = rechittools_->getEta(position, 0); //assume vertex at z=0
36  float phi = rechittools_->getPhi(position);
37  return std::make_pair(eta, phi);
38  };
39 
40  for (const auto& hit : *recHitsEE_)
41  hitEtaPhiCache_.push_back(makeEtaPhiPair(hit));
42  for (const auto& hit : *recHitsFH_)
43  hitEtaPhiCache_.push_back(makeEtaPhiPair(hit));
44  for (const auto& hit : *recHitsBH_)
45  hitEtaPhiCache_.push_back(makeEtaPhiPair(hit));
46 }
47 
49  if (!rechittools_)
50  throw cms::Exception("HGCalIsoCalculator::produceHGCalIso: rechittools not set");
51 
52  for (auto& r : isoringdeposits_)
53  r = 0;
54 
55  // make local temporaries to pass to the lambda
56  // avoids recomputing every iteration
57  const float seedEta = seed->eta();
58  const float seedPhi = seed->phi();
59  const std::vector<std::pair<DetId, float>>& seedHitsAndFractions = seed->hitsAndFractions();
60 
61  auto checkAndFill = [this, &seedEta, &seedPhi, &seedHitsAndFractions](const HGCRecHit& hit,
62  std::pair<float, float> etaphiVal) {
63  float deltar2 = reco::deltaR2(etaphiVal.first, etaphiVal.second, seedEta, seedPhi);
64  if (deltar2 > dr2_ || deltar2 < mindr2_)
65  return;
66 
67  unsigned int layer = rechittools_->getLayerWithOffset(hit.id());
68  if (layer >= nlayers_)
69  return;
70 
71  //do not consider hits associated to the photon cluster
72  if (std::none_of(seedHitsAndFractions.begin(), seedHitsAndFractions.end(), [&hit](const auto& seedhit) {
73  return hit.id() == seedhit.first;
74  })) {
75  const unsigned int ring = ringasso_.at(layer);
76  isoringdeposits_.at(ring) += hit.energy();
77  }
78  };
79 
80  // The cache order is EE,FH,BH, so we should loop over them the same here
81  auto itEtaPhiCache = hitEtaPhiCache_.cbegin();
82  for (const auto& hit : *recHitsEE_) {
83  checkAndFill(hit, *itEtaPhiCache);
84  itEtaPhiCache++;
85  }
86  for (const auto& hit : *recHitsFH_) {
87  checkAndFill(hit, *itEtaPhiCache);
88  itEtaPhiCache++;
89  }
90  for (const auto& hit : *recHitsBH_) {
91  checkAndFill(hit, *itEtaPhiCache);
92  itEtaPhiCache++;
93  }
94 }
95 
96 void HGCalIsoCalculator::setNRings(const size_t nrings) {
97  if (nrings > nlayers_)
98  throw std::logic_error("PhotonHGCalIsoCalculator::setNRings: max number of rings reached");
99 
100  ringasso_.clear();
101  isoringdeposits_.clear();
102  unsigned int separator = nlayers_ / nrings;
103  size_t counter = 0;
104  for (size_t i = 0; i < nlayers_ + 1; i++) {
105  ringasso_.push_back(counter);
106  //the last ring might be larger.
107  if (i && !(i % separator) && (int)counter < (int)nrings - 1) {
108  counter++;
109  }
110  }
111  isoringdeposits_.resize(nrings, 0);
112 }
113 
114 const float HGCalIsoCalculator::getIso(const unsigned int ring) const {
115  if (ring >= isoringdeposits_.size())
116  throw cms::Exception("HGCalIsoCalculator::getIso: ring index out of range");
117  return isoringdeposits_[ring];
118 }
counter
Definition: counter.py:1
mps_fire.i
i
Definition: mps_fire.py:355
funct::false
false
Definition: Factorize.h:34
hit::id
unsigned int id
Definition: SiStripHitEffFromCalibTree.cc:92
HGCalIsoCalculator::ringasso_
std::vector< unsigned int > ringasso_
Definition: HGCalIsoCalculator.h:73
HGCalIsoCalculator::setRecHits
void setRecHits(edm::Handle< HGCRecHitCollection > hitsEE, edm::Handle< HGCRecHitCollection > hitsFH, edm::Handle< HGCRecHitCollection > hitsBH)
fill - once per event
Definition: HGCalIsoCalculator.cc:18
HGCalIsoCalculator::getIso
const float getIso(const unsigned int ring) const
Definition: HGCalIsoCalculator.cc:114
HGCalIsoCalculator::HGCalIsoCalculator
HGCalIsoCalculator()
Definition: HGCalIsoCalculator.cc:11
edm::Handle
Definition: AssociativeIterator.h:50
mps_merge.separator
string separator
Definition: mps_merge.py:79
HGCalIsoCalculator::rechittools_
const hgcal::RecHitTools * rechittools_
Definition: HGCalIsoCalculator.h:77
deltaR.h
HGCalIsoCalculator::recHitsFH_
edm::Handle< HGCRecHitCollection > recHitsFH_
Definition: HGCalIsoCalculator.h:79
PVValHelper::eta
Definition: PVValidationHelpers.h:69
Point3DBase< float, GlobalTag >
HGCalIsoCalculator::mindr2_
float mindr2_
Definition: HGCalIsoCalculator.h:75
HGCalIsoCalculator::setNRings
void setNRings(const size_t nrings)
Definition: HGCalIsoCalculator.cc:96
HGCRecHit
Definition: HGCRecHit.h:14
HGCalIsoCalculator::isoringdeposits_
std::vector< float > isoringdeposits_
Definition: HGCalIsoCalculator.h:72
hgcal::RecHitTools::getPhi
float getPhi(const GlobalPoint &position) const
Definition: RecHitTools.cc:437
hgcal::RecHitTools::getLayerWithOffset
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:352
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
HGCalIsoCalculator::hitEtaPhiCache_
std::vector< std::pair< float, float > > hitEtaPhiCache_
Definition: HGCalIsoCalculator.h:81
HGCalIsoCalculator::nlayers_
unsigned int nlayers_
Definition: HGCalIsoCalculator.h:83
HGCalIsoCalculator::dr2_
float dr2_
Definition: HGCalIsoCalculator.h:75
edm::Ptr< CaloCluster >
alignCSCRings.r
r
Definition: alignCSCRings.py:93
DDAxes::phi
hgcal::RecHitTools::getEta
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:426
HGCalIsoCalculator::produceHGCalIso
void produceHGCalIso(const reco::CaloClusterPtr &seedCluster)
Definition: HGCalIsoCalculator.cc:48
Exception
Definition: hltDiff.cc:246
relativeConstraints.ring
ring
Definition: relativeConstraints.py:68
HGCalIsoCalculator.h
HGCalIsoCalculator::recHitsEE_
edm::Handle< HGCRecHitCollection > recHitsEE_
Definition: HGCalIsoCalculator.h:78
cms::Exception
Definition: Exception.h:70
HGCalIsoCalculator::~HGCalIsoCalculator
~HGCalIsoCalculator()
Definition: HGCalIsoCalculator.cc:16
hgcal::RecHitTools::getPosition
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:126
HGCalIsoCalculator::recHitsBH_
edm::Handle< HGCRecHitCollection > recHitsBH_
Definition: HGCalIsoCalculator.h:80
SurveyInfoScenario_cff.seed
seed
Definition: SurveyInfoScenario_cff.py:295
hit
Definition: SiStripHitEffFromCalibTree.cc:88