CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
FFTJetPFPileupCleaner Class Reference

#include <RecoJets/FFTJetProducers/plugins/FFTJetPFPileupCleaner.cc>

Inheritance diagram for FFTJetPFPileupCleaner:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 FFTJetPFPileupCleaner (const edm::ParameterSet &)
 
 ~FFTJetPFPileupCleaner ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Protected Member Functions

void beginJob ()
 
void endJob ()
 
void produce (edm::Event &, const edm::EventSetup &)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Private Member Functions

void buildRemovalMask ()
 
 FFTJetPFPileupCleaner ()
 
 FFTJetPFPileupCleaner (const FFTJetPFPileupCleaner &)
 
reco::VertexRef findSomeVertexWFakes (const edm::Handle< reco::VertexCollection > &vertices, const edm::Handle< reco::VertexCollection > &fakeVertices, const reco::PFCandidate &pfcand, bool *fromFakeSet) const
 
bool isAcceptableVtx (reco::VertexCollection::const_iterator iv) const
 
bool isRemovable (reco::PFCandidate::ParticleType ptype) const
 
FFTJetPFPileupCleaneroperator= (const FFTJetPFPileupCleaner &)
 
void setRemovalBit (reco::PFCandidate::ParticleType ptype, bool onOff)
 

Private Attributes

bool checkClosestZVertex
 
double etaMax
 
double etaMin
 
edm::InputTag FakePrimaryVertices
 
bool keepIfPVneighbor
 
edm::InputTag PFCandidates
 
unsigned removalMask
 
bool remove_e
 
bool remove_egamma_HF
 
bool remove_gamma
 
bool remove_h
 
bool remove_h0
 
bool remove_h_HF
 
bool remove_mu
 
bool remove_X
 
bool removeMainVertex
 
bool removeUnassociated
 
bool reverseRemovalDecision
 
bool useFakePrimaryVertex
 
double vertexNdofCut
 
double vertexZmaxCut
 
edm::InputTag Vertices
 
std::vector< std::pair< double,
unsigned > > 
zAssoc
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Detailed Description

Description: cleans up a collection of partice flow objects

Implementation: [Notes on implementation]

Definition at line 46 of file FFTJetPFPileupCleaner.cc.

Constructor & Destructor Documentation

FFTJetPFPileupCleaner::FFTJetPFPileupCleaner ( const edm::ParameterSet ps)
explicit

Definition at line 134 of file FFTJetPFPileupCleaner.cc.

References buildRemovalMask().

144  init_param(bool, remove_X ),
145  init_param(bool, remove_h ),
146  init_param(bool, remove_e ),
147  init_param(bool, remove_mu ),
148  init_param(bool, remove_gamma ),
149  init_param(bool, remove_h0 ),
150  init_param(bool, remove_h_HF ),
152  removalMask(0),
153  init_param(double, etaMin),
154  init_param(double, etaMax),
155  init_param(double, vertexNdofCut),
156  init_param(double, vertexZmaxCut)
157 {
159  produces<reco::PFCandidateCollection>();
160 }
#define init_param(type, varname)
FFTJetPFPileupCleaner::~FFTJetPFPileupCleaner ( )

Definition at line 163 of file FFTJetPFPileupCleaner.cc.

164 {
165 }
FFTJetPFPileupCleaner::FFTJetPFPileupCleaner ( )
private
FFTJetPFPileupCleaner::FFTJetPFPileupCleaner ( const FFTJetPFPileupCleaner )
private

Member Function Documentation

void FFTJetPFPileupCleaner::beginJob ( void  )
protectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 501 of file FFTJetPFPileupCleaner.cc.

502 {
503 }
void FFTJetPFPileupCleaner::buildRemovalMask ( )
private

Definition at line 487 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().

488 {
497 }
void setRemovalBit(reco::PFCandidate::ParticleType ptype, bool onOff)
void FFTJetPFPileupCleaner::endJob ( void  )
protectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 507 of file FFTJetPFPileupCleaner.cc.

508 {
509 }
reco::VertexRef FFTJetPFPileupCleaner::findSomeVertexWFakes ( const edm::Handle< reco::VertexCollection > &  vertices,
const edm::Handle< reco::VertexCollection > &  fakeVertices,
const reco::PFCandidate pfcand,
bool *  fromFakeSet 
) const
private

Definition at line 282 of file FFTJetPFPileupCleaner.cc.

References abs, checkClosestZVertex, first, getHLTprescales::index, isAcceptableVtx(), keepIfPVneighbor, python.multivaluedict::sort(), reco::PFCandidate::trackRef(), reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), useFakePrimaryVertex, reco::PFCandidate::vertex(), w(), detailsBasic3DVector::z, and zAssoc.

