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 
76  // The following, if true, will switch to an algorithm
77  // which takes a fake primary vertex into account
79 
80  // The following, if true, will cause association of a candidate
81  // with some vertex no matter what
83 
84  // The following switch will check if the primary vertex
85  // is a neighbor (in Z) of a track and will keep the
86  // track if this is so, even if it is not directly associated
87  // with the primary vertex. This switch is meaningful only if
88  // "checkClosestZVertex" is true.
90 
91  // The following, if true, will cause removal of candidates
92  // associated with the main vertex
94 
95  // The following will tell us to remove candidates not associated
96  // with any vertex
98 
99  // Do we want to reverse the decision?
101 
102  // Flags for removing things. See PFCandidate header
103  // for particle types.
104  bool remove_X; // undefined
105  bool remove_h; // charged hadron
106  bool remove_e; // electron
107  bool remove_mu; // muon
108  bool remove_gamma; // photon
109  bool remove_h0; // neutral hadron
110  bool remove_h_HF; // HF tower identified as a hadron
111  bool remove_egamma_HF; // HF tower identified as an EM particle
112 
113  // Mask for removing things
114  unsigned removalMask;
115 
116  // Min and max eta for keeping things
117  double etaMin;
118  double etaMax;
119 
120  // Cut for the vertex Ndof
122 
123  // Cut for the vertex Z
125 
126  // Vector for associating tracks with Z positions of the vertices
127  mutable std::vector<std::pair<double, unsigned> > zAssoc;
128 };
129 
130 //
131 // constructors and destructor
132 //
134  : init_param(edm::InputTag, PFCandidates),
135  init_param(edm::InputTag, Vertices),
136  init_param(edm::InputTag, FakePrimaryVertices),
137  init_param(bool, useFakePrimaryVertex),
138  init_param(bool, checkClosestZVertex),
139  init_param(bool, keepIfPVneighbor),
140  init_param(bool, removeMainVertex),
141  init_param(bool, removeUnassociated),
142  init_param(bool, reverseRemovalDecision),
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 ),
150  init_param(bool, remove_egamma_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 }
160 
161 
163 {
164 }
165 
166 
168  reco::VertexCollection::const_iterator iv) const
169 {
170  return !iv->isFake() &&
171  static_cast<double>(iv->ndof()) > vertexNdofCut &&
172  std::abs(iv->z()) < vertexZmaxCut;
173 }
174 
175 
176 // ------------ method called to produce the data ------------
178  edm::Event& iEvent, const edm::EventSetup& iSetup)
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 }
254 
255 
257  const reco::PFCandidate::ParticleType ptype) const
258 {
259  const unsigned shift = static_cast<unsigned>(ptype);
260  assert(shift < 32U);
261  return removalMask & (1U << shift);
262 }
263 
264 
266  const reco::PFCandidate::ParticleType ptype, const bool value)
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 }
276 
277 
278 // The following essentially follows the code in PFPileUp.cc,
279 // with added cut on ndof, vertex Z position, and iteration
280 // over fakes
282  const edm::Handle<reco::VertexCollection>& vertices,
283  const edm::Handle<reco::VertexCollection>& fakeVertices,
284  const reco::PFCandidate& pfcand, bool* fromFakeSet) const
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 }
484 
485 
487 {
496 }
497 
498 
499 // ------------ method called once each job just before starting event loop
501 {
502 }
503 
504 
505 // ------------ method called once each job just after ending the event loop
507 {
508 }
509 
510 
511 //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:43
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:123
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:429
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:116
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
bool isValid() const
Definition: HandleBase.h:76
void produce(edm::Event &, const edm::EventSetup &) override
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:390
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:38
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:37
static unsigned int const shift
T w() const
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
bool isAcceptableVtx(reco::VertexCollection::const_iterator iv) const
tuple PFCandidates