test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FFTJetPFPileupCleaner.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetProducers
4 // Class: FFTJetPFPileupCleaner
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Thu Jul 14 17:50:33 CDT 2011
16 //
17 //
18 #include <cmath>
19 #include <climits>
20 #include <utility>
21 #include <algorithm>
22 
23 // framework include files
30 
33 
36 
39 
40 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
41 
42 //
43 // class declaration
44 //
46 {
47 public:
50 
51 protected:
52  // methods
53  void beginJob() override;
54  void produce(edm::Event&, const edm::EventSetup&) override;
55  void endJob() override;
56 
57 private:
61 
63  void setRemovalBit(reco::PFCandidate::ParticleType ptype, bool onOff);
64  void buildRemovalMask();
65  bool isAcceptableVtx(reco::VertexCollection::const_iterator iv) const;
66 
69  const edm::Handle<reco::VertexCollection>& fakeVertices,
70  const reco::PFCandidate& pfcand, bool* fromFakeSet) const;
71 
75 
79 
80  // The following, if true, will switch to an algorithm
81  // which takes a fake primary vertex into account
83 
84  // The following, if true, will cause association of a candidate
85  // with some vertex no matter what
87 
88  // The following switch will check if the primary vertex
89  // is a neighbor (in Z) of a track and will keep the
90  // track if this is so, even if it is not directly associated
91  // with the primary vertex. This switch is meaningful only if
92  // "checkClosestZVertex" is true.
94 
95  // The following, if true, will cause removal of candidates
96  // associated with the main vertex
98 
99  // The following will tell us to remove candidates not associated
100  // with any vertex
102 
103  // Do we want to reverse the decision?
105 
106  // Flags for removing things. See PFCandidate header
107  // for particle types.
108  bool remove_X; // undefined
109  bool remove_h; // charged hadron
110  bool remove_e; // electron
111  bool remove_mu; // muon
112  bool remove_gamma; // photon
113  bool remove_h0; // neutral hadron
114  bool remove_h_HF; // HF tower identified as a hadron
115  bool remove_egamma_HF; // HF tower identified as an EM particle
116 
117  // Mask for removing things
118  unsigned removalMask;
119 
120  // Min and max eta for keeping things
121  double etaMin;
122  double etaMax;
123 
124  // Cut for the vertex Ndof
126 
127  // Cut for the vertex Z
129 
130  // Cut for the vertex rho
131  double vertexRhoCut;
132 
133  // Vector for associating tracks with Z positions of the vertices
134  mutable std::vector<std::pair<double, unsigned> > zAssoc;
135 };
136 
137 //
138 // constructors and destructor
139 //
143  init_param(edm::InputTag, FakePrimaryVertices),
144  init_param(bool, useFakePrimaryVertex),
145  init_param(bool, checkClosestZVertex),
146  init_param(bool, keepIfPVneighbor),
147  init_param(bool, removeMainVertex),
148  init_param(bool, removeUnassociated),
149  init_param(bool, reverseRemovalDecision),
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 ),
157  init_param(bool, remove_egamma_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 }
174 
175 
177 {
178 }
179 
180 
182  reco::VertexCollection::const_iterator iv) const
183 {
184  return !iv->isFake() &&
185  static_cast<double>(iv->ndof()) > vertexNdofCut &&
186  std::abs(iv->z()) < vertexZmaxCut &&
187  iv->position().rho() < vertexRhoCut;
188 }
189 
190 
191 // ------------ method called to produce the data ------------
193  edm::Event& iEvent, const edm::EventSetup& iSetup)
194 {
195  // get PFCandidates
196  std::auto_ptr<reco::PFCandidateCollection>
197  pOutput(new reco::PFCandidateCollection);
198 
200  iEvent.getByToken(PFCandidatesToken, pfCandidates);
201 
202  // get vertices
204  iEvent.getByToken(VerticesToken, vertices);
205 
208  {
209  iEvent.getByToken(FakePrimaryVerticesToken, fakeVertices);
210  if (!fakeVertices.isValid())
211  throw cms::Exception("FFTJetBadConfig")
212  << "ERROR in FFTJetPFPileupCleaner:"
213  " could not find fake vertices"
214  << std::endl;
215  }
216 
217  const unsigned ncand = pfCandidates->size();
218  for (unsigned i=0; i<ncand; ++i)
219  {
220  reco::PFCandidatePtr candptr(pfCandidates, i);
221  bool remove = false;
222 
223  if (isRemovable(candptr->particleId()))
224  {
225  bool fromFakeSet = false;
226  reco::VertexRef vertexref(findSomeVertexWFakes(vertices, fakeVertices,
227  *candptr, &fromFakeSet));
228  if (vertexref.isNull())
229  {
230  // Track is not associated with any vertex
231  // in any of the vertex sets
232  remove = removeUnassociated;
233  }
234  else if (vertexref.key() == 0 &&
235  (!useFakePrimaryVertex || fromFakeSet))
236  {
237  // Track is associated with the main primary vertex
238  // However, if we are using fake vertices, this only
239  // matters if the vertex comes from the fake set. If
240  // it comes from the real set, remove the track anyway
241  // because in the combined set the associated vertex
242  // would not be the main primary vertex.
243  remove = removeMainVertex;
244  }
245  else
246  remove = true;
247  }
248 
249  // Check the eta range
250  if (!remove)
251  {
252  const double eta = candptr->p4().Eta();
253  remove = eta < etaMin || eta > etaMax;
254  }
255 
256  // Should we remember this candidate?
258  remove = !remove;
259  if (!remove)
260  {
261  const reco::PFCandidate& cand = (*pfCandidates)[i];
262  pOutput->push_back(cand);
263  pOutput->back().setSourceCandidatePtr(candptr);
264  }
265  }
266 
267  iEvent.put(pOutput);
268 }
269 
270 
272  const reco::PFCandidate::ParticleType ptype) const
273 {
274  const unsigned shift = static_cast<unsigned>(ptype);
275  assert(shift < 32U);
276  return removalMask & (1U << shift);
277 }
278 
279 
281  const reco::PFCandidate::ParticleType ptype, const bool value)
282 {
283  const unsigned shift = static_cast<unsigned>(ptype);
284  assert(shift < 32U);
285  const unsigned mask = (1U << shift);
286  if (value)
287  removalMask |= mask;
288  else
289  removalMask &= ~mask;
290 }
291 
292 
293 // The following essentially follows the code in PFPileUp.cc,
294 // with added cut on ndof, vertex Z position, and iteration
295 // over fakes
298  const edm::Handle<reco::VertexCollection>& fakeVertices,
299  const reco::PFCandidate& pfcand, bool* fromFakeSet) const
300 {
301  typedef reco::VertexCollection::const_iterator IV;
303 
304  *fromFakeSet = false;
305  reco::TrackBaseRef trackBaseRef(pfcand.trackRef());
306 
307  size_t iVertex = 0;
308  unsigned nFoundVertex = 0;
309  const IV vertend(vertices->end());
310 
311  {
312  unsigned index = 0;
313  double bestweight = 0.0;
314  for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
315  if (isAcceptableVtx(iv))
316  {
317  const reco::Vertex& vtx = *iv;
318 
319  // loop on tracks in vertices
320  IT trackend(vtx.tracks_end());
321  for (IT iTrack=vtx.tracks_begin(); iTrack!=trackend; ++iTrack)
322  {
323  const reco::TrackBaseRef& baseRef = *iTrack;
324 
325  // one of the tracks in the vertex is the same as
326  // the track considered in the function
327  if (baseRef == trackBaseRef)
328  {
329  // select the vertex for which the track has the highest weight
330  const double w = vtx.trackWeight(baseRef);
331  if (w > bestweight)
332  {
333  bestweight=w;
334  iVertex=index;
335  nFoundVertex++;
336  }
337  }
338  }
339  }
340  }
341 
342  if (nFoundVertex > 0)
343  {
344  if (nFoundVertex != 1)
345  edm::LogWarning("TrackOnTwoVertex")
346  << "a track is shared by at least two vertices. "
347  << "Used to be an assert";
348 
349  // Check if we can re-associate this track with one
350  // of the fake vertices
352  {
353  const double ztrack = pfcand.vertex().z();
354  double dzmin = std::abs(ztrack - ((*vertices)[iVertex]).z());
355 
356  const IV fakeEnd(fakeVertices->end());
357  unsigned index = 0;
358  for (IV iv=fakeVertices->begin(); iv!=fakeEnd; ++iv, ++index)
359  if (isAcceptableVtx(iv))
360  {
361  const double dz = std::abs(ztrack - iv->z());
362  if (dz < dzmin)
363  {
364  dzmin = dz;
365  iVertex = index;
366  *fromFakeSet = true;
367  }
368  }
369  }
370 
371  return reco::VertexRef(*fromFakeSet ? fakeVertices : vertices, iVertex);
372  }
373 
374  // optional: as a secondary solution, associate the closest vertex in z
375  if (checkClosestZVertex)
376  {
377  const double ztrack = pfcand.vertex().z();
378  bool foundVertex = false;
379 
380  if (keepIfPVneighbor)
381  {
382  // Sort all vertices according to their Z coordinate
383  zAssoc.clear();
384  unsigned index = 0;
385  for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
386  if (isAcceptableVtx(iv))
387  zAssoc.push_back(std::pair<double,unsigned>(iv->z(), index));
388  const unsigned numRealVertices = index;
389 
390  // Mix the fake vertex collection into zAssoc.
391  // Note that we do not reset "index" before doing this.
393  {
394  const IV fakeEnd(fakeVertices->end());
395  for (IV iv=fakeVertices->begin(); iv!=fakeEnd; ++iv, ++index)
396  if (isAcceptableVtx(iv))
397  zAssoc.push_back(std::pair<double,unsigned>(iv->z(), index));
398  }
399 
400  // Check where the track z position fits into this sequence
401  if (!zAssoc.empty())
402  {
403  std::sort(zAssoc.begin(), zAssoc.end());
404  std::pair<double,unsigned> tPair(ztrack, UINT_MAX);
405  const unsigned iAbove = std::upper_bound(
406  zAssoc.begin(), zAssoc.end(), tPair) - zAssoc.begin();
407 
408  // Check whether one of the vertices with indices
409  // iAbove or (iAbove - 1) is a primary vertex.
410  // If so, return it. Otherwise return the one
411  // with closest distance to the track.
412  unsigned ich[2] = {0U, 0U};
413  unsigned nch = 1;
414  if (iAbove)
415  {
416  ich[0] = iAbove - 1U;
417  ich[1] = iAbove;
418  if (iAbove < zAssoc.size())
419  nch = 2;
420  }
421 
422  double dzmin = 1.0e100;
423  unsigned bestVertexNum = UINT_MAX;
424  for (unsigned icheck=0; icheck<nch; ++icheck)
425  {
426  const unsigned zAssocIndex = ich[icheck];
427  const unsigned vertexNum = zAssoc[zAssocIndex].second;
428 
429  if (vertexNum == numRealVertices ||
430  (!useFakePrimaryVertex && vertexNum == 0U))
431  {
432  bestVertexNum = vertexNum;
433  break;
434  }
435 
436  const double dz = std::abs(ztrack - zAssoc[zAssocIndex].first);
437  if (dz < dzmin)
438  {
439  dzmin = dz;
440  bestVertexNum = vertexNum;
441  }
442  }
443 
444  foundVertex = bestVertexNum < UINT_MAX;
445  if (foundVertex)
446  {
447  iVertex = bestVertexNum;
448  if (iVertex >= numRealVertices)
449  {
450  *fromFakeSet = true;
451  iVertex -= numRealVertices;
452  }
453  }
454  }
455  }
456  else
457  {
458  // This is a simple association algorithm (from PFPileUp)
459  // extended to take fake vertices into account
460  double dzmin = 1.0e100;
461  unsigned index = 0;
462  for (IV iv=vertices->begin(); iv!=vertend; ++iv, ++index)
463  if (isAcceptableVtx(iv))
464  {
465  const double dz = std::abs(ztrack - iv->z());
466  if (dz < dzmin)
467  {
468  dzmin = dz;
469  iVertex = index;
470  foundVertex = true;
471  }
472  }
473 
475  {
476  const IV fakeEnd(fakeVertices->end());
477  index = 0;
478  for (IV iv=fakeVertices->begin(); iv!=fakeEnd; ++iv, ++index)
479  if (isAcceptableVtx(iv))
480  {
481  const double dz = std::abs(ztrack - iv->z());
482  if (dz < dzmin)
483  {
484  dzmin = dz;
485  iVertex = index;
486  *fromFakeSet = true;
487  foundVertex = true;
488  }
489  }
490  }
491  }
492 
493  if (foundVertex)
494  return reco::VertexRef(*fromFakeSet ? fakeVertices : vertices, iVertex);
495  }
496 
497  return reco::VertexRef();
498 }
499 
500 
502 {
511 }
512 
513 
514 // ------------ method called once each job just before starting event loop
516 {
517 }
518 
519 
520 // ------------ method called once each job just after ending the event loop
522 {
523 }
524 
525 
526 //define this as a plug-in
int i
Definition: DBlmapReader.cc:9
void setRemovalBit(reco::PFCandidate::ParticleType ptype, bool onOff)
ParticleType
particle types
Definition: PFCandidate.h:44
const double w
Definition: UKUtility.cc:23
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
assert(m_qm.get())
edm::EDGetTokenT< reco::VertexCollection > FakePrimaryVerticesToken
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:125
key_type key() const
Accessor for product key.
Definition: Ref.h:264
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
bool isRemovable(reco::PFCandidate::ParticleType ptype) const
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:647
edm::EDGetTokenT< reco::VertexCollection > VerticesToken
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:73
std::vector< LinkConnSpec >::const_iterator IT
bool isValid() const
Definition: HandleBase.h:75
void produce(edm::Event &, const edm::EventSetup &) override
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
bool isNull() const
Checks for null.
Definition: Ref.h:249
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
std::vector< std::pair< double, unsigned > > zAssoc
edm::EDGetTokenT< reco::PFCandidateCollection > PFCandidatesToken
#define init_param(type, varname)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
FFTJetPFPileupCleaner & operator=(const FFTJetPFPileupCleaner &)
reco::VertexRef findSomeVertexWFakes(const edm::Handle< reco::VertexCollection > &vertices, const edm::Handle< reco::VertexCollection > &fakeVertices, const reco::PFCandidate &pfcand, bool *fromFakeSet) const
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
static unsigned int const shift
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
bool isAcceptableVtx(reco::VertexCollection::const_iterator iv) const
tuple PFCandidates