Referenced by produce().

286 {
287  typedef reco::VertexCollection::const_iterator IV;
289 
290  *fromFakeSet = false;
291  reco::TrackBaseRef trackBaseRef(pfcand.trackRef());
292 
293  size_t iVertex = 0;
294  unsigned nFoundVertex = 0;
295  const IV vertend(vertices->end());
296 
297  {
298  unsigned index = 0;
299  double bestweight = 0.0;
300  for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
301  if (isAcceptableVtx(iv))
302  {
303  const reco::Vertex& vtx = *iv;
304 
305  // loop on tracks in vertices
306  IT trackend(vtx.tracks_end());
307  for (IT iTrack=vtx.tracks_begin(); iTrack!=trackend; ++iTrack)
308  {
309  const reco::TrackBaseRef& baseRef = *iTrack;
310 
311  // one of the tracks in the vertex is the same as
312  // the track considered in the function
313  if (baseRef == trackBaseRef)
314  {
315  // select the vertex for which the track has the highest weight
316  const double w = vtx.trackWeight(baseRef);
317  if (w > bestweight)
318  {
319  bestweight=w;
320  iVertex=index;
321  nFoundVertex++;
322  }
323  }
324  }
325  }
326  }
327 
328  if (nFoundVertex > 0)
329  {
330  if (nFoundVertex != 1)
331  edm::LogWarning("TrackOnTwoVertex")
332  << "a track is shared by at least two vertices. "
333  << "Used to be an assert";
334 
335  // Check if we can re-associate this track with one
336  // of the fake vertices
338  {
339  const double ztrack = pfcand.vertex().z();
340  double dzmin = std::abs(ztrack - ((*vertices)[iVertex]).z());
341 
342  const IV fakeEnd(fakeVertices->end());
343  unsigned index = 0;
344  for (IV iv=fakeVertices->begin(); iv!=fakeEnd; ++iv, ++index)
345  if (isAcceptableVtx(iv))
346  {
347  const double dz = std::abs(ztrack - iv->z());
348  if (dz < dzmin)
349  {
350  dzmin = dz;
351  iVertex = index;
352  *fromFakeSet = true;
353  }
354  }
355  }
356 
357  return reco::VertexRef(*fromFakeSet ? fakeVertices : vertices, iVertex);
358  }
359 
360  // optional: as a secondary solution, associate the closest vertex in z
361  if (checkClosestZVertex)
362  {
363  const double ztrack = pfcand.vertex().z();
364  bool foundVertex = false;
365 
366  if (keepIfPVneighbor)
367  {
368  // Sort all vertices according to their Z coordinate
369  zAssoc.clear();
370  unsigned index = 0;
371  for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
372  if (isAcceptableVtx(iv))
373  zAssoc.push_back(std::pair<double,unsigned>(iv->z(), index));
374  const unsigned numRealVertices = index;
375 
376  // Mix the fake vertex collection into zAssoc.
377  // Note that we do not reset "index" before doing this.
379  {
380  const IV fakeEnd(fakeVertices->end());
381  for (IV iv=fakeVertices->begin(); iv!=fakeEnd; ++iv, ++index)
382  if (isAcceptableVtx(iv))
383  zAssoc.push_back(std::pair<double,unsigned>(iv->z(), index));
384  }
385 
386  // Check where the track z position fits into this sequence
387  if (!zAssoc.empty())
388  {
389  std::sort(zAssoc.begin(), zAssoc.end());
390  std::pair<double,unsigned> tPair(ztrack, UINT_MAX);
391  const unsigned iAbove = std::upper_bound(
392  zAssoc.begin(), zAssoc.end(), tPair) - zAssoc.begin();
393 
394  // Check whether one of the vertices with indices
395  // iAbove or (iAbove - 1) is a primary vertex.
396  // If so, return it. Otherwise return the one
397  // with closest distance to the track.
398  unsigned ich[2] = {0U, 0U};
399  unsigned nch = 1;
400  if (iAbove)
401  {
402  ich[0] = iAbove - 1U;
403  ich[1] = iAbove;
404  if (iAbove < zAssoc.size())
405  nch = 2;
406  }
407 
408  double dzmin = 1.0e100;
409  unsigned bestVertexNum = UINT_MAX;
410  for (unsigned icheck=0; icheck<nch; ++icheck)
411  {
412  const unsigned zAssocIndex = ich[icheck];
413  const unsigned vertexNum = zAssoc[zAssocIndex].second;
414 
415  if (vertexNum == numRealVertices ||
416  (!useFakePrimaryVertex && vertexNum == 0U))
417  {
418  bestVertexNum = vertexNum;
419  break;
420  }
421 
422  const double dz = std::abs(ztrack - zAssoc[zAssocIndex].first);
423  if (dz < dzmin)
424  {
425  dzmin = dz;
426  bestVertexNum = vertexNum;
427  }
428  }
429 
430  foundVertex = bestVertexNum < UINT_MAX;
431  if (foundVertex)
432  {
433  iVertex = bestVertexNum;
434  if (iVertex >= numRealVertices)
435  {
436  *fromFakeSet = true;
437  iVertex -= numRealVertices;
438  }
439  }
440  }
441  }
442  else
443  {
444  // This is a simple association algorithm (from PFPileUp)
445  // extended to take fake vertices into account
446  double dzmin = 1.0e100;
447  unsigned index = 0;
448  for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
449  if (isAcceptableVtx(iv))
450  {
451  const double dz = std::abs(ztrack - iv->z());
452  if (dz < dzmin)
453  {
454  dzmin = dz;
455  iVertex = index;
456  foundVertex = true;
457  }
458  }
459 
461  {
462  const IV fakeEnd(fakeVertices->end());
463  index = 0;
464  for (IV iv=fakeVertices->begin(); iv!=fakeEnd; ++iv, ++index)
465  if (isAcceptableVtx(iv))
466  {
467  const double dz = std::abs(ztrack - iv->z());
468  if (dz < dzmin)
469  {
470  dzmin = dz;
471  iVertex = index;
472  *fromFakeSet = true;
473  foundVertex = true;
474  }
475  }
476  }
477  }
478 
479  if (foundVertex)
480  return reco::VertexRef(*fromFakeSet ? fakeVertices : vertices, iVertex);
481  }
482 
483  return reco::VertexRef();
484 }
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
#define abs(x)
Definition: mlp_lapack.h:159
float float float z
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:349
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:563
std::vector< LinkConnSpec >::const_iterator IT
bool first
Definition: L1TdeRCT.cc:94
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
std::vector< std::pair< double, unsigned > > zAssoc
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
T w() const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
bool isAcceptableVtx(reco::VertexCollection::const_iterator iv) const
bool FFTJetPFPileupCleaner::isAcceptableVtx ( reco::VertexCollection::const_iterator  iv) const
private

