File indexing completed on 2024-04-06 12:24:55
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "RecoEgamma/EgammaIsolationAlgos/interface/EcalPFClusterIsolation.h"
0010
0011 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0012 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0013 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0014 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0015
0016 #include <DataFormats/Math/interface/deltaR.h>
0017
0018 template <typename T1>
0019 EcalPFClusterIsolation<T1>::EcalPFClusterIsolation(double drMax,
0020 double drVetoBarrel,
0021 double drVetoEndcap,
0022 double etaStripBarrel,
0023 double etaStripEndcap,
0024 double energyBarrel,
0025 double energyEndcap)
0026 : drMax_(drMax),
0027 drVetoBarrel_(drVetoBarrel),
0028 drVetoEndcap_(drVetoEndcap),
0029 etaStripBarrel_(etaStripBarrel),
0030 etaStripEndcap_(etaStripEndcap),
0031 energyBarrel_(energyBarrel),
0032 energyEndcap_(energyEndcap) {}
0033
0034 template <typename T1>
0035 EcalPFClusterIsolation<T1>::~EcalPFClusterIsolation() {}
0036
0037 template <typename T1>
0038 bool EcalPFClusterIsolation<T1>::computedRVeto(T1Ref candRef, reco::PFClusterRef pfclu) {
0039 float dR2 = deltaR2(candRef->eta(), candRef->phi(), pfclu->eta(), pfclu->phi());
0040 if (dR2 > (drMax_ * drMax_))
0041 return false;
0042
0043 if (candRef->superCluster().isNonnull()) {
0044
0045 for (reco::CaloCluster_iterator it = candRef->superCluster()->clustersBegin();
0046 it != candRef->superCluster()->clustersEnd();
0047 ++it) {
0048 if ((*it)->seed() == pfclu->seed()) {
0049 return false;
0050 }
0051 }
0052 }
0053
0054 return true;
0055 }
0056
0057 template <>
0058 bool EcalPFClusterIsolation<reco::RecoChargedCandidate>::computedRVeto(T1Ref candRef, reco::PFClusterRef pfclu) {
0059 float dR2 = deltaR2(candRef->eta(), candRef->phi(), pfclu->eta(), pfclu->phi());
0060 if (dR2 > (drMax_ * drMax_) || dR2 < drVeto2_)
0061 return false;
0062 else
0063 return true;
0064 }
0065
0066 template <typename T1>
0067 double EcalPFClusterIsolation<T1>::getSum(const T1 cand, edm::Handle<reco::PFClusterCollection> clusterHandle) {
0068 drVeto2_ = -1.;
0069 float etaStrip = -1;
0070
0071 if (std::abs(cand.eta()) < 1.479) {
0072 drVeto2_ = drVetoBarrel_ * drVetoBarrel_;
0073 etaStrip = etaStripBarrel_;
0074 } else {
0075 drVeto2_ = drVetoEndcap_ * drVetoEndcap_;
0076 etaStrip = etaStripEndcap_;
0077 }
0078
0079 float etSum = 0;
0080 for (size_t i = 0; i < clusterHandle->size(); i++) {
0081 reco::PFClusterRef pfclu(clusterHandle, i);
0082
0083 if (std::abs(cand.eta()) < 1.479) {
0084 if (std::abs(pfclu->pt()) < energyBarrel_)
0085 continue;
0086 } else {
0087 if (std::abs(pfclu->energy()) < energyEndcap_)
0088 continue;
0089 }
0090
0091 float dEta = std::abs(cand.eta() - pfclu->eta());
0092 if (dEta < etaStrip)
0093 continue;
0094 if (not computedRVeto(cand, pfclu))
0095 continue;
0096
0097 etSum += pfclu->pt();
0098 }
0099
0100 return etSum;
0101 }
0102
0103 template <typename T1>
0104 double EcalPFClusterIsolation<T1>::getSum(T1Ref ref, edm::Handle<std::vector<reco::PFCluster> > clusts) {
0105 return getSum(*ref, clusts);
0106 }
0107
0108 template <typename T1>
0109 bool EcalPFClusterIsolation<T1>::computedRVeto(T1 cand, reco::PFClusterRef pfclu) {
0110 float dR2 = deltaR2(cand.eta(), cand.phi(), pfclu->eta(), pfclu->phi());
0111 if (dR2 > (drMax_ * drMax_))
0112 return false;
0113
0114 if (cand.superCluster().isNonnull()) {
0115
0116 for (reco::CaloCluster_iterator it = cand.superCluster()->clustersBegin(); it != cand.superCluster()->clustersEnd();
0117 ++it) {
0118 if ((*it)->seed() == pfclu->seed()) {
0119 return false;
0120 }
0121 }
0122 }
0123
0124 return true;
0125 }
0126
0127 template class EcalPFClusterIsolation<reco::RecoEcalCandidate>;
0128 template class EcalPFClusterIsolation<reco::RecoChargedCandidate>;
0129 template class EcalPFClusterIsolation<reco::Photon>;
0130 template class EcalPFClusterIsolation<reco::GsfElectron>;