CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EgammaRecHitIsolation Class Reference

#include <EgammaRecHitIsolation.h>

Public Member Functions

void doFlagChecks (const std::vector< int > &v)
 
void doSeverityChecks (const EcalRecHitCollection *const recHits, const std::vector< int > &v)
 
 EgammaRecHitIsolation (double extRadius, double intRadius, double etaSlice, double etLow, double eLow, edm::ESHandle< CaloGeometry >, const EcalRecHitCollection &, const EcalSeverityLevelAlgo *, DetId::Detector detector)
 
double getEnergySum (const reco::Candidate *emObject, EcalPFRecHitThresholds const &thresholds) const
 
double getEnergySum (const reco::SuperCluster *emObject, EcalPFRecHitThresholds const &thresholds) const
 
double getEtSum (const reco::Candidate *emObject, EcalPFRecHitThresholds const &thresholds) const
 
double getEtSum (const reco::SuperCluster *emObject, EcalPFRecHitThresholds const &thresholds) const
 
void setUseNumCrystals (bool b=true)
 
void setVetoClustered (bool b=true)
 
 ~EgammaRecHitIsolation ()
 

Private Member Functions

double getSum_ (const reco::Candidate *, bool returnEt, const EcalPFRecHitThresholds *thresholds) const
 
double getSum_ (const reco::SuperCluster *, bool returnEt, const EcalPFRecHitThresholds *thresholds) const
 

Private Attributes

const EcalRecHitCollectioncaloHits_
 
const EcalRecHitCollectionecalBarHits_
 
double eLow_
 
double etaSlice_
 
double etLow_
 
double extRadius_
 
std::vector< int > flags_
 
double intRadius_
 
std::vector< int > severitiesexcl_
 
const EcalSeverityLevelAlgosevLevel_
 
const CaloSubdetectorGeometrysubdet_ [2]
 
edm::ESHandle< CaloGeometrytheCaloGeom_
 
bool useNumCrystals_
 
bool vetoClustered_
 

Detailed Description

Definition at line 28 of file EgammaRecHitIsolation.h.

Constructor & Destructor Documentation

◆ EgammaRecHitIsolation()

EgammaRecHitIsolation::EgammaRecHitIsolation ( double  extRadius,
double  intRadius,
double  etaSlice,
double  etLow,
double  eLow,
edm::ESHandle< CaloGeometry theCaloGeom,
const EcalRecHitCollection caloHits,
const EcalSeverityLevelAlgo sl,
DetId::Detector  detector 
)

Definition at line 32 of file EgammaRecHitIsolation.cc.

References DetId::Ecal, EcalBarrel, EcalEndcap, CaloGeometry::getSubdetectorGeometry(), edm::ESHandle< T >::product(), subdet_, and theCaloGeom_.

41  : // not used anymore, kept for compatibility
44  etaSlice_(etaSlice),
45  etLow_(etLow),
46  eLow_(eLow),
47  theCaloGeom_(theCaloGeom),
48  caloHits_(caloHits),
49  sevLevel_(sl),
50  useNumCrystals_(false),
51  vetoClustered_(false),
52  ecalBarHits_(nullptr),
53  //chStatus_(0),
54  severitiesexcl_(0),
55  //severityRecHitThreshold_(0),
56  //spId_(EcalSeverityLevelAlgo::kSwissCross),
57  //spIdThreshold_(0),
58  flags_(0) {
59  //set up the geometry and selector
60  const CaloGeometry* caloGeom = theCaloGeom_.product();
63 }
const CaloSubdetectorGeometry * subdet_[2]
edm::ESHandle< CaloGeometry > theCaloGeom_
std::vector< int > flags_
std::vector< int > severitiesexcl_
T const * product() const
Definition: ESHandle.h:86
const EcalRecHitCollection * ecalBarHits_
const EcalSeverityLevelAlgo * sevLevel_
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
const EcalRecHitCollection & caloHits_

◆ ~EgammaRecHitIsolation()

EgammaRecHitIsolation::~EgammaRecHitIsolation ( )

Definition at line 65 of file EgammaRecHitIsolation.cc.

65 {}

Member Function Documentation