Definition at line 168 of file FFTJetPFPileupCleaner.cc.

References abs, vertexNdofCut, and vertexZmaxCut.

Referenced by findSomeVertexWFakes().

170 {
171  return !iv->isFake() &&
172  static_cast<double>(iv->ndof()) > vertexNdofCut &&
173  std::abs(iv->z()) < vertexZmaxCut;
174 }
#define abs(x)
Definition: mlp_lapack.h:159
bool FFTJetPFPileupCleaner::isRemovable ( reco::PFCandidate::ParticleType  ptype) const
private

Definition at line 257 of file FFTJetPFPileupCleaner.cc.

References removalMask, and edm::shift.

Referenced by produce().

259 {
260  const unsigned shift = static_cast<unsigned>(ptype);
261  assert(shift < 32U);
262  return removalMask & (1U << shift);
263 }
static unsigned int const shift
FFTJetPFPileupCleaner& FFTJetPFPileupCleaner::operator= ( const FFTJetPFPileupCleaner )
private
void FFTJetPFPileupCleaner::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
protectedvirtual

Implements edm::EDProducer.

Definition at line 178 of file FFTJetPFPileupCleaner.cc.

References eta(), etaMax, FakePrimaryVertices, findSomeVertexWFakes(), edm::Event::getByLabel(), i, edm::Ref< C, T, F >::isNull(), isRemovable(), edm::HandleBase::isValid(), edm::Ref< C, T, F >::key(), reco::tau::pfCandidates(), PFCandidates, edm::Event::put(), removeMainVertex, removeUnassociated, reverseRemovalDecision, reco::PFCandidate::setSourceCandidatePtr(), useFakePrimaryVertex, and Vertices.

