CMS 3D CMS Logo

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 () override
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducer () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Protected Member Functions

void beginJob () override
 
void endJob () override
 
void produce (edm::Event &, const edm::EventSetup &) override
 
- 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 ()=delete
 
 FFTJetPFPileupCleaner (const FFTJetPFPileupCleaner &)=delete
 
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 &)=delete
 
void setRemovalBit (reco::PFCandidate::ParticleType ptype, bool onOff)
 

Private Attributes

bool checkClosestZVertex
 
double etaMax
 
double etaMin
 
edm::InputTag FakePrimaryVertices
 
edm::EDGetTokenT< reco::VertexCollectionFakePrimaryVerticesToken
 
bool keepIfPVneighbor
 
edm::InputTag PFCandidates
 
edm::EDGetTokenT< reco::PFCandidateCollectionPFCandidatesToken
 
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 vertexRhoCut
 
double vertexZmaxCut
 
edm::InputTag Vertices
 
edm::EDGetTokenT< reco::VertexCollectionVerticesToken
 
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
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 

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 140 of file FFTJetPFPileupCleaner.cc.

References buildRemovalMask(), FakePrimaryVertices, FakePrimaryVerticesToken, PFCandidates, PFCandidatesToken, useFakePrimaryVertex, Vertices, and VerticesToken.

150  init_param(bool, remove_X ),
151  init_param(bool, remove_h ),
152  init_param(bool, remove_e ),
153  init_param(bool, remove_mu ),
154  init_param(bool, remove_gamma ),
155  init_param(bool, remove_h0 ),
156  init_param(bool, remove_h_HF ),
158  removalMask(0),
159  init_param(double, etaMin),
160  init_param(double, etaMax),
161  init_param(double, vertexNdofCut),
162  init_param(double, vertexZmaxCut),
163  init_param(double, vertexRhoCut)
164 {
166 
167  PFCandidatesToken = consumes<reco::PFCandidateCollection>(PFCandidates);
168  VerticesToken = consumes<reco::VertexCollection>(Vertices);
170  FakePrimaryVerticesToken = consumes<reco::VertexCollection>(FakePrimaryVertices);
171 
172  produces<reco::PFCandidateCollection>();
173 }
edm::EDGetTokenT< reco::VertexCollection > FakePrimaryVerticesToken
edm::EDGetTokenT< reco::VertexCollection > VerticesToken
edm::EDGetTokenT< reco::PFCandidateCollection > PFCandidatesToken
#define init_param(type, varname)
FFTJetPFPileupCleaner::~FFTJetPFPileupCleaner ( )
override

Definition at line 176 of file FFTJetPFPileupCleaner.cc.

177 {
178 }
FFTJetPFPileupCleaner::FFTJetPFPileupCleaner ( )
privatedelete
FFTJetPFPileupCleaner::FFTJetPFPileupCleaner ( const FFTJetPFPileupCleaner )
privatedelete

Member Function Documentation

void FFTJetPFPileupCleaner::beginJob ( void  )
overrideprotectedvirtual

Reimplemented from edm::EDProducer.

Definition at line 514 of file FFTJetPFPileupCleaner.cc.

515 {
516 }
void FFTJetPFPileupCleaner::buildRemovalMask ( )
private

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

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

Reimplemented from edm::EDProducer.

Definition at line 520 of file FFTJetPFPileupCleaner.cc.

References DEFINE_FWK_MODULE.

521 {
522 }
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 295 of file FFTJetPFPileupCleaner.cc.

References funct::abs(), checkClosestZVertex, PVValHelper::dz, plotBeamSpotDB::first, isAcceptableVtx(), keepIfPVneighbor, reco::PFCandidate::trackRef(), reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), mitigatedMETSequence_cff::U, useFakePrimaryVertex, reco::PFCandidate::vertex(), badGlobalMuonTaggersAOD_cff::vtx, w, z, and zAssoc.

Referenced by produce().

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

Definition at line 181 of file FFTJetPFPileupCleaner.cc.

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

Referenced by findSomeVertexWFakes().

