CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PFPileUpAlgo Class Reference

#include <PFPileUpAlgo.h>

Public Member Functions

int chargedHadronVertex (const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand) const
 
const reco::PFCandidateCollectiongetPFCandidatesFromPU () const
 
const reco::PFCandidateCollectiongetPFCandidatesFromVtx () const
 
 PFPileUpAlgo ()
 
 PFPileUpAlgo (bool checkClosestZVertex, bool verbose=false)
 
void process (const reco::PFCandidateCollection &pfCandidates, const reco::VertexCollection &vertices, const edm::Handle< reco::PFCandidateCollection > *handle=0)
 
void setCheckClosestZVertex (bool val)
 
void setVerbose (bool verbose)
 
 ~PFPileUpAlgo ()
 

Private Attributes

bool checkClosestZVertex_
 use the closest z vertex if a track is not in a vertex More...
 
reco::PFCandidateCollection pfCandidatesFromPU_
 
reco::PFCandidateCollection pfCandidatesFromVtx_
 
bool verbose_
 verbose ? More...
 

Detailed Description

Definition at line 12 of file PFPileUpAlgo.h.

Constructor & Destructor Documentation

PFPileUpAlgo::PFPileUpAlgo ( )
inline

Definition at line 14 of file PFPileUpAlgo.h.

14 :checkClosestZVertex_(true), verbose_(false) {;}
bool verbose_
verbose ?
Definition: PFPileUpAlgo.h:45
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition: PFPileUpAlgo.h:41
PFPileUpAlgo::PFPileUpAlgo ( bool  checkClosestZVertex,
bool  verbose = false 
)
inline

Definition at line 16 of file PFPileUpAlgo.h.

16  :
17  checkClosestZVertex_(checkClosestZVertex), verbose_(verbose) {;}
bool verbose_
verbose ?
Definition: PFPileUpAlgo.h:45
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition: PFPileUpAlgo.h:41
PFPileUpAlgo::~PFPileUpAlgo ( )
inline

Definition at line 19 of file PFPileUpAlgo.h.

19 {;}

Member Function Documentation

int PFPileUpAlgo::chargedHadronVertex ( const reco::VertexCollection vertices,
const reco::PFCandidate pfcand 
) const

Definition at line 44 of file PFPileUpAlgo.cc.

References checkClosestZVertex_, getHLTprescales::index, reco::PFCandidate::trackRef(), reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), reco::PFCandidate::vertex(), and w().

Referenced by process().

44  {
45 
46 
47  reco::TrackBaseRef trackBaseRef( pfcand.trackRef() );
48 
49  size_t iVertex = 0;
50  unsigned index=0;
51  unsigned nFoundVertex = 0;
52  typedef reco::VertexCollection::const_iterator IV;
54  float bestweight=0;
55  for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) {
56 
57  const reco::Vertex& vtx = *iv;
58 
59  // loop on tracks in vertices
60  for(IT iTrack=vtx.tracks_begin();
61  iTrack!=vtx.tracks_end(); ++iTrack) {
62 
63  const reco::TrackBaseRef& baseRef = *iTrack;
64 
65  // one of the tracks in the vertex is the same as
66  // the track considered in the function
67  if(baseRef == trackBaseRef ) {
68  float w = vtx.trackWeight(baseRef);
69  //select the vertex for which the track has the highest weight
70  if (w > bestweight){
71  bestweight=w;
72  iVertex=index;
73  nFoundVertex++;
74  }
75  }
76  }
77  }
78 
79  if (nFoundVertex>0){
80  if (nFoundVertex!=1)
81  edm::LogWarning("TrackOnTwoVertex")<<"a track is shared by at least two verteces. Used to be an assert";
82  return iVertex;
83  }
84  // no vertex found with this track.
85 
86  // optional: as a secondary solution, associate the closest vertex in z
87  if ( checkClosestZVertex_ ) {
88 
89  double dzmin = 10000;
90  double ztrack = pfcand.vertex().z();
91  bool foundVertex = false;
92  index = 0;
93  for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) {
94 
95  double dz = fabs(ztrack - iv->z());
96  if(dz<dzmin) {
97  dzmin = dz;
98  iVertex = index;
99  foundVertex = true;
100  }
101  }
102 
103  if( foundVertex )
104  return iVertex;
105 
106  }
107 
108 
109  return -1 ;
110 }
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:339
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
Definition: PFCandidate.cc:546
std::vector< LinkConnSpec >::const_iterator IT
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition: PFPileUpAlgo.h:41
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
T w() const
const reco::PFCandidateCollection& PFPileUpAlgo::getPFCandidatesFromPU ( ) const
inline

