#include <RecoEgamma/PhotonIdentification/interface/PhotonIDAlgo.h>
Public Member Functions | |
void | baseSetup (const edm::ParameterSet &conf) |
double | calculateBasicClusterIso (const reco::Photon *photon, const edm::Event &iEvent, double RCone=0.4, double RConeInner=0, double etMin=0) |
double | calculateEcalRecHitIso (const reco::Photon *photon, const edm::Event &iEvent, const edm::EventSetup &iSetup, double RCone, double RConeInner, double etaSlice, double etMin) |
double | calculateHcalRecHitIso (const reco::Photon *photon, const edm::Event &iEvent, const edm::EventSetup &iSetup, double RCone, double RConeInner, double etaSlice, double etMin) |
double | calculateR9 (const reco::Photon *photon, const edm::Event &iEvent, const edm::EventSetup &iSetup) |
void | calculateTrackIso (const reco::Photon *photon, const edm::Event &e, double &trkCone, int &ntrkCone, double pTThresh=0, double RCone=.4, double RinnerCone=.1) |
void | classify (const reco::Photon *photon, bool &isEBPho, bool &isEEPho, bool &isEBGap, bool &isEEGap, bool &isEBEEGap) |
bool | isAlsoElectron (const reco::Photon *photon, const edm::Event &e) |
PhotonIDAlgo () | |
virtual | ~PhotonIDAlgo () |
Protected Attributes | |
std::string | barrelbasicclusterCollection_ |
std::string | barrelbasicclusterProducer_ |
std::string | barrelecalCollection_ |
std::string | barrelecalProducer_ |
std::string | barrelsuperclusterCollection_ |
std::string | barrelsuperclusterProducer_ |
std::string | endcapbasicclusterCollection_ |
std::string | endcapbasicclusterProducer_ |
std::string | endcapecalCollection_ |
std::string | endcapecalProducer_ |
std::string | endcapsuperclusterCollection_ |
std::string | endcapSuperClusterProducer_ |
edm::InputTag | gsfRecoInputTag_ |
std::string | hcalCollection_ |
std::string | hcalProducer_ |
std::vector< double > | moduleEtaBoundary_ |
double | modulePhiBoundary_ |
edm::InputTag | trackInputTag_ |
Definition at line 11 of file PhotonIDAlgo.h.
PhotonIDAlgo::PhotonIDAlgo | ( | ) | [inline] |
virtual PhotonIDAlgo::~PhotonIDAlgo | ( | ) | [inline, virtual] |
void PhotonIDAlgo::baseSetup | ( | const edm::ParameterSet & | conf | ) |
Definition at line 37 of file PhotonIDAlgo.cc.
References barrelecalCollection_, barrelecalProducer_, endcapecalCollection_, endcapecalProducer_, edm::ParameterSet::getParameter(), gsfRecoInputTag_, hcalCollection_, hcalProducer_, moduleEtaBoundary_, modulePhiBoundary_, and trackInputTag_.
Referenced by CutBasedPhotonIDAlgo::setup().
00037 { 00038 00039 00040 trackInputTag_ = conf.getParameter<edm::InputTag>("trackProducer"); 00041 00042 barrelecalCollection_ = conf.getParameter<std::string>("barrelEcalRecHitCollection"); 00043 barrelecalProducer_ = conf.getParameter<std::string>("barrelEcalRecHitProducer"); 00044 endcapecalCollection_ = conf.getParameter<std::string>("endcapEcalRecHitCollection"); 00045 endcapecalProducer_ = conf.getParameter<std::string>("endcapEcalRecHitProducer"); 00046 hcalCollection_ = conf.getParameter<std::string>("HcalRecHitCollection"); 00047 hcalProducer_ = conf.getParameter<std::string>("HcalRecHitProducer"); 00048 00049 gsfRecoInputTag_ = conf.getParameter<edm::InputTag>("GsfRecoCollection"); 00050 00051 modulePhiBoundary_ = conf.getParameter<double>("modulePhiBoundary"); 00052 moduleEtaBoundary_ = conf.getParameter<std::vector<double> >("moduleEtaBoundary"); 00053 00054 }
double PhotonIDAlgo::calculateBasicClusterIso | ( | const reco::Photon * | photon, | |
const edm::Event & | iEvent, | |||
double | RCone = 0.4 , |
|||
double | RConeInner = 0 , |
|||
double | etMin = 0 | |||
) |
Definition at line 137 of file PhotonIDAlgo.cc.
References barrelbasicclusterCollection_, barrelbasicclusterProducer_, barrelsuperclusterCollection_, barrelsuperclusterProducer_, endcapbasicclusterCollection_, endcapbasicclusterProducer_, endcapsuperclusterCollection_, endcapSuperClusterProducer_, edm::Event::getByLabel(), EgammaEcalIsolation::getEcalEtSum(), reco::Particle::p4(), and edm::Handle< T >::product().
00142 { 00143 00144 00145 edm::Handle<reco::BasicClusterCollection> basicClusterH; 00146 edm::Handle<reco::SuperClusterCollection> endcapSuperClusterH; 00147 00148 double peta = photon->p4().Eta(); 00149 if (fabs(peta) > 1.479){ 00150 iEvent.getByLabel(endcapbasicclusterProducer_,endcapbasicclusterCollection_,basicClusterH); 00151 iEvent.getByLabel(endcapSuperClusterProducer_,endcapsuperclusterCollection_,endcapSuperClusterH); 00152 } 00153 else{ 00154 iEvent.getByLabel(barrelbasicclusterProducer_,barrelbasicclusterCollection_,basicClusterH); 00155 iEvent.getByLabel(barrelsuperclusterProducer_,barrelsuperclusterCollection_,endcapSuperClusterH); 00156 } 00157 const reco::BasicClusterCollection* basicClusterCollection_ = basicClusterH.product(); 00158 const reco::SuperClusterCollection* endcapSuperClusterCollection_ = endcapSuperClusterH.product(); 00159 00160 double ecalIsol=0.; 00161 EgammaEcalIsolation phoIso(RCone,etMin, basicClusterCollection_, endcapSuperClusterCollection_); 00162 ecalIsol = phoIso.getEcalEtSum(photon); 00163 // delete phoIso; 00164 00165 return ecalIsol; 00166 00167 00168 }
double PhotonIDAlgo::calculateEcalRecHitIso | ( | const reco::Photon * | photon, | |
const edm::Event & | iEvent, | |||
const edm::EventSetup & | iSetup, | |||
double | RCone, | |||
double | RConeInner, | |||
double | etaSlice, | |||
double | etMin | |||
) |
Definition at line 170 of file PhotonIDAlgo.cc.
References barrelecalCollection_, barrelecalProducer_, DetId::Ecal, endcapecalCollection_, endcapecalProducer_, edm::EventSetup::get(), edm::Event::getByLabel(), EgammaRecHitIsolation::getEtSum(), edm::Handle< T >::product(), and reco::Photon::superCluster().
Referenced by CutBasedPhotonIDAlgo::calculate().
00176 { 00177 00178 00179 edm::Handle<EcalRecHitCollection> ecalhitsCollH; 00180 00181 double peta = photon->superCluster()->position().eta(); 00182 if (fabs(peta) > 1.479){ 00183 iEvent.getByLabel(endcapecalProducer_,endcapecalCollection_, ecalhitsCollH); 00184 } 00185 else{ 00186 iEvent.getByLabel(barrelecalProducer_,barrelecalCollection_, ecalhitsCollH); 00187 } 00188 const EcalRecHitCollection* rechitsCollection_ = ecalhitsCollH.product(); 00189 00190 std::auto_ptr<CaloRecHitMetaCollectionV> RecHits(0); 00191 RecHits = std::auto_ptr<CaloRecHitMetaCollectionV>(new EcalRecHitMetaCollection(*rechitsCollection_)); 00192 00193 edm::ESHandle<CaloGeometry> geoHandle; 00194 iSetup.get<CaloGeometryRecord>().get(geoHandle); 00195 double ecalIsol=0.; 00196 00197 00198 EgammaRecHitIsolation phoIso(RCone, 00199 RConeInner, 00200 etaSlice, 00201 etMin, 00202 geoHandle, 00203 &(*RecHits), 00204 DetId::Ecal); 00205 ecalIsol = phoIso.getEtSum(photon); 00206 // delete phoIso; 00207 00208 return ecalIsol; 00209 00210 00211 }
double PhotonIDAlgo::calculateHcalRecHitIso | ( | const reco::Photon * | photon, | |
const edm::Event & | iEvent, | |||
const edm::EventSetup & | iSetup, | |||
double | RCone, | |||
double | RConeInner, | |||
double | etaSlice, | |||
double | etMin | |||
) |
Definition at line 255 of file PhotonIDAlgo.cc.
References edm::EventSetup::get(), edm::Event::getByLabel(), EgammaRecHitIsolation::getEtSum(), DetId::Hcal, hcalCollection_, hcalProducer_, and edm::Handle< T >::product().
Referenced by CutBasedPhotonIDAlgo::calculate().
00261 { 00262 00263 00264 edm::Handle<HBHERecHitCollection> hcalhitsCollH; 00265 00266 iEvent.getByLabel(hcalProducer_,hcalCollection_, hcalhitsCollH); 00267 00268 const HBHERecHitCollection* rechitsCollection_ = hcalhitsCollH.product(); 00269 00270 std::auto_ptr<CaloRecHitMetaCollectionV> RecHits(0); 00271 RecHits = std::auto_ptr<CaloRecHitMetaCollectionV>(new HBHERecHitMetaCollection(*rechitsCollection_)); 00272 00273 edm::ESHandle<CaloGeometry> geoHandle; 00274 iSetup.get<CaloGeometryRecord>().get(geoHandle); 00275 double ecalIsol=0.; 00276 00277 00278 EgammaRecHitIsolation phoIso(RCone, 00279 RConeInner, 00280 etMin, 00281 etaSlice, 00282 geoHandle, 00283 &(*RecHits), 00284 DetId::Hcal); 00285 ecalIsol = phoIso.getEtSum(photon); 00286 // delete phoIso; 00287 00288 return ecalIsol; 00289 00290 00291 }
double PhotonIDAlgo::calculateR9 | ( | const reco::Photon * | photon, | |
const edm::Event & | iEvent, | |||
const edm::EventSetup & | iSetup | |||
) |
Definition at line 213 of file PhotonIDAlgo.cc.
References barrelecalCollection_, barrelecalProducer_, EcalClusterTools::e3x3(), endcapecalCollection_, endcapecalProducer_, edm::EventSetup::get(), edm::Ref< C, T, F >::get(), edm::Event::getByLabel(), edm::ESHandle< T >::product(), edm::Handle< T >::product(), r9, reco::SuperCluster::seed(), and reco::Photon::superCluster().
Referenced by CutBasedPhotonIDAlgo::calculate().
00216 { 00217 00218 00219 edm::Handle<EcalRecHitCollection> ecalhitsCollH; 00220 double peta = photon->superCluster()->position().eta(); 00221 edm::ESHandle<CaloGeometry> geoHandle; 00222 iSetup.get<CaloGeometryRecord>().get(geoHandle); 00223 //const CaloGeometry& geometry = *geoHandle; 00224 edm::ESHandle<CaloTopology> pTopology; 00225 iSetup.get<CaloTopologyRecord>().get(pTopology); 00226 const CaloTopology *topology = pTopology.product(); 00227 // const CaloSubdetectorGeometry *geometry_p; 00228 if (fabs(peta) > 1.479){ 00229 iEvent.getByLabel(endcapecalProducer_,endcapecalCollection_, ecalhitsCollH); 00230 } 00231 else{ 00232 iEvent.getByLabel(barrelecalProducer_,barrelecalCollection_, ecalhitsCollH); 00233 } 00234 const EcalRecHitCollection* rechitsCollection_ = ecalhitsCollH.product(); 00235 00236 reco::SuperClusterRef scref = photon->superCluster(); 00237 const reco::SuperCluster *sc = scref.get(); 00238 const reco::BasicClusterRef bcref = sc->seed(); 00239 const reco::BasicCluster *bc = bcref.get(); 00240 00241 if (fabs(peta) > 1.479){ 00242 00243 float e3x3 = EcalClusterTools::e3x3(*bc, rechitsCollection_, topology); 00244 double r9 = e3x3 / (photon->superCluster()->rawEnergy()); 00245 return r9; 00246 } 00247 else{ 00248 float e3x3 = EcalClusterTools::e3x3(*bc, rechitsCollection_, topology); 00249 double r9 = e3x3 / (photon->superCluster()->rawEnergy()); 00250 return r9; 00251 } 00252 00253 }
void PhotonIDAlgo::calculateTrackIso | ( | const reco::Photon * | photon, | |
const edm::Event & | e, | |||
double & | trkCone, | |||
int & | ntrkCone, | |||
double | pTThresh = 0 , |
|||
double | RCone = .4 , |
|||
double | RinnerCone = .1 | |||
) |
Definition at line 105 of file PhotonIDAlgo.cc.
References counter(), edm::Event::getByLabel(), PhotonTkIsolation::getNumberTracks(), PhotonTkIsolation::getPtTracks(), edm::Handle< T >::isValid(), edm::Handle< T >::product(), trackInputTag_, and tracks.
Referenced by CutBasedPhotonIDAlgo::calculate().
00111 { 00112 int counter =0; 00113 double ptSum =0.; 00114 00115 00116 //get the tracks 00117 edm::Handle<reco::TrackCollection> tracks; 00118 e.getByLabel(trackInputTag_,tracks); 00119 if(!tracks.isValid()) { 00120 return; 00121 } 00122 const reco::TrackCollection* trackCollection = tracks.product(); 00123 //Photon Eta and Phi. Hope these are correct. 00124 00125 00126 PhotonTkIsolation phoIso(RCone, RinnerCone, pTThresh, 2., trackCollection); 00127 counter = phoIso.getNumberTracks(photon); 00128 ptSum = phoIso.getPtTracks(photon); 00129 //delete phoIso; 00130 00131 ntrkCone = counter; 00132 trkCone = ptSum; 00133 }
void PhotonIDAlgo::classify | ( | const reco::Photon * | photon, | |
bool & | isEBPho, | |||
bool & | isEEPho, | |||
bool & | isEBGap, | |||
bool & | isEEGap, | |||
bool & | isEBEEGap | |||
) |
Definition at line 58 of file PhotonIDAlgo.cc.
References eta, i, moduleEtaBoundary_, modulePhiBoundary_, phi, Pi, and reco::Photon::superCluster().
Referenced by CutBasedPhotonIDAlgo::calculate().
00063 { 00064 00065 //Set fiducial flags for this photon. 00066 double eta = photon->superCluster()->position().eta(); 00067 double phi = photon->superCluster()->position().phi(); 00068 double feta = fabs(eta); 00069 00070 //Are you in the Ecal Endcap (EE)? 00071 if(feta>1.479) 00072 isEEPho = true; 00073 else 00074 isEBPho = true; 00075 00076 // Set fiducial flags (isEBGap, isEEGap... 00077 00078 //Are you in the gap between EE and Ecal Barrel (EB)? 00079 if (fabs(feta-1.479)<.1) isEBEEGap=true; 00080 00081 // Set isEBGap if photon is closer than "modulePhiBoundary_" (set in cfg) 00082 // to a phi module/supermodule boundary (same thing) 00083 if (phi < 0) phi += TMath::Pi()*2.; 00084 Float_t phiRelative = fmod( phi , 20*TMath::Pi()/180 ) - 10*TMath::Pi()/180; 00085 if ( fabs(phiRelative) < modulePhiBoundary_ ) isEBGap=true; 00086 00087 // Set isEBGap if photon is between specific eta values 00088 // in the "moduleEtaBoundary_" variable. 00089 // Loop over the vector of Eta boundaries given in the config file 00090 bool nearEtaBoundary = false; 00091 for (unsigned int i=0; i < moduleEtaBoundary_.size(); i+=2) { 00092 // Checks to see if it's between the 0th and 1st entry, the 2nd and 3rd entry...etc 00093 if ( (feta > moduleEtaBoundary_[i]) && (feta < moduleEtaBoundary_[i+1]) ) { 00094 //std::cout << "Photon between eta " << moduleEtaBoundary_[i] << " and " << moduleEtaBoundary_[i+1] << std::endl; 00095 nearEtaBoundary = true; 00096 break; 00097 } 00098 } 00099 00100 // If it's near an eta boundary and in the Barrel 00101 if (nearEtaBoundary) isEBGap=true; 00102 00103 }
bool PhotonIDAlgo::isAlsoElectron | ( | const reco::Photon * | photon, | |
const edm::Event & | e | |||
) |
Definition at line 295 of file PhotonIDAlgo.cc.
Referenced by CutBasedPhotonIDAlgo::calculate().
00296 { 00297 00298 //Currently some instability with GsfPixelMatchElectronCollection 00299 //causes this simple code to die horribly. Thus we simply return false 00300 //for now. 00301 00302 //Get MY supercluster position 00303 // std::cout << "Checking isAlsoElectron code: " << std::endl; 00304 // reco::SuperClusterRef sc = photon->superCluster(); 00305 // float PhoCaloE = sc.get()->energy(); 00306 00307 // math::XYZVector position(sc.get()->position().x(), 00308 // sc.get()->position().y(), 00309 // sc.get()->position().z()); 00310 00311 // std::cout << "Got supercluster position: Photon." << std::endl; 00312 //get the Gsf electrons: 00313 // edm::Handle<reco::PixelMatchGsfElectronCollection> pElectrons; 00314 // e.getByLabel(gsfRecoInputTag_, pElectrons); 00315 // std::cout << "Got GsfElectronCollection: " << std::endl; 00316 // float PhoCaloE=0; 00317 00318 // const reco::PixelMatchGsfElectronCollection *elec = pElectrons.product(); 00319 // for(reco::PixelMatchGsfElectronCollection::const_iterator gItr = elec->begin(); gItr != elec->end(); ++gItr){ 00320 00321 // std::cout << "Got Electron: " << std::endl; 00322 // float EleCaloE = gItr->caloEnergy(); 00323 // std::cout << "Energy: " << EleCaloE << std::endl; 00324 // std::cout << "Photon E: " << PhoCaloE << std::endl; 00325 // float dE = fabs(EleCaloE-PhoCaloE); 00326 // std::cout << "Made comparison. " << std::endl; 00327 00328 // if(dE < 0.0001) return true; 00329 00330 // } 00331 00332 return false; 00333 }
std::string PhotonIDAlgo::barrelbasicclusterCollection_ [protected] |
std::string PhotonIDAlgo::barrelbasicclusterProducer_ [protected] |
std::string PhotonIDAlgo::barrelecalCollection_ [protected] |
Definition at line 73 of file PhotonIDAlgo.h.
Referenced by baseSetup(), calculateEcalRecHitIso(), and calculateR9().
std::string PhotonIDAlgo::barrelecalProducer_ [protected] |
Definition at line 74 of file PhotonIDAlgo.h.
Referenced by baseSetup(), calculateEcalRecHitIso(), and calculateR9().
std::string PhotonIDAlgo::barrelsuperclusterCollection_ [protected] |
std::string PhotonIDAlgo::barrelsuperclusterProducer_ [protected] |
std::string PhotonIDAlgo::endcapbasicclusterCollection_ [protected] |
std::string PhotonIDAlgo::endcapbasicclusterProducer_ [protected] |
std::string PhotonIDAlgo::endcapecalCollection_ [protected] |
Definition at line 75 of file PhotonIDAlgo.h.
Referenced by baseSetup(), calculateEcalRecHitIso(), and calculateR9().
std::string PhotonIDAlgo::endcapecalProducer_ [protected] |
Definition at line 76 of file PhotonIDAlgo.h.
Referenced by baseSetup(), calculateEcalRecHitIso(), and calculateR9().
std::string PhotonIDAlgo::endcapsuperclusterCollection_ [protected] |
std::string PhotonIDAlgo::endcapSuperClusterProducer_ [protected] |
edm::InputTag PhotonIDAlgo::gsfRecoInputTag_ [protected] |
std::string PhotonIDAlgo::hcalCollection_ [protected] |
Definition at line 77 of file PhotonIDAlgo.h.
Referenced by baseSetup(), and calculateHcalRecHitIso().
std::string PhotonIDAlgo::hcalProducer_ [protected] |
Definition at line 78 of file PhotonIDAlgo.h.
Referenced by baseSetup(), and calculateHcalRecHitIso().
std::vector<double> PhotonIDAlgo::moduleEtaBoundary_ [protected] |
double PhotonIDAlgo::modulePhiBoundary_ [protected] |
edm::InputTag PhotonIDAlgo::trackInputTag_ [protected] |