#include <RecoEgamma/EgammaIsolationAlgos/interface/ElectronTkIsolation.h>
Public Member Functions | |
ElectronTkIsolation (double extRadius, double intRadius, double ptLow, double lip, const reco::TrackCollection *) | |
int | getNumberTracks (const reco::GsfElectron *) const |
double | getPtTracks (const reco::GsfElectron *) const |
~ElectronTkIsolation () | |
Private Member Functions | |
std::pair< int, double > | getIso (const reco::GsfElectron *) const |
Private Attributes | |
double | extRadius_ |
double | intRadius_ |
double | lip_ |
double | ptLow_ |
const reco::TrackCollection * | trackCollection_ |
Definition at line 24 of file ElectronTkIsolation.h.
ElectronTkIsolation::ElectronTkIsolation | ( | double | extRadius, | |
double | intRadius, | |||
double | ptLow, | |||
double | lip, | |||
const reco::TrackCollection * | trackCollection | |||
) |
Definition at line 27 of file ElectronTkIsolation.cc.
00031 : 00032 extRadius_(extRadius), 00033 intRadius_(intRadius), 00034 ptLow_(ptLow), 00035 lip_(lip), 00036 trackCollection_(trackCollection) 00037 { 00038 }
ElectronTkIsolation::~ElectronTkIsolation | ( | ) |
std::pair< int, double > ElectronTkIsolation::getIso | ( | const reco::GsfElectron * | electron | ) | const [private] |
Definition at line 45 of file ElectronTkIsolation.cc.
References counter(), extRadius_, reco::GsfElectron::gsfTrack(), intRadius_, lip_, and ptLow_.
Referenced by getNumberTracks(), and getPtTracks().
00046 { 00047 int counter =0 ; 00048 double ptSum =0.; 00049 //Take the electron track 00050 reco::GsfTrackRef tmpTrack = electron->gsfTrack() ; 00051 math::XYZVector tmpElectronMomentumAtVtx = (*tmpTrack).momentum () ; 00052 00053 for ( reco::TrackCollection::const_iterator itrTr = (*trackCollection_).begin() ; 00054 itrTr != (*trackCollection_).end() ; 00055 ++itrTr ) 00056 { 00057 math::XYZVector tmpTrackMomentumAtVtx = (*itrTr).momentum () ; 00058 double this_pt = (*itrTr).pt(); 00059 if ( this_pt < ptLow_ ) 00060 continue ; 00061 if (fabs( (*itrTr).dz() - (*tmpTrack).dz() ) > lip_ ) 00062 continue ; 00063 double dr = DeltaR(tmpTrackMomentumAtVtx,tmpElectronMomentumAtVtx) ; 00064 if ( fabs(dr) < extRadius_ && 00065 fabs(dr) >= intRadius_ ) 00066 { 00067 ++counter ; 00068 ptSum += this_pt; 00069 } 00070 }//end loop over tracks 00071 00072 std::pair<int,double> retval; 00073 retval.first = counter; 00074 retval.second = ptSum; 00075 00076 return retval; 00077 }
int ElectronTkIsolation::getNumberTracks | ( | const reco::GsfElectron * | electron | ) | const |
Definition at line 80 of file ElectronTkIsolation.cc.
References getIso().
Referenced by EgammaElectronTkNumIsolationProducer::produce().
00081 { 00082 //counter for the tracks in the isolation cone 00083 return getIso(electron).first ; 00084 }
double ElectronTkIsolation::getPtTracks | ( | const reco::GsfElectron * | electron | ) | const |
Definition at line 86 of file ElectronTkIsolation.cc.
References getIso().
Referenced by EgammaElectronTkIsolationProducer::produce().
00087 { 00088 return getIso(electron).second ; 00089 }
double ElectronTkIsolation::extRadius_ [private] |
double ElectronTkIsolation::intRadius_ [private] |
double ElectronTkIsolation::lip_ [private] |
double ElectronTkIsolation::ptLow_ [private] |
const reco::TrackCollection* ElectronTkIsolation::trackCollection_ [private] |
Definition at line 48 of file ElectronTkIsolation.h.