CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
PrimaryVertexAssignment Class Reference

#include <PrimaryVertexAssignment.h>

Public Types

enum  Quality {
  UsedInFit =7, PrimaryDz =6, PrimaryV0 =5, BTrack =4,
  Unused =3, OtherDz =2, NotReconstructedPrimary =1, Unassigned =0
}
 

Public Member Functions

std::pair< int, PrimaryVertexAssignment::QualitychargedHadronVertex (const reco::VertexCollection &vertices, const reco::TrackRef &trackRef, const reco::Track *track, float trackTime, float trackTimeResolution, const edm::View< reco::Candidate > &jets, const TransientTrackBuilder &builder) const
 
std::pair< int, PrimaryVertexAssignment::QualitychargedHadronVertex (const reco::VertexCollection &vertices, const reco::TrackRef &trackRef, float trackTime, float trackTimeResolution, const edm::View< reco::Candidate > &jets, const TransientTrackBuilder &builder) const
 
std::pair< int, PrimaryVertexAssignment::QualitychargedHadronVertex (const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand, const edm::View< reco::Candidate > &jets, const TransientTrackBuilder &builder) const
 
std::pair< int, PrimaryVertexAssignment::QualitychargedHadronVertex (const reco::VertexCollection &vertices, const reco::RecoChargedRefCandidate &chcand, const edm::ValueMap< float > *trackTimeTag, const edm::ValueMap< float > *trackTimeResoTag, const edm::View< reco::Candidate > &jets, const TransientTrackBuilder &builder) const
 
 PrimaryVertexAssignment (const edm::ParameterSet &iConfig)
 
 ~PrimaryVertexAssignment ()
 

Private Attributes

double maxDistanceToJetAxis_
 
double maxDtSigForPrimaryAssignment_
 
double maxDxyForJetAxisAssigment_
 
double maxDxyForNotReconstructedPrimary_
 
double maxDxySigForNotReconstructedPrimary_
 
double maxDzErrorForPrimaryAssignment_
 
double maxDzForJetAxisAssigment_
 
double maxDzForPrimaryAssignment_
 
double maxDzSigForPrimaryAssignment_
 
double maxJetDeltaR_
 
double minJetPt_
 
bool preferHighRanked_
 
bool useTiming_
 

Detailed Description

Definition at line 18 of file PrimaryVertexAssignment.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

PrimaryVertexAssignment::PrimaryVertexAssignment ( const edm::ParameterSet iConfig)
inline

Definition at line 22 of file PrimaryVertexAssignment.h.

22  :
23  maxDzSigForPrimaryAssignment_(iConfig.getParameter<double>("maxDzSigForPrimaryAssignment")),
24  maxDzForPrimaryAssignment_(iConfig.getParameter<double>("maxDzForPrimaryAssignment")),
25  maxDzErrorForPrimaryAssignment_(iConfig.getParameter<double>("maxDzErrorForPrimaryAssignment")),
26  maxDtSigForPrimaryAssignment_(iConfig.getParameter<double>("maxDtSigForPrimaryAssignment")),
27  maxJetDeltaR_(iConfig.getParameter<double>("maxJetDeltaR")),
28  minJetPt_(iConfig.getParameter<double>("minJetPt")),
29  maxDistanceToJetAxis_(iConfig.getParameter<double>("maxDistanceToJetAxis")),
30  maxDzForJetAxisAssigment_(iConfig.getParameter<double>("maxDzForJetAxisAssigment")),
31  maxDxyForJetAxisAssigment_(iConfig.getParameter<double>("maxDxyForJetAxisAssigment")),
32  maxDxySigForNotReconstructedPrimary_(iConfig.getParameter<double>("maxDxySigForNotReconstructedPrimary")),
33  maxDxyForNotReconstructedPrimary_(iConfig.getParameter<double>("maxDxyForNotReconstructedPrimary")),
34  useTiming_(iConfig.getParameter<bool>("useTiming")),
35  preferHighRanked_(iConfig.getParameter<bool>("preferHighRanked"))
36  {}
T getParameter(std::string const &) const
PrimaryVertexAssignment::~PrimaryVertexAssignment ( )
inline

Member Function Documentation

std::pair< int, PrimaryVertexAssignment::Quality > PrimaryVertexAssignment::chargedHadronVertex ( const reco::VertexCollection vertices,
const reco::TrackRef trackRef,
const reco::Track track,
float  trackTime,
float  trackTimeResolution,
const edm::View< reco::Candidate > &  jets,
const TransientTrackBuilder builder 
) const

Definition at line 11 of file PrimaryVertexAssignment.cc.

