Identifies pile-up candidates from a collection of PFCandidates, and produces the corresponding collection of PileUpCandidates. More...
#include <PFPileUp.h>
Public Member Functions | |
virtual void | beginJob () |
PFPileUp (const edm::ParameterSet &) | |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
~PFPileUp () | |
Private Member Functions | |
reco::VertexRef | chargedHadronVertex (const edm::Handle< reco::VertexCollection > &vertices, const reco::PFCandidate &pfcand) const |
Private Attributes | |
bool | checkClosestZVertex_ |
use the closest z vertex if a track is not in a vertex | |
bool | enable_ |
enable PFPileUp selection | |
edm::InputTag | inputTagPFCandidates_ |
PFCandidates to be analyzed. | |
edm::InputTag | inputTagVertices_ |
vertices | |
bool | verbose_ |
verbose ? |
Identifies pile-up candidates from a collection of PFCandidates, and produces the corresponding collection of PileUpCandidates.
Definition at line 31 of file PFPileUp.h.
PFPileUp::PFPileUp | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 19 of file PFPileUp.cc.
References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().
{ inputTagPFCandidates_ = iConfig.getParameter<InputTag>("PFCandidates"); inputTagVertices_ = iConfig.getParameter<InputTag>("Vertices"); enable_ = iConfig.getParameter<bool>("Enable"); verbose_ = iConfig.getUntrackedParameter<bool>("verbose",false); if ( iConfig.exists("checkClosestZVertex") ) { checkClosestZVertex_ = iConfig.getParameter<bool>("checkClosestZVertex"); } else { checkClosestZVertex_ = false; } produces<reco::PileUpPFCandidateCollection>(); }
PFPileUp::~PFPileUp | ( | ) |
Definition at line 48 of file PFPileUp.cc.
{ }
void PFPileUp::beginJob | ( | void | ) | [virtual] |
VertexRef PFPileUp::chargedHadronVertex | ( | const edm::Handle< reco::VertexCollection > & | vertices, |
const reco::PFCandidate & | pfcand | ||
) | const [private] |
Definition at line 111 of file PFPileUp.cc.
References getHLTprescales::index, reco::PFCandidate::trackRef(), reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), and reco::PFCandidate::vertex().
{ reco::TrackBaseRef trackBaseRef( pfcand.trackRef() ); size_t iVertex = 0; unsigned index=0; unsigned nFoundVertex = 0; typedef reco::VertexCollection::const_iterator IV; float bestweight=0; for(IV iv=vertices->begin(); iv!=vertices->end(); ++iv, ++index) { const reco::Vertex& vtx = *iv; typedef reco::Vertex::trackRef_iterator IT; // 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 float w = vtx.trackWeight(baseRef); if(baseRef == trackBaseRef ) { //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 VertexRef( vertices, 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 VertexRef( vertices, iVertex); } return VertexRef(); }
void PFPileUp::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 55 of file PFPileUp.cc.
References edm::Event::getByLabel(), h, i, edm::Ref< C, T, F >::isNull(), edm::Ref< C, T, F >::key(), reco::tau::pfCandidates(), and edm::Event::put().
{ // LogDebug("PFPileUp")<<"START event: "<<iEvent.id().event() // <<" in run "<<iEvent.id().run()<<endl; // get PFCandidates auto_ptr< reco::PileUpPFCandidateCollection > pOutput( new reco::PileUpPFCandidateCollection ); if(enable_) { Handle<PFCandidateCollection> pfCandidates; iEvent.getByLabel( inputTagPFCandidates_, pfCandidates); // get vertices Handle<VertexCollection> vertices; iEvent.getByLabel( inputTagVertices_, vertices); for( unsigned i=0; i<pfCandidates->size(); i++ ) { // const reco::PFCandidate& cand = (*pfCandidates)[i]; PFCandidatePtr candptr(pfCandidates, i); // PFCandidateRef pfcandref(pfCandidates,i); VertexRef vertexref; switch( candptr->particleId() ) { case PFCandidate::h: vertexref = chargedHadronVertex( vertices, *candptr ); break; default: continue; } // no associated vertex, or primary vertex // not pile-up if( vertexref.isNull() || vertexref.key()==0 ) continue; pOutput->push_back( PileUpPFCandidate( candptr, vertexref ) ); pOutput->back().setSourceCandidatePtr( candptr ); } } iEvent.put( pOutput ); }
bool PFPileUp::checkClosestZVertex_ [private] |
use the closest z vertex if a track is not in a vertex
Definition at line 62 of file PFPileUp.h.
bool PFPileUp::enable_ [private] |
enable PFPileUp selection
Definition at line 56 of file PFPileUp.h.
edm::InputTag PFPileUp::inputTagPFCandidates_ [private] |
PFCandidates to be analyzed.
Definition at line 50 of file PFPileUp.h.
edm::InputTag PFPileUp::inputTagVertices_ [private] |
vertices
Definition at line 53 of file PFPileUp.h.
bool PFPileUp::verbose_ [private] |
verbose ?
Definition at line 59 of file PFPileUp.h.