180 {
181  // get PFCandidates
182  std::auto_ptr<reco::PFCandidateCollection>
183  pOutput(new reco::PFCandidateCollection);
184 
186  iEvent.getByLabel(PFCandidates, pfCandidates);
187 
188  // get vertices
190  iEvent.getByLabel(Vertices, vertices);
191 
194  {
195  iEvent.getByLabel(FakePrimaryVertices, fakeVertices);
196  if (!fakeVertices.isValid())
197  throw cms::Exception("FFTJetBadConfig")
198  << "ERROR in FFTJetPFPileupCleaner:"
199  " could not find fake vertices"
200  << std::endl;
201  }
202 
203  const unsigned ncand = pfCandidates->size();
204  for (unsigned i=0; i<ncand; ++i)
205  {
206  reco::PFCandidatePtr candptr(pfCandidates, i);
207  bool remove = false;
208 
209  if (isRemovable(candptr->particleId()))
210  {
211  bool fromFakeSet = false;
212  reco::VertexRef vertexref(findSomeVertexWFakes(vertices, fakeVertices,
213  *candptr, &fromFakeSet));
214  if (vertexref.isNull())
215  {
216  // Track is not associated with any vertex
217  // in any of the vertex sets
218  remove = removeUnassociated;
219  }
220  else if (vertexref.key() == 0 &&
221  (!useFakePrimaryVertex || fromFakeSet))
222  {
223  // Track is associated with the main primary vertex
224  // However, if we are using fake vertices, this only
225  // matters if the vertex comes from the fake set. If
226  // it comes from the real set, remove the track anyway
227  // because in the combined set the associated vertex
228  // would not be the main primary vertex.
229  remove = removeMainVertex;
230  }
231  else
232  remove = true;
233  }
234 
235  // Check the eta range
236  if (!remove)
237  {
238  const double eta = candptr->p4().Eta();
239  remove = eta < etaMin || eta > etaMax;
240  }
241 
242  // Should we remember this candidate?
244  remove = !remove;
245  if (!remove)
246  {
247  const reco::PFCandidate& cand = (*pfCandidates)[i];
248  pOutput->push_back(cand);
249  pOutput->back().setSourceCandidatePtr(candptr);
250  }
251  }
252 
253  iEvent.put(pOutput);
254 }
int i
Definition: DBlmapReader.cc:9
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:115
T eta() const
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
bool isRemovable(reco::PFCandidate::ParticleType ptype) const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:35
reco::VertexRef findSomeVertexWFakes(const edm::Handle< reco::VertexCollection > &vertices, const edm::Handle< reco::VertexCollection > &fakeVertices, const reco::PFCandidate &pfcand, bool *fromFakeSet) const
void FFTJetPFPileupCleaner::setRemovalBit ( reco::PFCandidate::ParticleType  ptype,
bool  onOff 
)
private

Definition at line 266 of file FFTJetPFPileupCleaner.cc.

References removalMask, and edm::shift.

Referenced by buildRemovalMask().

268 {
269  const unsigned shift = static_cast<unsigned>(ptype);
270  assert(shift < 32U);
271  const unsigned mask = (1U << shift);
272  if (value)
273  removalMask |= mask;
274  else
275  removalMask &= ~mask;
276 }
static unsigned int const shift

Member Data Documentation

bool FFTJetPFPileupCleaner::checkClosestZVertex
private

Definition at line 83 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().

double FFTJetPFPileupCleaner::etaMax
private

Definition at line 119 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

double FFTJetPFPileupCleaner::etaMin
private

Definition at line 118 of file FFTJetPFPileupCleaner.cc.

edm::InputTag FFTJetPFPileupCleaner::FakePrimaryVertices
private

Definition at line 75 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::keepIfPVneighbor
private

Definition at line 90 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().

edm::InputTag FFTJetPFPileupCleaner::PFCandidates
private

Definition at line 73 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

unsigned FFTJetPFPileupCleaner::removalMask
private

Definition at line 115 of file FFTJetPFPileupCleaner.cc.

Referenced by isRemovable(), and setRemovalBit().

bool FFTJetPFPileupCleaner::remove_e
private

Definition at line 107 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_egamma_HF
private

Definition at line 112 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_gamma
private

Definition at line 109 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h
private

Definition at line 106 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h0
private

Definition at line 110 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h_HF
private

Definition at line 111 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_mu
private

Definition at line 108 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_X
private

Definition at line 105 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::removeMainVertex
private

Definition at line 94 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::removeUnassociated
private

Definition at line 98 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::reverseRemovalDecision
private

Definition at line 101 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::useFakePrimaryVertex
private

Definition at line 79 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes(), and produce().

double FFTJetPFPileupCleaner::vertexNdofCut
private

Definition at line 122 of file FFTJetPFPileupCleaner.cc.

Referenced by isAcceptableVtx().

double FFTJetPFPileupCleaner::vertexZmaxCut
private

Definition at line 125 of file FFTJetPFPileupCleaner.cc.

Referenced by isAcceptableVtx().

edm::InputTag FFTJetPFPileupCleaner::Vertices
private

Definition at line 74 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

std::vector<std::pair<double, unsigned> > FFTJetPFPileupCleaner::zAssoc
mutableprivate

Definition at line 128 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().