References funct::abs(), edm::View< T >::at(), edm::View< T >::begin(), BTrack, TransientTrackBuilder::build(), reco::deltaR(), boostedElectronIsolation_cff::deltaR, SoftLeptonByDistance_cfi::distance, dt, reco::TrackBase::dxy(), reco::TrackBase::dxyError(), PVValHelper::dz, reco::TrackBase::dz(), reco::TrackBase::dzError(), MillePedeFileConverter_cfg::e, edm::View< T >::end(), edm::isNotFinite(), electrons_cff::jetIdx, IPTools::jetTrackDistance(), SiStripPI::max, maxDistanceToJetAxis_, maxDtSigForPrimaryAssignment_, maxDxyForJetAxisAssigment_, maxDxyForNotReconstructedPrimary_, maxDxySigForNotReconstructedPrimary_, maxDzErrorForPrimaryAssignment_, maxDzForJetAxisAssigment_, maxDzForPrimaryAssignment_, maxDzSigForPrimaryAssignment_, maxJetDeltaR_, minJetPt_, NotReconstructedPrimary, or, OtherDz, preferHighRanked_, PrimaryDz, mathSSE::sqrt(), Unassigned, UsedInFit, useTiming_, badGlobalMuonTaggersAOD_cff::vtx, and w.

Referenced by chargedHadronVertex(), PrimaryVertexSorter< ParticlesCollection >::runAlgo(), and ~PrimaryVertexAssignment().

