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 
11 HGCalIsoCalculator::HGCalIsoCalculator():dr2_(0.15*0.15),mindr2_(0),rechittools_(nullptr),debug_(false),nlayers_(30){
12  setNRings(5);
13 }
14 
16 }
17 
22  ){
23  recHitsEE_=hitsEE;
24  recHitsFH_=hitsFH;
25  recHitsBH_=hitsBH;
26 
27  if(!rechittools_)
28  throw cms::Exception("HGCalIsoCalculator::produceHGCalIso: rechittools not set");
29 
30  hitEtaPhiCache_.clear();
31  hitEtaPhiCache_.reserve(recHitsEE_->size()+recHitsFH_->size()+recHitsBH_->size());
32 
33  // Since HGCal is not projective and the rechits don't cache any
34  // eta,phi, we make our own here
35  auto makeEtaPhiPair = [this](const auto& hit) {
37  float eta=rechittools_->getEta(position, 0); //assume vertex at z=0
38  float phi=rechittools_->getPhi(position);
39  return std::make_pair(eta, phi);
40  };
41 
42  for(const auto& hit : *recHitsEE_) hitEtaPhiCache_.push_back(makeEtaPhiPair(hit));
43  for(const auto& hit : *recHitsFH_) hitEtaPhiCache_.push_back(makeEtaPhiPair(hit));
44  for(const auto& hit : *recHitsBH_) hitEtaPhiCache_.push_back(makeEtaPhiPair(hit));
45 }
46 
48 
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, std::pair<float,float> etaphiVal) {
62  float deltar2=reco::deltaR2(etaphiVal.first,etaphiVal.second,seedEta,seedPhi);
63  if(deltar2>dr2_ || deltar2<mindr2_) return;
64 
65  unsigned int layer=rechittools_->getLayerWithOffset(hit.id());
66  if(layer>=nlayers_) return;
67 
68  //do not consider hits associated to the photon cluster
69  if( std::none_of(seedHitsAndFractions.begin(), seedHitsAndFractions.end(), [&hit](const auto& seedhit){return hit.id()==seedhit.first;}) ){
70  const unsigned int ring=ringasso_.at(layer);
71  isoringdeposits_.at(ring)+=hit.energy();
72  }
73  };
74 
75  // The cache order is EE,FH,BH, so we should loop over them the same here
76  auto itEtaPhiCache = hitEtaPhiCache_.cbegin();
77  for(const auto& hit : *recHitsEE_){
78  checkAndFill(hit, *itEtaPhiCache);
79  itEtaPhiCache++;
80  }
81  for(const auto& hit : *recHitsFH_){
82  checkAndFill(hit, *itEtaPhiCache);
83  itEtaPhiCache++;
84  }
85  for(const auto& hit : *recHitsBH_){
86  checkAndFill(hit, *itEtaPhiCache);
87  itEtaPhiCache++;
88  }
89 }
90 
91 void HGCalIsoCalculator::setNRings(const size_t nrings){
92  if(nrings>nlayers_)
93  throw std::logic_error("PhotonHGCalIsoCalculator::setNRings: max number of rings reached");
94 
95  ringasso_.clear();
96  isoringdeposits_.clear();
97  unsigned int separator=nlayers_/nrings;
98  size_t counter=0;
99  for(size_t i=0;i<nlayers_+1;i++){
100  ringasso_.push_back(counter);
101  //the last ring might be larger.
102  if(i && !(i%separator) && (int)counter<(int)nrings-1){
103  counter++;
104  }
105  }
106  isoringdeposits_.resize(nrings,0);
107 }
108 
109 const float HGCalIsoCalculator::getIso(const unsigned int ring)const{
110  if(ring>=isoringdeposits_.size())
111  throw cms::Exception("HGCalIsoCalculator::getIso: ring index out of range");
112  return isoringdeposits_[ring];
113 }
string separator
Definition: mps_merge.py:79
edm::Handle< HGCRecHitCollection > recHitsEE_
edm::Handle< HGCRecHitCollection > recHitsBH_
#define nullptr
void setNRings(const size_t nrings)
void setRecHits(edm::Handle< HGCRecHitCollection > hitsEE, edm::Handle< HGCRecHitCollection > hitsFH, edm::Handle< HGCRecHitCollection > hitsBH)
fill - once per event
std::vector< unsigned int > ringasso_
const hgcal::RecHitTools * rechittools_
float getEta(const GlobalPoint &position, const float &vertex_z=0.) const
Definition: RecHitTools.cc:340
float getPhi(const GlobalPoint &position) const
Definition: RecHitTools.cc:351
std::vector< std::pair< float, float > > hitEtaPhiCache_
void produceHGCalIso(const reco::CaloClusterPtr &seedCluster)
edm::Handle< HGCRecHitCollection > recHitsFH_
std::vector< float > isoringdeposits_
unsigned int id
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:283
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:112
static int position[264][3]
Definition: ReadPGInfo.cc:509
const float getIso(const unsigned int ring) const