183 {
184  return !iv->isFake() &&
185  static_cast<double>(iv->ndof()) > vertexNdofCut &&
186  std::abs(iv->z()) < vertexZmaxCut &&
187  iv->position().rho() < vertexRhoCut;
188 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool FFTJetPFPileupCleaner::isRemovable ( reco::PFCandidate::ParticleType  ptype) const
private

Definition at line 270 of file FFTJetPFPileupCleaner.cc.

References removalMask, edm::shift, and mitigatedMETSequence_cff::U.

Referenced by produce().

272 {
273  const unsigned shift = static_cast<unsigned>(ptype);
274  assert(shift < 32U);
275  return removalMask & (1U << shift);
276 }
static unsigned int const shift
FFTJetPFPileupCleaner& FFTJetPFPileupCleaner::operator= ( const FFTJetPFPileupCleaner )
privatedelete
void FFTJetPFPileupCleaner::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprotected

Definition at line 192 of file FFTJetPFPileupCleaner.cc.

References PVValHelper::eta, etaMax, FakePrimaryVerticesToken, findSomeVertexWFakes(), edm::Event::getByToken(), mps_fire::i, edm::Ref< C, T, F >::isNull(), isRemovable(), edm::HandleBase::isValid(), edm::Ref< C, T, F >::key(), eostools::move(), slimmedMuons_cfi::pfCandidates, PFCandidatesToken, edm::Event::put(), removeMainVertex, removeUnassociated, reverseRemovalDecision, useFakePrimaryVertex, electrons_cff::vertices, and VerticesToken.

194 {
195  // get PFCandidates
196  auto pOutput = std::make_unique<reco::PFCandidateCollection>();
197 
199  iEvent.getByToken(PFCandidatesToken, pfCandidates);
200 
201  // get vertices
203  iEvent.getByToken(VerticesToken, vertices);
204 
207  {
208  iEvent.getByToken(FakePrimaryVerticesToken, fakeVertices);
209  if (!fakeVertices.isValid())
210  throw cms::Exception("FFTJetBadConfig")
211  << "ERROR in FFTJetPFPileupCleaner:"
212  " could not find fake vertices"
213  << std::endl;
214  }
215 
216  const unsigned ncand = pfCandidates->size();
217  for (unsigned i=0; i<ncand; ++i)
218  {
219  reco::PFCandidatePtr candptr(pfCandidates, i);
220  bool remove = false;
221 
222  if (isRemovable(candptr->particleId()))
223  {
224  bool fromFakeSet = false;
225  reco::VertexRef vertexref(findSomeVertexWFakes(vertices, fakeVertices,
226  *candptr, &fromFakeSet));
227  if (vertexref.isNull())
228  {
229  // Track is not associated with any vertex
230  // in any of the vertex sets
231  remove = removeUnassociated;
232  }
233  else if (vertexref.key() == 0 &&
234  (!useFakePrimaryVertex || fromFakeSet))
235  {
236  // Track is associated with the main primary vertex
237  // However, if we are using fake vertices, this only
238  // matters if the vertex comes from the fake set. If
239  // it comes from the real set, remove the track anyway
240  // because in the combined set the associated vertex
241  // would not be the main primary vertex.
242  remove = removeMainVertex;
243  }
244  else
245  remove = true;
246  }
247 
248  // Check the eta range
249  if (!remove)
250  {
251  const double eta = candptr->p4().Eta();
252  remove = eta < etaMin || eta > etaMax;
253  }
254 
255  // Should we remember this candidate?
257  remove = !remove;
258  if (!remove)
259  {
260  const reco::PFCandidate& cand = (*pfCandidates)[i];
261  pOutput->push_back(cand);
262  pOutput->back().setSourceCandidatePtr(candptr);
263  }
264  }
265 
266  iEvent.put(std::move(pOutput));
267 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
edm::EDGetTokenT< reco::VertexCollection > FakePrimaryVerticesToken
bool isRemovable(reco::PFCandidate::ParticleType ptype) const
edm::EDGetTokenT< reco::VertexCollection > VerticesToken
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< reco::PFCandidateCollection > PFCandidatesToken
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
reco::VertexRef findSomeVertexWFakes(const edm::Handle< reco::VertexCollection > &vertices, const edm::Handle< reco::VertexCollection > &fakeVertices, const reco::PFCandidate &pfcand, bool *fromFakeSet) const
def move(src, dest)
Definition: eostools.py:510
void FFTJetPFPileupCleaner::setRemovalBit ( reco::PFCandidate::ParticleType  ptype,
bool  onOff 
)
private

Definition at line 279 of file FFTJetPFPileupCleaner.cc.

References RecoTauDiscriminantConfiguration::mask, removalMask, edm::shift, and mitigatedMETSequence_cff::U.

Referenced by buildRemovalMask().

281 {
282  const unsigned shift = static_cast<unsigned>(ptype);
283  assert(shift < 32U);
284  const unsigned mask = (1U << shift);
285  if (value)
286  removalMask |= mask;
287  else
288  removalMask &= ~mask;
289 }
Definition: value.py:1
static unsigned int const shift

Member Data Documentation

bool FFTJetPFPileupCleaner::checkClosestZVertex
private

Definition at line 86 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().

double FFTJetPFPileupCleaner::etaMax
private

Definition at line 122 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

double FFTJetPFPileupCleaner::etaMin
private

Definition at line 121 of file FFTJetPFPileupCleaner.cc.

edm::InputTag FFTJetPFPileupCleaner::FakePrimaryVertices
private

Definition at line 74 of file FFTJetPFPileupCleaner.cc.

Referenced by FFTJetPFPileupCleaner().

edm::EDGetTokenT<reco::VertexCollection> FFTJetPFPileupCleaner::FakePrimaryVerticesToken
private

Definition at line 78 of file FFTJetPFPileupCleaner.cc.

Referenced by FFTJetPFPileupCleaner(), and produce().

bool FFTJetPFPileupCleaner::keepIfPVneighbor
private

Definition at line 93 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().

edm::InputTag FFTJetPFPileupCleaner::PFCandidates
private

Definition at line 72 of file FFTJetPFPileupCleaner.cc.

Referenced by FFTJetPFPileupCleaner().

edm::EDGetTokenT<reco::PFCandidateCollection> FFTJetPFPileupCleaner::PFCandidatesToken
private

Definition at line 76 of file FFTJetPFPileupCleaner.cc.

Referenced by FFTJetPFPileupCleaner(), and produce().

unsigned FFTJetPFPileupCleaner::removalMask
private

Definition at line 118 of file FFTJetPFPileupCleaner.cc.

Referenced by isRemovable(), and setRemovalBit().

bool FFTJetPFPileupCleaner::remove_e
private

Definition at line 110 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_egamma_HF
private

Definition at line 115 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_gamma
private

Definition at line 112 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h
private

Definition at line 109 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h0
private

Definition at line 113 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_h_HF
private

Definition at line 114 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_mu
private

Definition at line 111 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::remove_X
private

Definition at line 108 of file FFTJetPFPileupCleaner.cc.

Referenced by buildRemovalMask().

bool FFTJetPFPileupCleaner::removeMainVertex
private

Definition at line 97 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::removeUnassociated
private

Definition at line 101 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::reverseRemovalDecision
private

Definition at line 104 of file FFTJetPFPileupCleaner.cc.

Referenced by produce().

bool FFTJetPFPileupCleaner::useFakePrimaryVertex
private

Definition at line 82 of file FFTJetPFPileupCleaner.cc.

Referenced by FFTJetPFPileupCleaner(), findSomeVertexWFakes(), and produce().

double FFTJetPFPileupCleaner::vertexNdofCut
private

Definition at line 125 of file FFTJetPFPileupCleaner.cc.

Referenced by isAcceptableVtx().

double FFTJetPFPileupCleaner::vertexRhoCut
private

Definition at line 131 of file FFTJetPFPileupCleaner.cc.

Referenced by isAcceptableVtx().

double FFTJetPFPileupCleaner::vertexZmaxCut
private

Definition at line 128 of file FFTJetPFPileupCleaner.cc.

Referenced by isAcceptableVtx().

edm::InputTag FFTJetPFPileupCleaner::Vertices
private

Definition at line 73 of file FFTJetPFPileupCleaner.cc.

Referenced by FFTJetPFPileupCleaner().

edm::EDGetTokenT<reco::VertexCollection> FFTJetPFPileupCleaner::VerticesToken
private

Definition at line 77 of file FFTJetPFPileupCleaner.cc.

Referenced by FFTJetPFPileupCleaner(), and produce().

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

Definition at line 134 of file FFTJetPFPileupCleaner.cc.

Referenced by findSomeVertexWFakes().