◆ doFlagChecks()

void EgammaRecHitIsolation::doFlagChecks ( const std::vector< int > &  v)
inline

Definition at line 64 of file EgammaRecHitIsolation.h.

References flags_, jetUpdater_cfi::sort, and findQualityFiles::v.

Referenced by PhotonIsolationCalculator::calculateEcalRecHitIso().

64  {
65  flags_.clear();
66  flags_.insert(flags_.begin(), v.begin(), v.end());
67  std::sort(flags_.begin(), flags_.end());
68  }
std::vector< int > flags_

◆ doSeverityChecks()

void EgammaRecHitIsolation::doSeverityChecks ( const EcalRecHitCollection *const  recHits,
const std::vector< int > &  v 
)
inline

◆ getEnergySum() [1/2]

double EgammaRecHitIsolation::getEnergySum ( const reco::Candidate emObject,
EcalPFRecHitThresholds const &  thresholds 
) const
inline

Definition at line 44 of file EgammaRecHitIsolation.h.

References getSum_(), and particleFlowZeroSuppressionECAL_cff::thresholds.

44  {
45  return getSum_(emObject, false, &thresholds);
46  }
double getSum_(const reco::Candidate *, bool returnEt, const EcalPFRecHitThresholds *thresholds) const

◆ getEnergySum() [2/2]

double EgammaRecHitIsolation::getEnergySum ( const reco::SuperCluster emObject,
EcalPFRecHitThresholds const &  thresholds 
) const
inline

Definition at line 51 of file EgammaRecHitIsolation.h.

References getSum_(), and particleFlowZeroSuppressionECAL_cff::thresholds.

51  {
52  return getSum_(emObject, false, &thresholds);
53  }
double getSum_(const reco::Candidate *, bool returnEt, const EcalPFRecHitThresholds *thresholds) const

◆ getEtSum() [1/2]

double EgammaRecHitIsolation::getEtSum ( const reco::Candidate emObject,
EcalPFRecHitThresholds const &  thresholds 
) const
inline

Definition at line 41 of file EgammaRecHitIsolation.h.

References getSum_(), and particleFlowZeroSuppressionECAL_cff::thresholds.

Referenced by PhotonIsolationCalculator::calculateEcalRecHitIso(), and GsfElectronAlgo::createElectron().

41  {
42  return getSum_(emObject, true, &thresholds);
43  }
double getSum_(const reco::Candidate *, bool returnEt, const EcalPFRecHitThresholds *thresholds) const

◆ getEtSum() [2/2]

double EgammaRecHitIsolation::getEtSum ( const reco::SuperCluster emObject,
EcalPFRecHitThresholds const &  thresholds 
) const
inline

Definition at line 48 of file EgammaRecHitIsolation.h.

References getSum_(), and particleFlowZeroSuppressionECAL_cff::thresholds.

48  {
49  return getSum_(emObject, true, &thresholds);
50  }
double getSum_(const reco::Candidate *, bool returnEt, const EcalPFRecHitThresholds *thresholds) const

◆ getSum_() [1/2]

double EgammaRecHitIsolation::getSum_ ( const reco::Candidate emObject,
bool  returnEt,
const EcalPFRecHitThresholds thresholds 
) const
private

Definition at line 67 of file EgammaRecHitIsolation.cc.

References caloHits_, reco::deltaPhi(), ecalBarHits_, eLow_, edm::SortedCollection< T, SORT >::empty(), edm::SortedCollection< T, SORT >::end(), hcalRecHitTable_cff::energy, CastorDataFrameFilter_impl::energySum(), EgHLTOffHistBins_cfi::et, PVValHelper::eta, PV3DBase< T, PVType, FrameType >::eta(), etaSlice_, etLow_, extRadius_, spr::find(), edm::SortedCollection< T, SORT >::find(), flags_, reco::Candidate::get(), edm::Ref< C, T, F >::get(), CaloSubdetectorGeometry::getCells(), CaloGeometry::getGeometry(), mps_fire::i, intRadius_, dqmiolumiharvest::j, EcalRecHit::kGood, phi, diffTwoXMLs::r2, severitiesexcl_, EcalSeverityLevelAlgo::severityLevel(), sevLevel_, mathSSE::sqrt(), subdet_, theCaloGeom_, particleFlowZeroSuppressionECAL_cff::thresholds, useNumCrystals_, and vetoClustered_.

