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  dEdXWeights_(iConfig.getParameter<std::vector<double> >("dEdXWeights"))
10 {
11  isoHelper_.setDeltaR(iConfig.getParameter<double>("isoDeltaR"));
12  isoHelper_.setNRings(iConfig.getParameter<unsigned int>("isoNRings"));
13  isoHelper_.setMinDeltaR(iConfig.getParameter<double>("isoDeltaRmin"));
14 
19  debug_ = iConfig.getUntrackedParameter<bool>("debug", false);
20 }
21 
23  edm::Handle<HGCRecHitCollection> recHitHandleEE;
24  iEvent.getByToken(recHitsEE_, recHitHandleEE);
25  edm::Handle<HGCRecHitCollection> recHitHandleFH;
26  iEvent.getByToken(recHitsFH_, recHitHandleFH);
27  edm::Handle<HGCRecHitCollection> recHitHandleBH;
28  iEvent.getByToken(recHitsBH_, recHitHandleBH);
29 
33  pcaHelper_.fillHitMap(*recHitHandleEE,*recHitHandleFH,*recHitHandleBH);
34  isoHelper_.setRecHits(recHitHandleEE, recHitHandleFH, recHitHandleBH);
35 }
36 
37 void HGCalEgammaIDHelper::computeHGCAL(const reco::Photon & thePhoton, float radius) {
38  if (thePhoton.isEB()) {
39  if (debug_) 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
53  pcaHelper_.computePCA(radius);
54  // second computation within cylinder, halo hits included
55  pcaHelper_.computePCA(radius);
57 
58  // isolation
59  isoHelper_.produceHGCalIso(thePhoton.superCluster()->seed());
60 }
61 
63  if (theElectron.isEB()) {
64  if (debug_) std::cout << "The electron is in the barrel" <<std::endl;
65  pcaHelper_.clear();
66  return;
67  }
68 
69  pcaHelper_.storeRecHits(*theElectron.electronCluster());
70  if (debug_)
71  std::cout << " Stored the hits belonging to the electronCluster " << std::endl;
72 
73  // initial computation, no radius cut, but halo hits not taken
74  if (debug_)
75  std::cout << " Calling PCA initial computation" << std::endl;
77  // first computation within cylinder, halo hits included
78  pcaHelper_.computePCA(radius);
79  // second computation within cylinder, halo hits included
80  pcaHelper_.computePCA(radius);
83 }
void computeShowerWidth(float radius, bool withHalo=true)
CaloClusterPtr electronCluster() const
Definition: GsfElectron.h:248
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void setMinDeltaR(const float dr)
void computePCA(float radius, bool withHalo=true)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
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
hgcal::RecHitTools recHitTools_
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
void setDeltaR(const float dr)
bool isEB() const
Definition: GsfElectron.h:356
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:73
void storeRecHits(const reco::CaloCluster &theCluster)
int iEvent
Definition: GenABIO.cc:224
void setRecHitTools(const hgcal::RecHitTools *recHitTools)
void produceHGCalIso(const reco::CaloClusterPtr &seedCluster)
edm::InputTag bhRecHitInputTag_
void setRecHitTools(const hgcal::RecHitTools *recHitTools)
std::vector< double > dEdXWeights_
hgcal::EGammaPCAHelper pcaHelper_
void computeHGCAL(const reco::Photon &thePhoton, float radius)
bool isEB() const
Definition: Photon.h:121
void eventInit(const edm::Event &iEvent, const edm::EventSetup &iSetup)
HLT enums.
HGCalIsoCalculator isoHelper_
void fillHitMap(const HGCRecHitCollection &HGCEERecHits, const HGCRecHitCollection &HGCFHRecHits, const HGCRecHitCollection &HGCBHRecHits)
to compute from inside - once per event
edm::InputTag eeRecHitInputTag_