17  {
18 
19  typedef reco::VertexCollection::const_iterator IV;
21 
22  int iVertex = -1;
23  size_t index=0;
24  float bestweight=0;
25  for( auto const & vtx : vertices) {
26  float w = vtx.trackWeight(trackRef);
27  if (w > bestweight){
28  bestweight=w;
29  iVertex=index;
30  }
31  index++;
32  }
33 
34  bool useTime = useTiming_;
35  if (edm::isNotFinite(time) || timeReso<1e-6) {
36  useTime = false;
37  time = 0.;
38  timeReso = -1.;
39  }
40 
41  if (preferHighRanked_) {
42  for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv) {
43  int ivtx = iv - vertices.begin();
44  if (iVertex == ivtx) return std::pair<int,PrimaryVertexAssignment::Quality>(ivtx,PrimaryVertexAssignment::UsedInFit);
45 
46  double dz = std::abs(track->dz(iv->position()));
47  double dt = std::abs(time-iv->t());
48 
49  bool useTimeVtx = useTime && iv->tError()>0.;
50 
51  if ((dz < maxDzForPrimaryAssignment_ or dz/track->dzError() < maxDzSigForPrimaryAssignment_ ) and (!useTimeVtx or dt/timeReso < maxDtSigForPrimaryAssignment_)) {
52  return std::pair<int,PrimaryVertexAssignment::Quality>(ivtx,PrimaryVertexAssignment::PrimaryDz);
53  }
54  }
55  }
56 
57 
58  if(iVertex >= 0 ) return std::pair<int,PrimaryVertexAssignment::Quality>(iVertex,PrimaryVertexAssignment::UsedInFit);
59 
60  double distmin = std::numeric_limits<double>::max();
61  double dzmin = std::numeric_limits<double>::max();
62  double dtmin = std::numeric_limits<double>::max();
63  int vtxIdMinSignif = -1;
64  for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv) {
65  double dz = std::abs(track->dz(iv->position()));
66  double dt = std::abs(time-iv->t());
67 
68  double dzsig = dz/track->dzError();
69  double dist = dzsig*dzsig;
70 
71  bool useTimeVtx = useTime && iv->tError()>0.;
72  if (useTime && !useTimeVtx) {
73  dt = timeReso;
74  }
75 
76  if (useTime) {
77  double dtsig = dt/timeReso;
78 
79  dist += dtsig*dtsig;
80  }
81  if(dist<distmin) {
82  distmin = dist;
83  dzmin = dz;
84  dtmin = dt;
85  vtxIdMinSignif = iv-vertices.begin();
86  }
87  }
88 
89  // first use "closest in Z" with tight cuts (targetting primary particles)
90  const float add_cov = vtxIdMinSignif >= 0 ? vertices[vtxIdMinSignif].covariance(2,2) : 0.f;
91  const float dzE=sqrt(track->dzError()*track->dzError()+add_cov);
92  if(vtxIdMinSignif>=0 and
93  (dzmin < maxDzForPrimaryAssignment_ and dzmin/dzE < maxDzSigForPrimaryAssignment_ and track->dzError()<maxDzErrorForPrimaryAssignment_) and
94  (!useTime or dtmin/timeReso < maxDtSigForPrimaryAssignment_) )
95  {
96  iVertex=vtxIdMinSignif;
97  }
98 
99  if(iVertex >= 0 ) return std::pair<int,PrimaryVertexAssignment::Quality>(iVertex,PrimaryVertexAssignment::PrimaryDz);
100 
101 
102 
103  // if track not assigned yet, it could be a b-decay secondary , use jet axis dist criterion
104  // first find the closest jet within maxJetDeltaR_
105  int jetIdx = -1;
106  double minDeltaR = maxJetDeltaR_;
107  for(edm::View<reco::Candidate>::const_iterator ij=jets.begin(); ij!=jets.end(); ++ij)
108  {
109  if( ij->pt() < minJetPt_ ) continue; // skip jets below the jet Pt threshold
110 
111  double deltaR = reco::deltaR( *ij, *track );
112  if( deltaR < minDeltaR and track->dzError()<maxDzErrorForPrimaryAssignment_ )
113  {
114  minDeltaR = deltaR;
115  jetIdx = std::distance(jets.begin(), ij);
116  }
117  }
118  // if jet found
119  if( jetIdx!=-1 )
120  {
121  reco::TransientTrack transientTrack = builder.build(*track);
122  GlobalVector direction(jets.at(jetIdx).px(), jets.at(jetIdx).py(), jets.at(jetIdx).pz());
123  // find the vertex with the smallest distanceToJetAxis that is still within maxDistaneToJetAxis_
124  int vtxIdx = -1;
125  double minDistanceToJetAxis = maxDistanceToJetAxis_;
126  for(IV iv=vertices.begin(); iv!=vertices.end(); ++iv)
127  {
128  // only check for vertices that are close enough in Z and for tracks that have not too high dXY
129  if(std::abs(track->dz(iv->position())) > maxDzForJetAxisAssigment_ || std::abs(track->dxy(iv->position())) > maxDxyForJetAxisAssigment_)
130  continue;
131 
132  double distanceToJetAxis = IPTools::jetTrackDistance(transientTrack, direction, *iv).second.value();
133  if( distanceToJetAxis < minDistanceToJetAxis )
134  {
135  minDistanceToJetAxis = distanceToJetAxis;
136  vtxIdx = std::distance(vertices.begin(), iv);
137  }
138  }
139  if( vtxIdx>=0 )
140  {
141  iVertex=vtxIdx;
142  }
143  }
144  if(iVertex >= 0 )
145  return std::pair<int,PrimaryVertexAssignment::Quality>(iVertex,PrimaryVertexAssignment::BTrack);
146 
147  // if the track is not compatible with other PVs but is compatible with the BeamSpot, we may simply have not reco'ed the PV!
148  // we still point it to the closest in Z, but flag it as possible orphan-primary
149  if(!vertices.empty() && std::abs(track->dxy(vertices[0].position()))<maxDxyForNotReconstructedPrimary_ && std::abs(track->dxy(vertices[0].position())/track->dxyError())<maxDxySigForNotReconstructedPrimary_)
150  return std::pair<int,PrimaryVertexAssignment::Quality>(vtxIdMinSignif,PrimaryVertexAssignment::NotReconstructedPrimary);
151 
152  //FIXME: here we could better handle V0s and NucInt
153 
154  // all other tracks could be non-B secondaries and we just attach them with closest Z
155  if(vtxIdMinSignif>=0)
156  return std::pair<int,PrimaryVertexAssignment::Quality>(vtxIdMinSignif,PrimaryVertexAssignment::OtherDz);
157  //If for some reason even the dz failed (when?) we consider the track not assigned
158  return std::pair<int,PrimaryVertexAssignment::Quality>(-1,PrimaryVertexAssignment::Unassigned);
159 }
float dt
Definition: AMPTWrapper.h:126
const double w
Definition: UKUtility.cc:23
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
double dxyError() const
error on dxy
Definition: TrackBase.h:841
reco::TransientTrack build(const reco::Track *p) const
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:200
const_iterator begin() const
T sqrt(T t)
Definition: SSEVec.h:18
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< LinkConnSpec >::const_iterator IT
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:642
double dzError() const
error on dz
Definition: TrackBase.h:859
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
const_reference at(size_type pos) const
const_iterator end() const
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:624
std::pair<int,PrimaryVertexAssignment::Quality> PrimaryVertexAssignment::chargedHadronVertex ( const reco::VertexCollection vertices,
const reco::TrackRef trackRef,
float  trackTime,
float  trackTimeResolution,
const edm::View< reco::Candidate > &  jets,
const TransientTrackBuilder builder 
) const
inline

Definition at line 48 of file PrimaryVertexAssignment.h.

References chargedHadronVertex().

54  {
55  return chargedHadronVertex(vertices,trackRef,&(*trackRef),trackTime,trackTimeResolution,jets,builder);
56  }
std::pair< int, PrimaryVertexAssignment::Quality > chargedHadronVertex(const reco::VertexCollection &vertices, const reco::TrackRef &trackRef, const reco::Track *track, float trackTime, float trackTimeResolution, const edm::View< reco::Candidate > &jets, const TransientTrackBuilder &builder) const
std::pair<int,PrimaryVertexAssignment::Quality> PrimaryVertexAssignment::chargedHadronVertex ( const reco::VertexCollection vertices,
const reco::PFCandidate pfcand,
const edm::View< reco::Candidate > &  jets,
const TransientTrackBuilder builder 
) const
inline