Referenced by getEnergySum(), and getEtSum().

69  {
70  double energySum = 0.;
71  if (!caloHits_.empty()) {
72  //Take the SC position
74  math::XYZPoint const& theCaloPosition = sc.get()->position();
75  GlobalPoint pclu(theCaloPosition.x(), theCaloPosition.y(), theCaloPosition.z());
76  float etaclus = pclu.eta();
77  float phiclus = pclu.phi();
78  float r2 = intRadius_ * intRadius_;
79 
80  std::vector<std::pair<DetId, float> >::const_iterator rhIt;
81 
82  for (int subdetnr = 0; subdetnr <= 1; subdetnr++) { // look in barrel and endcap
83  if (nullptr == subdet_[subdetnr])
84  continue;
85 
87  subdet_[subdetnr]->getCells(pclu, extRadius_); // select cells around cluster
89 
90  for (CaloSubdetectorGeometry::DetIdSet::const_iterator i = chosen.begin(); i != chosen.end();
91  ++i) { //loop selected cells
92  j = caloHits_.find(*i); // find selected cell among rechits
93  if (j != caloHits_.end()) { // add rechit only if available
94  auto cell = theCaloGeom_->getGeometry(*i);
95  float eta = cell->etaPos();
96  float phi = cell->phiPos();
97  float etaDiff = eta - etaclus;
98  float phiDiff = reco::deltaPhi(phi, phiclus);
99  float energy = j->energy();
100 
101  float rhThres = (thresholds != nullptr) ? (*thresholds)[j->detid()] : 0.f;
102  if (energy <= rhThres)
103  continue;
104 
105  if (useNumCrystals_) {
106  if (fabs(etaclus) < 1.479) { // Barrel num crystals, crystal width = 0.0174
107  if (fabs(etaDiff) < 0.0174 * etaSlice_)
108  continue;
109  //if (sqrt(etaDiff*etaDiff + phiDiff*phiDiff) < 0.0174*intRadius_)
110  //continue;
111  if ((etaDiff * etaDiff + phiDiff * phiDiff) < 0.00030276 * r2)
112  continue;
113  } else { // Endcap num crystals, crystal width = 0.00864*fabs(sinh(eta))
114  if (fabs(etaDiff) < 0.00864 * fabs(sinh(eta)) * etaSlice_)
115  continue;
116  //if (sqrt(etaDiff*etaDiff + phiDiff*phiDiff) < 0.00864*fabs(sinh(eta))*intRadius_)
117  // continue;
118  if ((etaDiff * etaDiff + phiDiff * phiDiff) < (0.000037325 * (cosh(2 * eta) - 1) * r2))
119  continue;
120  }
121  } else {
122  if (fabs(etaDiff) < etaSlice_)
123  continue; // jurassic strip cut
124  if (etaDiff * etaDiff + phiDiff * phiDiff < r2)
125  continue; // jurassic exclusion cone cut
126  }
127  //Check if RecHit is in SC
128  if (vetoClustered_) {
129  //Loop over basic clusters:
130  bool isClustered = false;
131  for (reco::CaloCluster_iterator bcIt = sc->clustersBegin(); bcIt != sc->clustersEnd(); ++bcIt) {
132  for (rhIt = (*bcIt)->hitsAndFractions().begin(); rhIt != (*bcIt)->hitsAndFractions().end(); ++rhIt) {
133  if (rhIt->first == *i)
134  isClustered = true;
135  if (isClustered)
136  break;
137  }
138 
139  if (isClustered)
140  break;
141  } //end loop over basic clusters
142 
143  if (isClustered)
144  continue;
145  } //end if removeClustered
146 
147  //std::cout << "detid " << j->detid() << std::endl;
148  int severityFlag = ecalBarHits_ == nullptr ? -1 : sevLevel_->severityLevel(j->detid(), *ecalBarHits_);
149  std::vector<int>::const_iterator sit =
150  std::find(severitiesexcl_.begin(), severitiesexcl_.end(), severityFlag);
151 
152  if (sit != severitiesexcl_.end())
153  continue;
154 
155  // new rechit flag checks
156  //std::vector<int>::const_iterator vit = std::find(flags_.begin(),
157  // flags_.end(),
158  // j->recoFlag());
159  //if (vit != flags_.end())
160  // continue;
161  if (!j->checkFlag(EcalRecHit::kGood)) {
162  if (j->checkFlags(flags_)) {
163  continue;
164  }
165  }
166 
167  float et = energy * std::sqrt(cell->getPosition().perp2() / cell->getPosition().mag2());
168  if (et > etLow_ && energy > eLow_) { //Changed energy --> fabs(energy) - now changed back to energy
169  if (returnEt)
170  energySum += et;
171  else
172  energySum += energy;
173  }
174 
175  } //End if not end of list
176  } //End loop over rechits
177  } //End loop over barrel/endcap
178  } //End if caloHits_
179 
180  return energySum;
181 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
const CaloSubdetectorGeometry * subdet_[2]
edm::ESHandle< CaloGeometry > theCaloGeom_
std::vector< int > flags_
T get() const
get a component
Definition: Candidate.h:221
T eta() const
Definition: PV3DBase.h:73
std::vector< EcalRecHit >::const_iterator const_iterator
std::vector< int > severitiesexcl_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id) const
Evaluate status from id use channelStatus from DB.
virtual DetIdSet getCells(const GlobalPoint &r, double dR) const
Get a list of all cells within a dR of the given cell.
T sqrt(T t)
Definition: SSEVec.h:19
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:60
const EcalRecHitCollection * ecalBarHits_
const_iterator end() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
iterator find(key_type k)
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
const EcalSeverityLevelAlgo * sevLevel_
double energySum(const DataFrame &df, int fs, int ls)
const EcalRecHitCollection & caloHits_

