CMS 3D CMS Logo

ElectronHcalHelper.cc
Go to the documentation of this file.
5 
6 using namespace reco;
7 
9  if (cfg_.hOverEConeSize == 0. and !cfg_.onlyBehindCluster) {
10  return;
11  }
12 
13  caloGeometryToken_ = cc.esConsumes();
14  hcalTopologyToken_ = cc.esConsumes();
15  hcalChannelQualityToken_ = cc.esConsumes(edm::ESInputTag("", "withTopo"));
16  hcalSevLvlComputerToken_ = cc.esConsumes();
17  towerMapToken_ = cc.esConsumes();
18 }
19 
21  if (cfg_.hOverEConeSize == 0. and !cfg_.onlyBehindCluster) {
22  return;
23  }
24 
29 
30  if (cfg_.onlyBehindCluster) {
31  hcalIso_ = std::make_unique<EgammaHcalIsolation>(EgammaHcalIsolation::InclusionRule::isBehindClusterSeed,
32  0.,
34  0.,
35  cfg_.eThresHB,
36  EgammaHcalIsolation::arrayHB{{0., 0., 0., 0.}},
38  cfg_.eThresHE,
39  EgammaHcalIsolation::arrayHE{{0., 0., 0., 0., 0., 0., 0.}},
41  evt.get(cfg_.hbheRecHits),
42  eventSetup.getHandle(caloGeometryToken_),
43  eventSetup.getHandle(hcalTopologyToken_),
46  towerMap_);
47  } else {
48  hcalIso_ = std::make_unique<EgammaHcalIsolation>(EgammaHcalIsolation::InclusionRule::withinConeAroundCluster,
51  0.,
52  cfg_.eThresHB,
53  EgammaHcalIsolation::arrayHB{{0., 0., 0., 0.}},
55  cfg_.eThresHE,
56  EgammaHcalIsolation::arrayHE{{0., 0., 0., 0., 0., 0., 0.}},
58  evt.get(cfg_.hbheRecHits),
59  eventSetup.getHandle(caloGeometryToken_),
60  eventSetup.getHandle(hcalTopologyToken_),
63  towerMap_);
64  }
65 }
66 
68  return (cfg_.checkHcalStatus)
70  : true;
71 }
72 
73 double ElectronHcalHelper::hcalESum(const SuperCluster& sc, int depth, const HcalPFCuts* hcalCuts) const {
74  return (cfg_.onlyBehindCluster) ? hcalIso_->getHcalESumBc(&sc, depth, hcalCuts)
75  : (cfg_.hOverEConeSize > 0.) ? hcalIso_->getHcalESum(&sc, depth, hcalCuts)
76  : 0.;
77 }
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > hcalTopologyToken_
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:344
bool hasActiveHcal(const reco::SuperCluster &sc) const
edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd > hcalSevLvlComputerToken_
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
HcalChannelQuality const * hcalChannelQuality_
const Configuration cfg_
ElectronHcalHelper(const Configuration &cfg, edm::ConsumesCollector &&cc)
edm::EDGetTokenT< HBHERecHitCollection > hbheRecHits
std::unique_ptr< EgammaHcalIsolation > hcalIso_
HcalSeverityLevelComputer const * hcalSevLvlComputer_
void beginEvent(const edm::Event &evt, const edm::EventSetup &eventSetup)
bool hasActiveHcal(std::vector< CaloTowerDetId > const &towers, CaloTowerConstituentsMap const &towerMap, HcalChannelQuality const &hcalQuality, HcalTopology const &hcalTopology)
edm::ESGetToken< CaloTowerConstituentsMap, CaloGeometryRecord > towerMapToken_
edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcd > hcalChannelQualityToken_
auto hcalTowersBehindClusters(const reco::SuperCluster &sc) const
HcalTopology const * hcalTopology_
fixed size matrix
EgammaHcalIsolation::arrayHB eThresHB
EgammaHcalIsolation::arrayHE eThresHE
double hcalESum(const reco::SuperCluster &, int depth, const HcalPFCuts *hcalCuts) const
std::array< double, 4 > arrayHB
CaloTowerConstituentsMap const * towerMap_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
std::array< double, 7 > arrayHE