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 ()
 
ModuleDescription const & moduleDescription () const
 
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 ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Protected Member Functions

void beginJob () override
 
void endJob () override
 
void produce (edm::Event &, const edm::EventSetup &) override
 
- 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
 
- 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 45 of file FFTJetPFPileupCleaner.cc.

Constructor & Destructor Documentation

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

Definition at line 133 of file FFTJetPFPileupCleaner.cc.

References buildRemovalMask().

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

Definition at line 162 of file FFTJetPFPileupCleaner.cc.

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

Member Function Documentation

void FFTJetPFPileupCleaner::beginJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 500 of file FFTJetPFPileupCleaner.cc.

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

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

487 {
496 }
void setRemovalBit(reco::PFCandidate::ParticleType ptype, bool onOff)
void FFTJetPFPileupCleaner::endJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 506 of file FFTJetPFPileupCleaner.cc.

507 {
508 }
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 281 of file FFTJetPFPileupCleaner.cc.

References funct::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().

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

Definition at line 167 of file FFTJetPFPileupCleaner.cc.

References funct::abs(), vertexNdofCut, and vertexZmaxCut.

Referenced by findSomeVertexWFakes().

169 {
170  return !iv->isFake() &&
171  static_cast<double>(iv->ndof()) > vertexNdofCut &&
172  std::abs(iv->z()) < vertexZmaxCut;
173 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool FFTJetPFPileupCleaner::isRemovable ( reco::PFCandidate::ParticleType  ptype) const
private

Definition at line 256 of file FFTJetPFPileupCleaner.cc.

References removalMask, and edm::shift.

Referenced by produce().

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

Implements edm::EDProducer.

Definition at line 177 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.

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

References removalMask, and edm::shift.

Referenced by buildRemovalMask().

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

Member Data Documentation

bool FFTJetPFPileupCleaner::checkClosestZVertex
private

Definition at line 82 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().

double FFTJetPFPileupCleaner::etaMax
private

Definition at line 118 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

double FFTJetPFPileupCleaner::etaMin
private

Definition at line 117 of file FFTJetPFPileupCleaner.cc.

edm::InputTag FFTJetPFPileupCleaner::FakePrimaryVertices
private

Definition at line 74 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::keepIfPVneighbor
private

Definition at line 89 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().

edm::InputTag FFTJetPFPileupCleaner::PFCandidates
private

Definition at line 72 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

unsigned FFTJetPFPileupCleaner::removalMask
private

Definition at line 114 of file FFTJetPFPileupCleaner.cc.

Referenced by isRemovable(), and setRemovalBit().

bool FFTJetPFPileupCleaner::remove_e
private

Definition at line 106 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_egamma_HF
private

Definition at line 111 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_gamma
private

Definition at line 108 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h
private

Definition at line 105 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h0
private

Definition at line 109 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h_HF
private

Definition at line 110 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_mu
private

Definition at line 107 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_X
private

Definition at line 104 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::removeMainVertex
private

Definition at line 93 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::removeUnassociated
private

Definition at line 97 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::reverseRemovalDecision
private

Definition at line 100 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::useFakePrimaryVertex
private

Definition at line 78 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes(), and produce().

double FFTJetPFPileupCleaner::vertexNdofCut
private

Definition at line 121 of file FFTJetPFPileupCleaner.cc.

Referenced by isAcceptableVtx().

double FFTJetPFPileupCleaner::vertexZmaxCut
private

Definition at line 124 of file FFTJetPFPileupCleaner.cc.

Referenced by isAcceptableVtx().

edm::InputTag FFTJetPFPileupCleaner::Vertices
private

Definition at line 73 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

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

Definition at line 127 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().