#include <RecoJets/FFTJetProducers/plugins/FFTJetPFPileupCleaner.cc>
Description: cleans up a collection of partice flow objects
Implementation: [Notes on implementation]
Definition at line 42 of file FFTJetPFPileupCleaner.cc.
FFTJetPFPileupCleaner::FFTJetPFPileupCleaner | ( | const edm::ParameterSet & | ps | ) | [explicit] |
Definition at line 110 of file FFTJetPFPileupCleaner.cc.
References buildRemovalMask().
: init_param(edm::InputTag, PFCandidates), init_param(edm::InputTag, Vertices), init_param(bool, checkClosestZVertex), init_param(bool, removeMainVertex), init_param(bool, removeUnassociated), init_param(bool, reverseRemovalDecision), init_param(bool, remove_X ), init_param(bool, remove_h ), init_param(bool, remove_e ), init_param(bool, remove_mu ), init_param(bool, remove_gamma ), init_param(bool, remove_h0 ), init_param(bool, remove_h_HF ), init_param(bool, remove_egamma_HF), removalMask(0), init_param(double, etaMin), init_param(double, etaMax), init_param(double, vertexNdofCut) { buildRemovalMask(); produces<reco::PFCandidateCollection>(); }
FFTJetPFPileupCleaner::~FFTJetPFPileupCleaner | ( | ) |
Definition at line 135 of file FFTJetPFPileupCleaner.cc.
{ }
FFTJetPFPileupCleaner::FFTJetPFPileupCleaner | ( | ) | [private] |
FFTJetPFPileupCleaner::FFTJetPFPileupCleaner | ( | const FFTJetPFPileupCleaner & | ) | [private] |
void FFTJetPFPileupCleaner::beginJob | ( | void | ) | [protected, virtual] |
void FFTJetPFPileupCleaner::buildRemovalMask | ( | ) | [private] |
Definition at line 302 of file FFTJetPFPileupCleaner.cc.
References reco::PFCandidate::e, reco::PFCandidate::egamma_HF, reco::PFCandidate::gamma, reco::PFCandidate::h, reco::PFCandidate::h0, reco::PFCandidate::h_HF, reco::PFCandidate::mu, remove_e, remove_egamma_HF, remove_gamma, remove_h, remove_h0, remove_h_HF, remove_mu, remove_X, setRemovalBit(), and reco::PFCandidate::X.
Referenced by FFTJetPFPileupCleaner().
{ setRemovalBit(reco::PFCandidate::X, remove_X ); setRemovalBit(reco::PFCandidate::h, remove_h ); setRemovalBit(reco::PFCandidate::e, remove_e ); setRemovalBit(reco::PFCandidate::mu, remove_mu ); setRemovalBit(reco::PFCandidate::gamma, remove_gamma ); setRemovalBit(reco::PFCandidate::h0, remove_h0 ); setRemovalBit(reco::PFCandidate::h_HF, remove_h_HF ); setRemovalBit(reco::PFCandidate::egamma_HF, remove_egamma_HF); }
void FFTJetPFPileupCleaner::endJob | ( | void | ) | [protected, virtual] |
reco::VertexRef FFTJetPFPileupCleaner::findSomeVertex | ( | const edm::Handle< reco::VertexCollection > & | vertices, |
const reco::PFCandidate & | pfcand | ||
) | const [private] |
Definition at line 218 of file FFTJetPFPileupCleaner.cc.
References checkClosestZVertex, getHLTprescales::index, reco::PFCandidate::trackRef(), reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), reco::PFCandidate::vertex(), vertexNdofCut, and w().
Referenced by produce().
{ typedef reco::VertexCollection::const_iterator IV; typedef reco::Vertex::trackRef_iterator IT; reco::TrackBaseRef trackBaseRef(pfcand.trackRef()); size_t iVertex = 0; unsigned index = 0; unsigned nFoundVertex = 0; double bestweight = 0; const IV vertend(vertices->end()); for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index) { const double ndof = iv->ndof(); if (!iv->isFake() && ndof > vertexNdofCut) { const reco::Vertex& vtx = *iv; // loop on tracks in vertices IT trackend(vtx.tracks_end()); for (IT iTrack=vtx.tracks_begin(); iTrack!=trackend; ++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) { // select the vertex for which the track has the highest weight const double w = vtx.trackWeight(baseRef); 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 vertices. " << "Used to be an assert"; return reco::VertexRef(vertices, iVertex); } // 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!=vertend; ++iv, ++index) { const double ndof = iv->ndof(); if (!iv->isFake() && ndof > vertexNdofCut) { const double dz = fabs(ztrack - iv->z()); if (dz < dzmin) { dzmin = dz; iVertex = index; foundVertex = true; } } } if (foundVertex) return reco::VertexRef(vertices, iVertex); } return reco::VertexRef(); }
bool FFTJetPFPileupCleaner::isRemovable | ( | reco::PFCandidate::ParticleType | ptype | ) | const [private] |
Definition at line 194 of file FFTJetPFPileupCleaner.cc.
References removalMask, and edm::shift.
Referenced by produce().
{ const unsigned shift = static_cast<unsigned>(ptype); assert(shift < 32U); return removalMask & (1U << shift); }
FFTJetPFPileupCleaner& FFTJetPFPileupCleaner::operator= | ( | const FFTJetPFPileupCleaner & | ) | [private] |
void FFTJetPFPileupCleaner::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [protected, virtual] |
Implements edm::EDProducer.
Definition at line 141 of file FFTJetPFPileupCleaner.cc.
References eta(), etaMax, findSomeVertex(), edm::Event::getByLabel(), i, edm::Ref< C, T, F >::isNull(), isRemovable(), edm::Ref< C, T, F >::key(), PFCandidates, reco::tau::pfCandidates(), edm::Event::put(), removeMainVertex, removeUnassociated, reverseRemovalDecision, reco::PFCandidate::setSourceCandidatePtr(), and Vertices.
{ // get PFCandidates std::auto_ptr<reco::PFCandidateCollection> pOutput(new reco::PFCandidateCollection); edm::Handle<reco::PFCandidateCollection> pfCandidates; iEvent.getByLabel(PFCandidates, pfCandidates); // get vertices edm::Handle<reco::VertexCollection> vertices; iEvent.getByLabel(Vertices, vertices); const unsigned ncand = pfCandidates->size(); for (unsigned i=0; i<ncand; ++i) { reco::PFCandidatePtr candptr(pfCandidates, i); bool remove = false; if (isRemovable(candptr->particleId())) { reco::VertexRef vertexref(findSomeVertex(vertices, *candptr)); if (vertexref.isNull()) remove = removeUnassociated; else if (vertexref.key() == 0) remove = removeMainVertex; else remove = true; } // Check the eta range if (!remove) { const double eta = candptr->p4().Eta(); remove = eta < etaMin || eta > etaMax; } // Should we remember this candidate? if (reverseRemovalDecision) remove = !remove; if (!remove) { const reco::PFCandidate& cand = (*pfCandidates)[i]; pOutput->push_back(cand); pOutput->back().setSourceCandidatePtr(candptr); } } iEvent.put(pOutput); }
void FFTJetPFPileupCleaner::setRemovalBit | ( | reco::PFCandidate::ParticleType | ptype, |
bool | onOff | ||
) | [private] |
Definition at line 203 of file FFTJetPFPileupCleaner.cc.
References removalMask, and edm::shift.
Referenced by buildRemovalMask().
{ const unsigned shift = static_cast<unsigned>(ptype); assert(shift < 32U); const unsigned mask = (1U << shift); if (value) removalMask |= mask; else removalMask &= ~mask; }
bool FFTJetPFPileupCleaner::checkClosestZVertex [private] |
Definition at line 72 of file FFTJetPFPileupCleaner.cc.
Referenced by findSomeVertex().
double FFTJetPFPileupCleaner::etaMax [private] |
Definition at line 101 of file FFTJetPFPileupCleaner.cc.
Referenced by produce().
double FFTJetPFPileupCleaner::etaMin [private] |
Definition at line 100 of file FFTJetPFPileupCleaner.cc.
Definition at line 67 of file FFTJetPFPileupCleaner.cc.
Referenced by produce().
unsigned FFTJetPFPileupCleaner::removalMask [private] |
Definition at line 97 of file FFTJetPFPileupCleaner.cc.
Referenced by isRemovable(), and setRemovalBit().
bool FFTJetPFPileupCleaner::remove_e [private] |
Definition at line 89 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::remove_egamma_HF [private] |
Definition at line 94 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::remove_gamma [private] |
Definition at line 91 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::remove_h [private] |
Definition at line 88 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::remove_h0 [private] |
Definition at line 92 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::remove_h_HF [private] |
Definition at line 93 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::remove_mu [private] |
Definition at line 90 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::remove_X [private] |
Definition at line 87 of file FFTJetPFPileupCleaner.cc.
Referenced by buildRemovalMask().
bool FFTJetPFPileupCleaner::removeMainVertex [private] |
Definition at line 76 of file FFTJetPFPileupCleaner.cc.
Referenced by produce().
bool FFTJetPFPileupCleaner::removeUnassociated [private] |
Definition at line 80 of file FFTJetPFPileupCleaner.cc.
Referenced by produce().
bool FFTJetPFPileupCleaner::reverseRemovalDecision [private] |
Definition at line 83 of file FFTJetPFPileupCleaner.cc.
Referenced by produce().
double FFTJetPFPileupCleaner::vertexNdofCut [private] |
Definition at line 104 of file FFTJetPFPileupCleaner.cc.
Referenced by findSomeVertex().
edm::InputTag FFTJetPFPileupCleaner::Vertices [private] |
Definition at line 68 of file FFTJetPFPileupCleaner.cc.
Referenced by produce().