CMS 3D CMS Logo

PrimaryVertexAssignment.cc
Go to the documentation of this file.
8 
9 std::pair<int, PrimaryVertexAssignment::Quality> PrimaryVertexAssignment::chargedHadronVertex(
11  const reco::TrackRef& trackRef,
12  const reco::Track* track,
13  float time,
14  float timeReso, // <0 if timing not available for this object
16  const TransientTrackBuilder& builder) const {
17  typedef reco::VertexCollection::const_iterator IV;
19 
20  int iVertex = -1;
21  size_t index = 0;
22  float bestweight = 0;
23  for (auto const& vtx : vertices) {
24  float w = vtx.trackWeight(trackRef);
25  if (w > bestweight) {
26  bestweight = w;
27  iVertex = index;
28  }
29  index++;
30  }
31 
32  bool useTime = useTiming_;
33  if (edm::isNotFinite(time) || timeReso < 1e-6) {
34  useTime = false;
35  time = 0.;
36  timeReso = -1.;
37  }
38 
39  if (preferHighRanked_) {
40  for (IV iv = vertices.begin(); iv != vertices.end(); ++iv) {
41  int ivtx = iv - vertices.begin();
42  if (iVertex == ivtx)
43  return std::pair<int, PrimaryVertexAssignment::Quality>(ivtx, PrimaryVertexAssignment::UsedInFit);
44 
45  double dz = std::abs(track->dz(iv->position()));
46  double dt = std::abs(time - iv->t());
47 
48  bool useTimeVtx = useTime && iv->tError() > 0.;
49 
50  if ((dz < maxDzForPrimaryAssignment_ or dz / track->dzError() < maxDzSigForPrimaryAssignment_) and
51  (!useTimeVtx or dt / timeReso < maxDtSigForPrimaryAssignment_)) {
52  return std::pair<int, PrimaryVertexAssignment::Quality>(ivtx, PrimaryVertexAssignment::PrimaryDz);
53  }
54  }
55  }
56 
57  if (iVertex >= 0)
58  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
94  track->dzError() < maxDzErrorForPrimaryAssignment_) and
95  (!useTime or dtmin / timeReso < maxDtSigForPrimaryAssignment_)) {
96  iVertex = vtxIdMinSignif;
97  }
98 
99  if (iVertex >= 0)
100  return std::pair<int, PrimaryVertexAssignment::Quality>(iVertex, PrimaryVertexAssignment::PrimaryDz);
101 
102  // if track not assigned yet, it could be a b-decay secondary , use jet axis dist criterion
103  // first find the closest jet within maxJetDeltaR_
104  int jetIdx = -1;
105  double minDeltaR = maxJetDeltaR_;
106  for (edm::View<reco::Candidate>::const_iterator ij = jets.begin(); ij != jets.end(); ++ij) {
107  if (ij->pt() < minJetPt_)
108  continue; // skip jets below the jet Pt threshold
109 
110  double deltaR = reco::deltaR(*ij, *track);
111  if (deltaR < minDeltaR and track->dzError() < maxDzErrorForPrimaryAssignment_) {
112  minDeltaR = deltaR;
113  jetIdx = std::distance(jets.begin(), ij);
114  }
115  }
116  // if jet found
117  if (jetIdx != -1) {
118  reco::TransientTrack transientTrack = builder.build(*track);
119  GlobalVector direction(jets.at(jetIdx).px(), jets.at(jetIdx).py(), jets.at(jetIdx).pz());
120  // find the vertex with the smallest distanceToJetAxis that is still within maxDistaneToJetAxis_
121  int vtxIdx = -1;
122  double minDistanceToJetAxis = maxDistanceToJetAxis_;
123  for (IV iv = vertices.begin(); iv != vertices.end(); ++iv) {
124  // only check for vertices that are close enough in Z and for tracks that have not too high dXY
125  if (std::abs(track->dz(iv->position())) > maxDzForJetAxisAssigment_ ||
126  std::abs(track->dxy(iv->position())) > maxDxyForJetAxisAssigment_)
127  continue;
128 
129  double distanceToJetAxis = IPTools::jetTrackDistance(transientTrack, direction, *iv).second.value();
130  if (distanceToJetAxis < minDistanceToJetAxis) {
131  minDistanceToJetAxis = distanceToJetAxis;
132  vtxIdx = std::distance(vertices.begin(), iv);
133  }
134  }
135  if (vtxIdx >= 0) {
136  iVertex = vtxIdx;
137  }
138  }
139  if (iVertex >= 0)
140  return std::pair<int, PrimaryVertexAssignment::Quality>(iVertex, PrimaryVertexAssignment::BTrack);
141 
142  // if the track is not compatible with other PVs but is compatible with the BeamSpot, we may simply have not reco'ed the PV!
143  // we still point it to the closest in Z, but flag it as possible orphan-primary
144  if (!vertices.empty() && std::abs(track->dxy(vertices[0].position())) < maxDxyForNotReconstructedPrimary_ &&
145  std::abs(track->dxy(vertices[0].position()) / track->dxyError()) < maxDxySigForNotReconstructedPrimary_)
146  return std::pair<int, PrimaryVertexAssignment::Quality>(vtxIdMinSignif,
148 
149  //FIXME: here we could better handle V0s and NucInt
150 
151  // all other tracks could be non-B secondaries and we just attach them with closest Z
152  if (vtxIdMinSignif >= 0)
153  return std::pair<int, PrimaryVertexAssignment::Quality>(vtxIdMinSignif, PrimaryVertexAssignment::OtherDz);
154  //If for some reason even the dz failed (when?) we consider the track not assigned
155  return std::pair<int, PrimaryVertexAssignment::Quality>(-1, PrimaryVertexAssignment::Unassigned);
156 }
Vector3DBase
Definition: Vector3DBase.h:8
reco::Vertex::trackRef_iterator
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38
PrimaryVertexAssignment::preferHighRanked_
bool preferHighRanked_
Definition: PrimaryVertexAssignment.h:114
PrimaryVertexAssignment::OtherDz
Definition: PrimaryVertexAssignment.h:26
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
PFCandidate.h
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
PrimaryVertexAssignment::useTiming_
bool useTiming_
Definition: PrimaryVertexAssignment.h:113
PrimaryVertexAssignment::PrimaryDz
Definition: PrimaryVertexAssignment.h:22
PrimaryVertexAssignment::maxDtSigForPrimaryAssignment_
double maxDtSigForPrimaryAssignment_
Definition: PrimaryVertexAssignment.h:105
HLT_2018_cff.distance
distance
Definition: HLT_2018_cff.py:6417
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
PrimaryVertexAssignment::maxDzForPrimaryAssignment_
double maxDzForPrimaryAssignment_
Definition: PrimaryVertexAssignment.h:103
PrimaryVertexAssignment::BTrack
Definition: PrimaryVertexAssignment.h:24
PrimaryVertexAssignment::maxDxySigForNotReconstructedPrimary_
double maxDxySigForNotReconstructedPrimary_
Definition: PrimaryVertexAssignment.h:111
edm::Ref< TrackCollection >
deltaR.h
dt
float dt
Definition: AMPTWrapper.h:136
PrimaryVertexAssignment::maxDistanceToJetAxis_
double maxDistanceToJetAxis_
Definition: PrimaryVertexAssignment.h:108
PrimaryVertexAssignment::maxDzForJetAxisAssigment_
double maxDzForJetAxisAssigment_
Definition: PrimaryVertexAssignment.h:109
PrimaryVertexAssignment::maxDzErrorForPrimaryAssignment_
double maxDzErrorForPrimaryAssignment_
Definition: PrimaryVertexAssignment.h:104
w
const double w
Definition: UKUtility.cc:23
HLT_2018_cff.minDeltaR
minDeltaR
Definition: HLT_2018_cff.py:24818
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
PrimaryVertexAssignment::Unassigned
Definition: PrimaryVertexAssignment.h:28
reco::Track
Definition: Track.h:27
PrimaryVertexAssignment::maxJetDeltaR_
double maxJetDeltaR_
Definition: PrimaryVertexAssignment.h:106
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
badGlobalMuonTaggersAOD_cff.vtx
vtx
Definition: badGlobalMuonTaggersAOD_cff.py:5
Vertex.h
IPTools::jetTrackDistance
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:206
edm::View
Definition: CaloClusterFwd.h:14
PrimaryVertexAssignment::maxDzSigForPrimaryAssignment_
double maxDzSigForPrimaryAssignment_
Definition: PrimaryVertexAssignment.h:102
TransientTrackBuilder.h
PrimaryVertexAssignment::chargedHadronVertex
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
Definition: PrimaryVertexAssignment.cc:9
electrons_cff.jetIdx
jetIdx
Definition: electrons_cff.py:356
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
PrimaryVertexAssignment::maxDxyForNotReconstructedPrimary_
double maxDxyForNotReconstructedPrimary_
Definition: PrimaryVertexAssignment.h:112
PrimaryVertexAssignment::maxDxyForJetAxisAssigment_
double maxDxyForJetAxisAssigment_
Definition: PrimaryVertexAssignment.h:110
PrimaryVertexAssignment::minJetPt_
double minJetPt_
Definition: PrimaryVertexAssignment.h:107
TransientTrackBuilder
Definition: TransientTrackBuilder.h:16
PrimaryVertexAssignment::NotReconstructedPrimary
Definition: PrimaryVertexAssignment.h:27
IT
std::vector< LinkConnSpec >::const_iterator IT
Definition: TriggerBoardSpec.cc:5
PrimaryVertexAssignment.h
IPTools.h
reco::TransientTrack
Definition: TransientTrack.h:19
isFinite.h
PVValHelper::dz
Definition: PVValidationHelpers.h:50
or
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
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
TransientTrackBuilder::build
reco::TransientTrack build(const reco::Track *p) const
Definition: TransientTrackBuilder.cc:20
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PrimaryVertexAssignment::UsedInFit
Definition: PrimaryVertexAssignment.h:21
ntuplemaker.time
time
Definition: ntuplemaker.py:310
pwdgSkimBPark_cfi.vertices
vertices
Definition: pwdgSkimBPark_cfi.py:7
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37