#include <ElectronSqPtTkIsolation.h>
Public Member Functions | |
ElectronSqPtTkIsolation (double extRadius, double intRadius, double ptLow, double lip, const reco::TrackCollection *) | |
int | getNumberTracks (const reco::GsfElectron *) const |
double | getPtTracks (const reco::GsfElectron *) const |
~ElectronSqPtTkIsolation () | |
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 17 of file ElectronSqPtTkIsolation.h.
ElectronSqPtTkIsolation::ElectronSqPtTkIsolation | ( | double | extRadius, |
double | intRadius, | ||
double | ptLow, | ||
double | lip, | ||
const reco::TrackCollection * | trackCollection | ||
) |
Definition at line 20 of file ElectronSqPtTkIsolation.cc.
: extRadius_(extRadius), intRadius_(intRadius), ptLow_(ptLow), lip_(lip), trackCollection_(trackCollection) { }
ElectronSqPtTkIsolation::~ElectronSqPtTkIsolation | ( | ) |
Definition at line 33 of file ElectronSqPtTkIsolation.cc.
{ }
std::pair< int, double > ElectronSqPtTkIsolation::getIso | ( | const reco::GsfElectron * | electron | ) | const [private] |
Definition at line 38 of file ElectronSqPtTkIsolation.cc.
References extRadius_, reco::GsfElectron::gsfTrack(), intRadius_, lip_, and ptLow_.
Referenced by getNumberTracks(), and getPtTracks().
{ int counter =0 ; double ptSum =0.; //Take the electron track reco::GsfTrackRef tmpTrack = electron->gsfTrack() ; math::XYZVector tmpElectronMomentumAtVtx = (*tmpTrack).momentum () ; for ( reco::TrackCollection::const_iterator itrTr = (*trackCollection_).begin() ; itrTr != (*trackCollection_).end() ; ++itrTr ) { math::XYZVector tmpTrackMomentumAtVtx = (*itrTr).momentum () ; double this_pt = (*itrTr).pt(); if ( this_pt < ptLow_ ) continue ; if (fabs( (*itrTr).dz() - (*tmpTrack).dz() ) > lip_ ) continue ; double dr = DeltaR(tmpTrackMomentumAtVtx,tmpElectronMomentumAtVtx) ; if ( fabs(dr) < extRadius_ && fabs(dr) >= intRadius_ ) { ++counter ; ptSum += this_pt*this_pt; // sum of squared pT } }//end loop over tracks std::pair<int,double> retval; retval.first = counter; retval.second = ptSum; return retval; }
int ElectronSqPtTkIsolation::getNumberTracks | ( | const reco::GsfElectron * | electron | ) | const |
Definition at line 73 of file ElectronSqPtTkIsolation.cc.
References getIso().
{ //counter for the tracks in the isolation cone return getIso(electron).first ; }
double ElectronSqPtTkIsolation::getPtTracks | ( | const reco::GsfElectron * | electron | ) | const |
Definition at line 79 of file ElectronSqPtTkIsolation.cc.
References getIso().
Referenced by ElectronSqPtTkIsolationProducer::produce().
{ return getIso(electron).second ; }
double ElectronSqPtTkIsolation::extRadius_ [private] |
Definition at line 36 of file ElectronSqPtTkIsolation.h.
Referenced by getIso().
double ElectronSqPtTkIsolation::intRadius_ [private] |
Definition at line 37 of file ElectronSqPtTkIsolation.h.
Referenced by getIso().
double ElectronSqPtTkIsolation::lip_ [private] |
Definition at line 39 of file ElectronSqPtTkIsolation.h.
Referenced by getIso().
double ElectronSqPtTkIsolation::ptLow_ [private] |
Definition at line 38 of file ElectronSqPtTkIsolation.h.
Referenced by getIso().
const reco::TrackCollection* ElectronSqPtTkIsolation::trackCollection_ [private] |
Definition at line 41 of file ElectronSqPtTkIsolation.h.