#include <PFPileUpAlgo.h>
Public Member Functions | |
int | chargedHadronVertex (const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand) const |
const reco::PFCandidateCollection & | getPFCandidatesFromPU () const |
const reco::PFCandidateCollection & | getPFCandidatesFromVtx () const |
PFPileUpAlgo (bool checkClosestZVertex, bool verbose=false) | |
PFPileUpAlgo () | |
void | process (const reco::PFCandidateCollection &pfCandidates, const reco::VertexCollection &vertices, const edm::Handle< reco::PFCandidateCollection > *handle=0) |
void | setCheckClosestZVertex (bool val) |
void | setVerbose (bool verbose) |
~PFPileUpAlgo () | |
Private Attributes | |
bool | checkClosestZVertex_ |
use the closest z vertex if a track is not in a vertex | |
reco::PFCandidateCollection | pfCandidatesFromPU_ |
reco::PFCandidateCollection | pfCandidatesFromVtx_ |
bool | verbose_ |
verbose ? |
Definition at line 12 of file PFPileUpAlgo.h.
PFPileUpAlgo::PFPileUpAlgo | ( | ) | [inline] |
Definition at line 14 of file PFPileUpAlgo.h.
:checkClosestZVertex_(true), verbose_(false) {;}
PFPileUpAlgo::PFPileUpAlgo | ( | bool | checkClosestZVertex, |
bool | verbose = false |
||
) | [inline] |
Definition at line 16 of file PFPileUpAlgo.h.
: checkClosestZVertex_(checkClosestZVertex), verbose_(verbose) {;}
PFPileUpAlgo::~PFPileUpAlgo | ( | ) | [inline] |
Definition at line 19 of file PFPileUpAlgo.h.
{;}
int PFPileUpAlgo::chargedHadronVertex | ( | const reco::VertexCollection & | vertices, |
const reco::PFCandidate & | pfcand | ||
) | const |
Definition at line 44 of file PFPileUpAlgo.cc.
References checkClosestZVertex_, getHLTprescales::index, reco::PFCandidate::trackRef(), reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), reco::PFCandidate::vertex(), and w().
Referenced by process().
{ reco::TrackBaseRef trackBaseRef( pfcand.trackRef() ); size_t iVertex = 0; unsigned index=0; unsigned nFoundVertex = 0; typedef reco::VertexCollection::const_iterator IV; typedef reco::Vertex::trackRef_iterator IT; float bestweight=0; for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) { const reco::Vertex& vtx = *iv; // loop on tracks in vertices for(IT iTrack=vtx.tracks_begin(); iTrack!=vtx.tracks_end(); ++iTrack) { const reco::TrackBaseRef& baseRef = *iTrack; // one of the tracks in the vertex is the same as // the track considered in the function if(baseRef == trackBaseRef ) { float w = vtx.trackWeight(baseRef); //select the vertex for which the track has the highest weight if (w > bestweight){ bestweight=w; iVertex=index; nFoundVertex++; } } } } if (nFoundVertex>0){ if (nFoundVertex!=1) edm::LogWarning("TrackOnTwoVertex")<<"a track is shared by at least two verteces. Used to be an assert"; return iVertex; } // no vertex found with this track. // optional: as a secondary solution, associate the closest vertex in z if ( checkClosestZVertex_ ) { double dzmin = 10000; double ztrack = pfcand.vertex().z(); bool foundVertex = false; index = 0; for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) { double dz = fabs(ztrack - iv->z()); if(dz<dzmin) { dzmin = dz; iVertex = index; foundVertex = true; } } if( foundVertex ) return iVertex; } return -1 ; }
const reco::PFCandidateCollection& PFPileUpAlgo::getPFCandidatesFromPU | ( | ) | const [inline] |
Definition at line 30 of file PFPileUpAlgo.h.
References pfCandidatesFromPU_.
{return pfCandidatesFromPU_;}
const reco::PFCandidateCollection& PFPileUpAlgo::getPFCandidatesFromVtx | ( | ) | const [inline] |
Definition at line 32 of file PFPileUpAlgo.h.
References pfCandidatesFromVtx_.
{return pfCandidatesFromVtx_;}
void PFPileUpAlgo::process | ( | const reco::PFCandidateCollection & | pfCandidates, |
const reco::VertexCollection & | vertices, | ||
const edm::Handle< reco::PFCandidateCollection > * | handle = 0 |
||
) |
Definition at line 5 of file PFPileUpAlgo.cc.
References chargedHadronVertex(), reco::PFCandidate::h, i, reco::PFCandidate::particleId(), pfCandidatesFromPU_, and pfCandidatesFromVtx_.
{ pfCandidatesFromVtx_.clear(); pfCandidatesFromPU_.clear(); for( unsigned i=0; i<pfCandidates.size(); i++ ) { const reco::PFCandidate& cand = pfCandidates[i]; int ivertex; switch( cand.particleId() ) { case reco::PFCandidate::h: ivertex = chargedHadronVertex( vertices, cand ); break; default: continue; } // no associated vertex, or primary vertex // not pile-up if( ivertex == -1 || ivertex == 0 ) { pfCandidatesFromVtx_.push_back(cand); pfCandidatesFromVtx_.back().setSourceCandidatePtr( reco::PFCandidatePtr(*handle,i) ); } else { // associated to a vertex pfCandidatesFromPU_.push_back(cand); if(handle) { pfCandidatesFromPU_.back().setSourceCandidatePtr( reco::PFCandidatePtr(*handle,i) ); } } } }
void PFPileUpAlgo::setCheckClosestZVertex | ( | bool | val | ) | [inline] |
Definition at line 28 of file PFPileUpAlgo.h.
References checkClosestZVertex_.
{ checkClosestZVertex_ = val;}
void PFPileUpAlgo::setVerbose | ( | bool | verbose | ) | [inline] |
Definition at line 26 of file PFPileUpAlgo.h.
References validate_alignment_devdb10_cfg::verbose, and verbose_.
bool PFPileUpAlgo::checkClosestZVertex_ [private] |
use the closest z vertex if a track is not in a vertex
Definition at line 41 of file PFPileUpAlgo.h.
Referenced by chargedHadronVertex(), and setCheckClosestZVertex().
Definition at line 48 of file PFPileUpAlgo.h.
Referenced by getPFCandidatesFromPU(), and process().
Definition at line 47 of file PFPileUpAlgo.h.
Referenced by getPFCandidatesFromVtx(), and process().
bool PFPileUpAlgo::verbose_ [private] |