Definition at line 30 of file PFPileUpAlgo.h.

References pfCandidatesFromPU_.

30 {return pfCandidatesFromPU_;}
reco::PFCandidateCollection pfCandidatesFromPU_
Definition: PFPileUpAlgo.h:48
const reco::PFCandidateCollection& PFPileUpAlgo::getPFCandidatesFromVtx ( ) const
inline

Definition at line 32 of file PFPileUpAlgo.h.

References pfCandidatesFromVtx_.

32 {return pfCandidatesFromVtx_;}
reco::PFCandidateCollection pfCandidatesFromVtx_
Definition: PFPileUpAlgo.h:47
void PFPileUpAlgo::process ( const reco::PFCandidateCollection pfCandidates,
const reco::VertexCollection vertices,
const edm::Handle< reco::PFCandidateCollection > *  handle = 0 
)

Definition at line 5 of file PFPileUpAlgo.cc.

References chargedHadronVertex(), reco::PFCandidate::h, i, reco::PFCandidate::particleId(), pfCandidatesFromPU_, and pfCandidatesFromVtx_.

Referenced by Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::dumpPython(), ConfigBuilder.ConfigBuilder.PrintAllModules::leave(), Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::open(), Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::outputEventContent(), Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::setProcess(), and Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::setProperty().

7  {
8 
9  pfCandidatesFromVtx_.clear();
10  pfCandidatesFromPU_.clear();
11 
12  for( unsigned i=0; i<pfCandidates.size(); i++ ) {
13 
14  const reco::PFCandidate& cand = pfCandidates[i];
15 
16  int ivertex;
17 
18  switch( cand.particleId() ) {
20  ivertex = chargedHadronVertex( vertices, cand );
21  break;
22  default:
23  continue;
24  }
25 
26  // no associated vertex, or primary vertex
27  // not pile-up
28  if( ivertex == -1 ||
29  ivertex == 0 ) {
30  pfCandidatesFromVtx_.push_back(cand);
31  pfCandidatesFromVtx_.back().setSourceCandidatePtr( reco::PFCandidatePtr(*handle,i) );
32  } else {
33  // associated to a vertex
34  pfCandidatesFromPU_.push_back(cand);
35  if(handle) {
36  pfCandidatesFromPU_.back().setSourceCandidatePtr( reco::PFCandidatePtr(*handle,i) );
37  }
38  }
39  }
40 }
int i
Definition: DBlmapReader.cc:9
reco::PFCandidateCollection pfCandidatesFromVtx_
Definition: PFPileUpAlgo.h:47
reco::PFCandidateCollection pfCandidatesFromPU_
Definition: PFPileUpAlgo.h:48
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
int chargedHadronVertex(const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand) const
Definition: PFPileUpAlgo.cc:44
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
virtual ParticleType particleId() const
Definition: PFCandidate.h:324
void PFPileUpAlgo::setCheckClosestZVertex ( bool  val)
inline

Definition at line 28 of file PFPileUpAlgo.h.

References checkClosestZVertex_.

28 { checkClosestZVertex_ = val;}
bool checkClosestZVertex_
use the closest z vertex if a track is not in a vertex
Definition: PFPileUpAlgo.h:41
void PFPileUpAlgo::setVerbose ( bool  verbose)
inline

Member Data Documentation

bool PFPileUpAlgo::checkClosestZVertex_
private

use the closest z vertex if a track is not in a vertex

Definition at line 41 of file PFPileUpAlgo.h.

Referenced by chargedHadronVertex(), and setCheckClosestZVertex().

reco::PFCandidateCollection PFPileUpAlgo::pfCandidatesFromPU_
private

Definition at line 48 of file PFPileUpAlgo.h.

Referenced by getPFCandidatesFromPU(), and process().

reco::PFCandidateCollection PFPileUpAlgo::pfCandidatesFromVtx_
private

Definition at line 47 of file PFPileUpAlgo.h.

Referenced by getPFCandidatesFromVtx(), and process().

bool PFPileUpAlgo::verbose_
private

verbose ?

Definition at line 45 of file PFPileUpAlgo.h.

Referenced by setVerbose().