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 // $Id: FFTJetPFPileupCleaner.cc,v 1.5 2012/06/22 07:23:21 igv Exp $
17 //
18 //
19 #include <cmath>
20 #include <climits>
21 #include <utility>
22 #include <algorithm>
23 
24 // framework include files
31 
34 
37 
40 
41 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
42 
43 //
44 // class declaration
45 //
47 {
48 public:
51 
52 protected:
53  // methods
54  void beginJob();
55  void produce(edm::Event&, const edm::EventSetup&);
56  void endJob();
57 
58 private:
62 
64  void setRemovalBit(reco::PFCandidate::ParticleType ptype, bool onOff);
65  void buildRemovalMask();
66  bool isAcceptableVtx(reco::VertexCollection::const_iterator iv) const;
67 
70  const edm::Handle<reco::VertexCollection>& fakeVertices,
71  const reco::PFCandidate& pfcand, bool* fromFakeSet) const;
72 
76 
77  // The following, if true, will switch to an algorithm
78  // which takes a fake primary vertex into account
80 
81  // The following, if true, will cause association of a candidate
82  // with some vertex no matter what
84 
85  // The following switch will check if the primary vertex
86  // is a neighbor (in Z) of a track and will keep the
87  // track if this is so, even if it is not directly associated
88  // with the primary vertex. This switch is meaningful only if
89  // "checkClosestZVertex" is true.
91 
92  // The following, if true, will cause removal of candidates
93  // associated with the main vertex
95 
96  // The following will tell us to remove candidates not associated
97  // with any vertex
99 
100  // Do we want to reverse the decision?
102 
103  // Flags for removing things. See PFCandidate header
104  // for particle types.
105  bool remove_X; // undefined
106  bool remove_h; // charged hadron
107  bool remove_e; // electron
108  bool remove_mu; // muon
109  bool remove_gamma; // photon
110  bool remove_h0; // neutral hadron
111  bool remove_h_HF; // HF tower identified as a hadron
112  bool remove_egamma_HF; // HF tower identified as an EM particle
113 
114  // Mask for removing things
115  unsigned removalMask;
116 
117  // Min and max eta for keeping things
118  double etaMin;
119  double etaMax;
120 
121  // Cut for the vertex Ndof
123 
124  // Cut for the vertex Z
126 
127  // Vector for associating tracks with Z positions of the vertices
128  mutable std::vector<std::pair<double, unsigned> > zAssoc;
129 };
130 
131 //
132 // constructors and destructor
133 //
135  : init_param(edm::InputTag, PFCandidates),
136  init_param(edm::InputTag, Vertices),
137  init_param(edm::InputTag, FakePrimaryVertices),
138  init_param(bool, useFakePrimaryVertex),
139  init_param(bool, checkClosestZVertex),
140  init_param(bool, keepIfPVneighbor),
141  init_param(bool, removeMainVertex),
142  init_param(bool, removeUnassociated),
143  init_param(bool, reverseRemovalDecision),
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 ),
151  init_param(bool, remove_egamma_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 }
161 
162 
164 {
165 }
166 
167 
169  reco::VertexCollection::const_iterator iv) const
170 {
171  return !iv->isFake() &&
172  static_cast<double>(iv->ndof()) > vertexNdofCut &&
173  std::abs(iv->z()) < vertexZmaxCut;
174 }
175 
176 
177 // ------------ method called to produce the data ------------
179  edm::Event& iEvent, const edm::EventSetup& iSetup)
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 }
255 
256 
258  const reco::PFCandidate::ParticleType ptype) const
259 {
260  const unsigned shift = static_cast<unsigned>(ptype);
261  assert(shift < 32U);
262  return removalMask & (1U << shift);
263 }
264 
265 
267  const reco::PFCandidate::ParticleType ptype, const bool value)
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 }
277 
278 
279 // The following essentially follows the code in PFPileUp.cc,
280 // with added cut on ndof, vertex Z position, and iteration
281 // over fakes
283  const edm::Handle<reco::VertexCollection>& vertices,
284  const edm::Handle<reco::VertexCollection>& fakeVertices,
285  const reco::PFCandidate& pfcand, bool* fromFakeSet) const
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 }
485 
486 
488 {
497 }
498 
499 
500 // ------------ method called once each job just before starting event loop
502 {
503 }
504 
505 
506 // ------------ method called once each job just after ending the event loop
508 {
509 }
510 
511 
512 //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:40
void produce(edm::Event &, const edm::EventSetup &)
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
#define abs(x)
Definition: mlp_lapack.h:159
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:115
T eta() const
float float float z
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:349
bool isRemovable(reco::PFCandidate::ParticleType ptype) const
int iEvent
Definition: GenABIO.cc:243
bool isNull() const
Checks for null.
Definition: Ref.h:247
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
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
bool isValid() const
Definition: HandleBase.h:76
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
std::vector< std::pair< double, unsigned > > zAssoc
key_type key() const
Accessor for product key.
Definition: Ref.h:266
#define init_param(type, varname)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:35
FFTJetPFPileupCleaner & operator=(const FFTJetPFPileupCleaner &)
std::vector< VertexType > Vertices
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:38
static unsigned int const shift
T w() const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
bool isAcceptableVtx(reco::VertexCollection::const_iterator iv) const
tuple PFCandidates