◆ getSum_() [2/2]

double EgammaRecHitIsolation::getSum_ ( const reco::SuperCluster sc,
bool  returnEt,
const EcalPFRecHitThresholds thresholds 
) const
private

Definition at line 183 of file EgammaRecHitIsolation.cc.

References caloHits_, reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), reco::deltaPhi(), ecalBarHits_, eLow_, edm::SortedCollection< T, SORT >::empty(), edm::SortedCollection< T, SORT >::end(), hcalRecHitTable_cff::energy, CastorDataFrameFilter_impl::energySum(), EgHLTOffHistBins_cfi::et, PVValHelper::eta, PV3DBase< T, PVType, FrameType >::eta(), etaSlice_, etLow_, extRadius_, spr::find(), edm::SortedCollection< T, SORT >::find(), flags_, CaloSubdetectorGeometry::getCells(), mps_fire::i, intRadius_, dqmiolumiharvest::j, EcalRecHit::kGood, phi, reco::CaloCluster::position(), position, edm::ESHandle< T >::product(), diffTwoXMLs::r2, severitiesexcl_, EcalSeverityLevelAlgo::severityLevel(), sevLevel_, subdet_, theCaloGeom_, particleFlowZeroSuppressionECAL_cff::thresholds, useNumCrystals_, and vetoClustered_.