Definition at line 58 of file PrimaryVertexAssignment.h.

References chargedHadronVertex(), reco::PFCandidate::gsfTrackRef(), edm::Ref< C, T, F >::isNull(), reco::PFCandidate::isTimeValid(), fwrapper::jets, ntuplemaker::time, reco::PFCandidate::time(), reco::PFCandidate::timeError(), muonGEMDigis_cfi::timeResolution, reco::PFCandidate::trackRef(), Unassigned, and useTiming_.

61  {
62  float time = 0, timeResolution = -1;
63  if (useTiming_ && pfcand.isTimeValid()) {
64  time = pfcand.time(); timeResolution = pfcand.timeError();
65  }
66  if(pfcand.gsfTrackRef().isNull())
67  {
68  if(pfcand.trackRef().isNull())
69  return std::pair<int,PrimaryVertexAssignment::Quality>(-1,PrimaryVertexAssignment::Unassigned);
70  else
72  }
73  return chargedHadronVertex(vertices,reco::TrackRef(),&(*pfcand.gsfTrackRef()),time,timeResolution,jets,builder);
74  }
float time() const
Definition: PFCandidate.h:421
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:442
bool isTimeValid() const
do we have a valid time information
Definition: PFCandidate.h:419
vector< PseudoJet > jets
std::pair< int, PrimaryVertexAssignment::Quality > chargedHadronVertex(const reco::VertexCollection &vertices, const reco::TrackRef &trackRef, const reco::Track *track, float trackTime, float trackTimeResolution, const edm::View< reco::Candidate > &jets, const TransientTrackBuilder &builder) const
bool isNull() const
Checks for null.
Definition: Ref.h:248
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:480
float timeError() const
Definition: PFCandidate.h:423
std::pair<int,PrimaryVertexAssignment::Quality> PrimaryVertexAssignment::chargedHadronVertex ( const reco::VertexCollection vertices,
const reco::RecoChargedRefCandidate chcand,
const edm::ValueMap< float > *  trackTimeTag,
const edm::ValueMap< float > *  trackTimeResoTag,
const edm::View< reco::Candidate > &  jets,
const TransientTrackBuilder builder 
) const
inline

Definition at line 75 of file PrimaryVertexAssignment.h.

References chargedHadronVertex(), edm::Ref< C, T, F >::isNull(), fwrapper::jets, ntuplemaker::time, muonGEMDigis_cfi::timeResolution, reco::RecoChargedRefCandidate::track(), Unassigned, and useTiming_.

80  {
81  float time = 0, timeResolution = -1;
82  if (useTiming_) {
83  time = (*trackTimeTag)[chcand.track()];
84  timeResolution = (*trackTimeResoTag)[chcand.track()];
85  }
86  if(chcand.track().isNull())
87  return std::pair<int,PrimaryVertexAssignment::Quality>(-1,PrimaryVertexAssignment::Unassigned);
88  return chargedHadronVertex(vertices,chcand.track(),time,timeResolution,jets,builder);
89  }
vector< PseudoJet > jets
std::pair< int, PrimaryVertexAssignment::Quality > chargedHadronVertex(const reco::VertexCollection &vertices, const reco::TrackRef &trackRef, const reco::Track *track, float trackTime, float trackTimeResolution, const edm::View< reco::Candidate > &jets, const TransientTrackBuilder &builder) const
bool isNull() const
Checks for null.
Definition: Ref.h:248

Member Data Documentation

double PrimaryVertexAssignment::maxDistanceToJetAxis_
private

Definition at line 99 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDtSigForPrimaryAssignment_
private

Definition at line 96 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDxyForJetAxisAssigment_
private

Definition at line 101 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDxyForNotReconstructedPrimary_
private

Definition at line 103 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDxySigForNotReconstructedPrimary_
private

Definition at line 102 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDzErrorForPrimaryAssignment_
private

Definition at line 95 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDzForJetAxisAssigment_
private

Definition at line 100 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDzForPrimaryAssignment_
private

Definition at line 94 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxDzSigForPrimaryAssignment_
private

Definition at line 93 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::maxJetDeltaR_
private

Definition at line 97 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

double PrimaryVertexAssignment::minJetPt_
private

Definition at line 98 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

bool PrimaryVertexAssignment::preferHighRanked_
private

Definition at line 105 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().

bool PrimaryVertexAssignment::useTiming_
private

Definition at line 104 of file PrimaryVertexAssignment.h.

Referenced by chargedHadronVertex().