CMS 3D CMS Logo

HGCalEgammaIDHelper.cc
Go to the documentation of this file.
2 
3 #include <iostream>
4 
6  : eeRecHitInputTag_(iConfig.getParameter<edm::InputTag>("EERecHits")),
7  fhRecHitInputTag_(iConfig.getParameter<edm::InputTag>("FHRecHits")),
8  bhRecHitInputTag_(iConfig.getParameter<edm::InputTag>("BHRecHits")),
9  hitMapInputTag_(iConfig.getParameter<edm::InputTag>("hitMapTag")),
10  dEdXWeights_(iConfig.getParameter<std::vector<double>>("dEdXWeights")) {
11  isoHelper_.setDeltaR(iConfig.getParameter<double>("isoDeltaR"));
12  isoHelper_.setNRings(iConfig.getParameter<unsigned int>("isoNRings"));
13  isoHelper_.setMinDeltaR(iConfig.getParameter<double>("isoDeltaRmin"));
14 
18  hitMap_ = iC.consumes<std::unordered_map<DetId, const HGCRecHit*>>(hitMapInputTag_);
19  caloGeometry_ = iC.esConsumes();
21  debug_ = iConfig.getUntrackedParameter<bool>("debug", false);
22 }
23 
25  auto recHitHandleEE = iEvent.getHandle(recHitsEE_);
26  auto recHitHandleFH = iEvent.getHandle(recHitsFH_);
27  auto recHitHandleBH = iEvent.getHandle(recHitsBH_);
28 
33  isoHelper_.setRecHits(recHitHandleEE, recHitHandleFH, recHitHandleBH);
34 }
35 
36 void HGCalEgammaIDHelper::computeHGCAL(const reco::Photon& thePhoton, float radius) {
37  if (thePhoton.isEB()) {
38  if (debug_)
39  std::cout << "The photon is in the barrel" << std::endl;
40  pcaHelper_.clear();
41  return;
42  }
43 
44  pcaHelper_.storeRecHits(*thePhoton.superCluster()->seed());
45  if (debug_)
46  std::cout << " Stored the hits belonging to the photon superCluster seed " << std::endl;
47 
48  // initial computation, no radius cut, but halo hits not taken
49  if (debug_)
50  std::cout << " Calling PCA initial computation" << std::endl;
52  // first computation within cylinder, halo hits included
54  // second computation within cylinder, halo hits included
57 
58  // isolation
59  isoHelper_.produceHGCalIso(thePhoton.superCluster()->seed());
60 }
61 
63  if (theElectron.isEB()) {
64  if (debug_)
65  std::cout << "The electron is in the barrel" << std::endl;
66  pcaHelper_.clear();
67  return;
68  }
69 
70  pcaHelper_.storeRecHits(*theElectron.electronCluster());
71  if (debug_)
72  std::cout << " Stored the hits belonging to the electronCluster " << std::endl;
73 
74  // initial computation, no radius cut, but halo hits not taken
75  if (debug_)
76  std::cout << " Calling PCA initial computation" << std::endl;
78  // first computation within cylinder, halo hits included
80  // second computation within cylinder, halo hits included
84 }
void computeShowerWidth(float radius, bool withHalo=true)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void setMinDeltaR(const float dr)
void computePCA(float radius, bool withHalo=true)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
edm::InputTag fhRecHitInputTag_
void setNRings(const size_t nrings)
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
void setdEdXWeights(const std::vector< double > &dEdX)
void setRecHits(edm::Handle< HGCRecHitCollection > hitsEE, edm::Handle< HGCRecHitCollection > hitsFH, edm::Handle< HGCRecHitCollection > hitsBH)
fill - once per event
bool isEB() const
Definition: Photon.h:126
hgcal::RecHitTools recHitTools_
void setDeltaR(const float dr)
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
bool isEB() const
Definition: GsfElectron.h:328
T getUntrackedParameter(std::string const &, T const &) const
void storeRecHits(const reco::CaloCluster &theCluster)
int iEvent
Definition: GenABIO.cc:224
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
void setRecHitTools(const hgcal::RecHitTools *recHitTools)
void produceHGCalIso(const reco::CaloClusterPtr &seedCluster)
edm::InputTag bhRecHitInputTag_
void setRecHitTools(const hgcal::RecHitTools *recHitTools)
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometry_
edm::InputTag hitMapInputTag_
std::vector< double > dEdXWeights_
hgcal::EGammaPCAHelper pcaHelper_
void computeHGCAL(const reco::Photon &thePhoton, float radius)
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79
void eventInit(const edm::Event &iEvent, const edm::EventSetup &iSetup)
HLT enums.
CaloClusterPtr electronCluster() const
Definition: GsfElectron.h:220
HGCalIsoCalculator isoHelper_
void setHitMap(const std::unordered_map< DetId, const HGCRecHit *> *hitMap)
to set once per event
edm::EDGetTokenT< std::unordered_map< DetId, const HGCRecHit * > > hitMap_
edm::InputTag eeRecHitInputTag_