185  {
186  double energySum = 0.;
187  if (!caloHits_.empty()) {
188  //Take the SC position
189 
190  const math::XYZPoint& theCaloPosition = sc->position();
191  GlobalPoint pclu(theCaloPosition.x(), theCaloPosition.y(), theCaloPosition.z());
192  double etaclus = pclu.eta();
193  double phiclus = pclu.phi();
194  double r2 = intRadius_ * intRadius_;
195 
196  std::vector<std::pair<DetId, float> >::const_iterator rhIt;
197 
198  for (int subdetnr = 0; subdetnr <= 1; subdetnr++) { // look in barrel and endcap
199  if (nullptr == subdet_[subdetnr])
200  continue;
202  subdet_[subdetnr]->getCells(pclu, extRadius_); // select cells around cluster
204  for (CaloSubdetectorGeometry::DetIdSet::const_iterator i = chosen.begin(); i != chosen.end();
205  ++i) { //loop selected cells
206 
207  j = caloHits_.find(*i); // find selected cell among rechits
208  if (j != caloHits_.end()) { // add rechit only if available
209  const GlobalPoint& position = (theCaloGeom_.product())->getPosition(*i);
210  double eta = position.eta();
211  double phi = position.phi();
212  double etaDiff = eta - etaclus;
213  double phiDiff = reco::deltaPhi(phi, phiclus);
214  double energy = j->energy();
215 
216  float rhThres = (thresholds != nullptr) ? (*thresholds)[j->detid()] : 0.f;
217  if (energy <= rhThres)
218  continue;
219 
220  if (useNumCrystals_) {
221  if (fabs(etaclus) < 1.479) { // Barrel num crystals, crystal width = 0.0174
222  if (fabs(etaDiff) < 0.0174 * etaSlice_)
223  continue;
224  // if ( sqrt(etaDiff*etaDiff + phiDiff*phiDiff) < 0.0174*intRadius_) continue;
225  if ((etaDiff * etaDiff + phiDiff * phiDiff) < 0.00030276 * r2)
226  continue;
227  } else { // Endcap num crystals, crystal width = 0.00864*fabs(sinh(eta))
228  if (fabs(etaDiff) < 0.00864 * fabs(sinh(eta)) * etaSlice_)
229  continue;
230  // if ( sqrt(etaDiff*etaDiff + phiDiff*phiDiff) < 0.00864*fabs(sinh(eta))*intRadius_) continue;
231  if ((etaDiff * etaDiff + phiDiff * phiDiff) < (0.000037325 * (cosh(2 * eta) - 1) * r2))
232  continue;
233  }
234  } else {
235  if (fabs(etaDiff) < etaSlice_)
236  continue; // jurassic strip cut
237  if (etaDiff * etaDiff + phiDiff * phiDiff < r2)
238  continue; // jurassic exclusion cone cut
239  }
240 
241  //Check if RecHit is in SC
242  if (vetoClustered_) {
243  //Loop over basic clusters:
244  bool isClustered = false;
245  for (reco::CaloCluster_iterator bcIt = sc->clustersBegin(); bcIt != sc->clustersEnd(); ++bcIt) {
246  for (rhIt = (*bcIt)->hitsAndFractions().begin(); rhIt != (*bcIt)->hitsAndFractions().end(); ++rhIt) {
247  if (rhIt->first == *i)
248  isClustered = true;
249  if (isClustered)
250  break;
251  }
252  if (isClustered)
253  break;
254  } //end loop over basic clusters
255 
256  if (isClustered)
257  continue;
258  } //end if removeClustered
259 
260  int severityFlag = sevLevel_->severityLevel(j->detid(), *ecalBarHits_);
261  std::vector<int>::const_iterator sit =
262  std::find(severitiesexcl_.begin(), severitiesexcl_.end(), severityFlag);
263 
264  if (sit != severitiesexcl_.end())
265  continue;
266 
267  // new rechit flag checks
268  //std::vector<int>::const_iterator vit = std::find(flags_.begin(),
269  // flags_.end(),
270  // j->recoFlag());
271  //if (vit != flags_.end())
272  // continue;
273  if (!j->checkFlag(EcalRecHit::kGood)) {
274  if (j->checkFlags(flags_)) {
275  continue;
276  }
277  }
278 
279  double et = energy * position.perp() / position.mag();
280  if (et > etLow_ && energy > eLow_) { //Changed energy --> fabs(energy) -- then changed into energy
281  if (returnEt)
282  energySum += et;
283  else
284  energySum += energy;
285  }
286 
287  } //End if not end of list
288  } //End loop over rechits
289  } //End loop over barrel/endcap
290  } //End if caloHits_
291 
292  return energySum;
293 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:153
const CaloSubdetectorGeometry * subdet_[2]
edm::ESHandle< CaloGeometry > theCaloGeom_
std::vector< int > flags_
T eta() const
Definition: PV3DBase.h:73
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:88
std::vector< EcalRecHit >::const_iterator const_iterator
std::vector< int > severitiesexcl_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:91
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id) const
Evaluate status from id use channelStatus from DB.
T const * product() const
Definition: ESHandle.h:86
virtual DetIdSet getCells(const GlobalPoint &r, double dR) const
Get a list of all cells within a dR of the given cell.
const EcalRecHitCollection * ecalBarHits_
const_iterator end() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
iterator find(key_type k)
static int position[264][3]
Definition: ReadPGInfo.cc:289
const EcalSeverityLevelAlgo * sevLevel_
double energySum(const DataFrame &df, int fs, int ls)
const EcalRecHitCollection & caloHits_

