Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:56

0001 #include "RecoEgamma/EgammaIsolationAlgos/interface/HcalPFClusterIsolation.h"
0002 
0003 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0004 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0005 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0006 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0007 
0008 #include <DataFormats/Math/interface/deltaR.h>
0009 
0010 template <typename T1>
0011 HcalPFClusterIsolation<T1>::HcalPFClusterIsolation(double drMax,
0012                                                    double drVetoBarrel,
0013                                                    double drVetoEndcap,
0014                                                    double etaStripBarrel,
0015                                                    double etaStripEndcap,
0016                                                    double energyBarrel,
0017                                                    double energyEndcap,
0018                                                    bool useEt)
0019     : drMax_(drMax),
0020       drVetoBarrel_(drVetoBarrel),
0021       drVetoEndcap_(drVetoEndcap),
0022       etaStripBarrel_(etaStripBarrel),
0023       etaStripEndcap_(etaStripEndcap),
0024       energyBarrel_(energyBarrel),
0025       energyEndcap_(energyEndcap),
0026       useEt_(useEt) {}
0027 
0028 template <typename T1>
0029 HcalPFClusterIsolation<T1>::~HcalPFClusterIsolation() {}
0030 
0031 template <typename T1>
0032 double HcalPFClusterIsolation<T1>::getSum(const T1 cand,
0033                                           const std::vector<edm::Handle<reco::PFClusterCollection>>& clusterHandles) {
0034   double etSum = 0.;
0035   double candAbsEta = std::abs(cand.eta());
0036 
0037   float etaStrip = 0;
0038   float dRVeto = 0;
0039   if (candAbsEta < 1.479) {
0040     dRVeto = drVetoBarrel_;
0041     etaStrip = etaStripBarrel_;
0042   } else {
0043     dRVeto = drVetoEndcap_;
0044     etaStrip = etaStripEndcap_;
0045   }
0046 
0047   for (unsigned int nHandle = 0; nHandle < clusterHandles.size(); nHandle++) {
0048     for (unsigned i = 0; i < clusterHandles[nHandle]->size(); i++) {
0049       const reco::PFClusterRef pfclu(clusterHandles[nHandle], i);
0050 
0051       if (candAbsEta < 1.479) {
0052         if (std::abs(pfclu->pt()) < energyBarrel_)
0053           continue;
0054       } else {
0055         if (std::abs(pfclu->energy()) < energyEndcap_)
0056           continue;
0057       }
0058 
0059       float dEta = std::abs(cand.eta() - pfclu->eta());
0060       if (dEta < etaStrip)
0061         continue;
0062 
0063       float dR2 = deltaR2(cand.eta(), cand.phi(), pfclu->eta(), pfclu->phi());
0064       if (dR2 > (drMax_ * drMax_) || dR2 < (dRVeto * dRVeto))
0065         continue;
0066 
0067       if (useEt_)
0068         etSum += pfclu->pt();
0069       else
0070         etSum += pfclu->energy();
0071     }
0072   }
0073 
0074   return etSum;
0075 }
0076 
0077 template <typename T1>
0078 double HcalPFClusterIsolation<T1>::getSum(T1Ref ref,
0079                                           const std::vector<edm::Handle<reco::PFClusterCollection>>& clusterHandles) {
0080   return getSum(*ref, clusterHandles);
0081 }
0082 
0083 template class HcalPFClusterIsolation<reco::RecoEcalCandidate>;
0084 template class HcalPFClusterIsolation<reco::RecoChargedCandidate>;
0085 template class HcalPFClusterIsolation<reco::Photon>;
0086 template class HcalPFClusterIsolation<reco::GsfElectron>;