◆ setUseNumCrystals()

void EgammaRecHitIsolation::setUseNumCrystals ( bool  b = true)
inline

Definition at line 55 of file EgammaRecHitIsolation.h.

References b, and useNumCrystals_.

Referenced by PhotonIsolationCalculator::calculateEcalRecHitIso().

55 { useNumCrystals_ = b; }
double b
Definition: hdecay.h:120

◆ setVetoClustered()

void EgammaRecHitIsolation::setVetoClustered ( bool  b = true)
inline

Definition at line 56 of file EgammaRecHitIsolation.h.

References b, and vetoClustered_.

Referenced by PhotonIsolationCalculator::calculateEcalRecHitIso().

56 { vetoClustered_ = b; }
double b
Definition: hdecay.h:120

Member Data Documentation

◆ caloHits_

const EcalRecHitCollection& EgammaRecHitIsolation::caloHits_
private

Definition at line 84 of file EgammaRecHitIsolation.h.

Referenced by getSum_().

◆ ecalBarHits_

const EcalRecHitCollection* EgammaRecHitIsolation::ecalBarHits_
private

Definition at line 89 of file EgammaRecHitIsolation.h.

Referenced by doSeverityChecks(), and getSum_().

◆ eLow_

double EgammaRecHitIsolation::eLow_
private

Definition at line 81 of file EgammaRecHitIsolation.h.

Referenced by getSum_().

◆ etaSlice_

double EgammaRecHitIsolation::etaSlice_
private

Definition at line 79 of file EgammaRecHitIsolation.h.

Referenced by getSum_().

◆ etLow_

double EgammaRecHitIsolation::etLow_
private

Definition at line 80 of file EgammaRecHitIsolation.h.

Referenced by getSum_().

◆ extRadius_

double EgammaRecHitIsolation::extRadius_
private

Definition at line 77 of file EgammaRecHitIsolation.h.

Referenced by getSum_().

◆ flags_

std::vector<int> EgammaRecHitIsolation::flags_
private

Definition at line 91 of file EgammaRecHitIsolation.h.

Referenced by doFlagChecks(), and getSum_().

◆ intRadius_

double EgammaRecHitIsolation::intRadius_
private

Definition at line 78 of file EgammaRecHitIsolation.h.

Referenced by getSum_().

◆ severitiesexcl_

std::vector<int> EgammaRecHitIsolation::severitiesexcl_
private

Definition at line 90 of file EgammaRecHitIsolation.h.

Referenced by doSeverityChecks(), and getSum_().

◆ sevLevel_

const EcalSeverityLevelAlgo* EgammaRecHitIsolation::sevLevel_
private

Definition at line 85 of file EgammaRecHitIsolation.h.

Referenced by getSum_().

◆ subdet_

const CaloSubdetectorGeometry* EgammaRecHitIsolation::subdet_[2]
private

Definition at line 93 of file EgammaRecHitIsolation.h.

Referenced by EgammaRecHitIsolation(), and getSum_().

◆ theCaloGeom_

edm::ESHandle<CaloGeometry> EgammaRecHitIsolation::theCaloGeom_
private

Definition at line 83 of file EgammaRecHitIsolation.h.

Referenced by EgammaRecHitIsolation(), and getSum_().

◆ useNumCrystals_

bool EgammaRecHitIsolation::useNumCrystals_
private

Definition at line 87 of file EgammaRecHitIsolation.h.

Referenced by getSum_(), and setUseNumCrystals().

◆ vetoClustered_

bool EgammaRecHitIsolation::vetoClustered_
private

Definition at line 88 of file EgammaRecHitIsolation.h.

Referenced by getSum_(), and setVetoClustered().