CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
TrackVertexArbitration< VTX > Class Template Reference

#include <TrackVertexArbitration.h>

Public Member Functions

 TrackVertexArbitration (const edm::ParameterSet &params)
 
std::vector< VTX > trackVertexArbitrator (edm::Handle< reco::BeamSpot > &beamSpot, const reco::Vertex &pv, std::vector< reco::TransientTrack > &selectedTracks, std::vector< VTX > &secondaryVertices)
 

Private Member Functions

bool trackFilterArbitrator (const reco::TransientTrack &track) const
 

Private Attributes

edm::InputTag beamSpotCollection
 
double distCut
 
double dLenFraction
 
double dRCut
 
double fitterRatio
 
double fitterSigmacut
 
double fitterTini
 
double maxTimeSignificance
 
edm::InputTag primaryVertexCollection
 
edm::InputTag secondaryVertexCollection
 
double sigCut
 
double sign
 
edm::InputTag trackCollection
 
int trackMinLayers
 
int trackMinPixels
 
double trackMinPt
 

Detailed Description

template<class VTX>
class TrackVertexArbitration< VTX >

Definition at line 54 of file TrackVertexArbitration.h.

Constructor & Destructor Documentation

template<class VTX >
TrackVertexArbitration< VTX >::TrackVertexArbitration ( const edm::ParameterSet params)

Definition at line 86 of file TrackVertexArbitration.h.

References TrackVertexArbitration< VTX >::dRCut, and TrackVertexArbitration< VTX >::sign.

87  : primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
88  secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
89  trackCollection(params.getParameter<edm::InputTag>("tracks")),
90  beamSpotCollection(params.getParameter<edm::InputTag>("beamSpot")),
91  dRCut(params.getParameter<double>("dRCut")),
92  distCut(params.getParameter<double>("distCut")),
93  sigCut(params.getParameter<double>("sigCut")),
94  dLenFraction(params.getParameter<double>("dLenFraction")),
95  fitterSigmacut(params.getParameter<double>("fitterSigmacut")),
96  fitterTini(params.getParameter<double>("fitterTini")),
97  fitterRatio(params.getParameter<double>("fitterRatio")),
98  trackMinLayers(params.getParameter<int32_t>("trackMinLayers")),
99  trackMinPt(params.getParameter<double>("trackMinPt")),
100  trackMinPixels(params.getParameter<int32_t>("trackMinPixels")),
101  maxTimeSignificance(params.getParameter<double>("maxTimeSignificance")) {
102  sign = 1.0;
103  if (dRCut < 0) {
104  sign = -1.0;
105  }
106  dRCut *= dRCut;
107 }
edm::InputTag secondaryVertexCollection
T getParameter(std::string const &) const
Definition: ParameterSet.h:303

Member Function Documentation

template<class VTX >
bool TrackVertexArbitration< VTX >::trackFilterArbitrator ( const reco::TransientTrack track) const
private

Definition at line 109 of file TrackVertexArbitration.h.

References reco::TrackBase::hitPattern(), reco::TransientTrack::isValid(), reco::HitPattern::numberOfValidPixelHits(), reco::TrackBase::pt(), reco::TransientTrack::track(), reco::HitPattern::trackerLayersWithMeasurement(), HLT_FULL_cff::trackMinLayers, HLT_FULL_cff::trackMinPixels, and HLT_FULL_cff::trackMinPt.

109  {
110  if (!track.isValid())
111  return false;
113  return false;
114  if (track.track().pt() < trackMinPt)
115  return false;
117  return false;
118 
119  return true;
120 }
bool isValid() const
Make the ReferenceCountingProxy method to check validity public.
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:521
double pt() const
track transverse momentum
Definition: TrackBase.h:637
const Track & track() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
int numberOfValidPixelHits() const
Definition: HitPattern.h:825
template<class VTX >
std::vector< VTX > TrackVertexArbitration< VTX >::trackVertexArbitrator ( edm::Handle< reco::BeamSpot > &  beamSpot,
const reco::Vertex pv,
std::vector< reco::TransientTrack > &  selectedTracks,
std::vector< VTX > &  secondaryVertices 
)

Definition at line 141 of file TrackVertexArbitration.h.

References funct::abs(), IPTools::absoluteImpactParameter3D(), RecoVertex::convertError(), RecoVertex::convertPos(), gather_cfg::cout, svhelper::cov33(), reco::Vertex::covariance(), HLT_FULL_cff::deltaR, reco::deltaR2(), HLT_FULL_cff::distCut, HLT_FULL_cff::dLenFraction, HLT_FULL_cff::dRCut, reco::TransientTrack::dtErrorExt(), reco::TransientTrack::field(), HLT_FULL_cff::fitterRatio, HLT_FULL_cff::fitterSigmacut, HLT_FULL_cff::fitterTini, reco::TransientTrack::impactPointState(), edm::isFinite(), TransientVertex::isValid(), HLT_FULL_cff::maxTimeSignificance, reco::Vertex::position(), funct::pow(), dt_dqm_sourceclient_common_cff::reco, TrackCollections2monitor_cff::selectedTracks, reco::TransientTrack::setBeamSpot(), HLT_FULL_cff::sigCut, jetcorrextractor::sign(), Measurement1D::significance(), mathSSE::sqrt(), reco::TransientTrack::timeExt(), reco::TransientTrack::track(), trackWeight(), groupFilesInBlocks::tt, svhelper::updateVertexTime(), Measurement1D::value(), AdaptiveVertexFitter::vertex(), and w.

144  {
145  using namespace reco;
146 
147  //std::cout << "PV: " << pv.position() << std::endl;
148  VertexDistance3D dist;
154 
155  std::vector<VTX> recoVertices;
156  VertexDistance3D vdist;
157  GlobalPoint ppv(pv.position().x(), pv.position().y(), pv.position().z());
158 
159  std::unordered_map<unsigned int, Measurement1D> cachedIP;
160  for (typename std::vector<VTX>::const_iterator sv = secondaryVertices.begin(); sv != secondaryVertices.end(); ++sv) {
161  const double svTime(sv->t()), svTimeCov(svhelper::cov33(*sv));
162  const bool svHasTime = svTimeCov > 0.;
163 
164  GlobalPoint ssv(sv->position().x(), sv->position().y(), sv->position().z());
165  GlobalVector flightDir = ssv - ppv;
166  // std::cout << "Vertex : " << sv-secondaryVertices->begin() << " " << sv->position() << std::endl;
167  Measurement1D dlen =
168  vdist.distance(pv, VertexState(RecoVertex::convertPos(sv->position()), RecoVertex::convertError(sv->error())));
169  std::vector<reco::TransientTrack> selTracks;
170  for (unsigned int itrack = 0; itrack < selectedTracks.size(); itrack++) {
171  TransientTrack &tt = (selectedTracks)[itrack];
172  if (!trackFilterArbitrator(tt))
173  continue;
174  tt.setBeamSpot(*beamSpot);
175  float w = trackWeight(*sv, tt);
176  Measurement1D ipv;
177  if (cachedIP.count(itrack)) {
178  ipv = cachedIP[itrack];
179  } else {
180  std::pair<bool, Measurement1D> ipvp = IPTools::absoluteImpactParameter3D(tt, pv);
181  cachedIP[itrack] = ipvp.second;
182  ipv = ipvp.second;
183  }
184 
185  AnalyticalImpactPointExtrapolator extrapolator(tt.field());
187  extrapolator.extrapolate(tt.impactPointState(), RecoVertex::convertPos(sv->position()));
188  if (!tsos.isValid())
189  continue;
190  GlobalPoint refPoint = tsos.globalPosition();
191  GlobalError refPointErr = tsos.cartesianError().position();
192  Measurement1D isv =
193  dist.distance(VertexState(RecoVertex::convertPos(sv->position()), RecoVertex::convertError(sv->error())),
194  VertexState(refPoint, refPointErr));
195 
196  float dR = Geom::deltaR2(((sign > 0) ? flightDir : -flightDir),
197  tt.track()); //.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi());
198 
199  double timeSig = 0.;
200  if (svHasTime && edm::isFinite(tt.timeExt())) {
201  double tError = std::sqrt(std::pow(tt.dtErrorExt(), 2) + svTimeCov);
202  timeSig = std::abs(tt.timeExt() - svTime) / tError;
203  }
204 
205  if (w > 0 || (isv.significance() < sigCut && isv.value() < distCut && isv.value() < dlen.value() * dLenFraction &&
206  timeSig < maxTimeSignificance)) {
207  if ((isv.value() < ipv.value()) && isv.value() < distCut && isv.value() < dlen.value() * dLenFraction &&
208  dR < dRCut && timeSig < maxTimeSignificance) {
209 #ifdef VTXDEBUG
210  if (w > 0.5)
211  std::cout << " = ";
212  else
213  std::cout << " + ";
214 #endif
215  selTracks.push_back(tt);
216  } else {
217 #ifdef VTXDEBUG
218  if (w > 0.5 && isv.value() > ipv.value())
219  std::cout << " - ";
220  else
221  std::cout << " ";
222 #endif
223  //add also the tracks used in previous fitting that are still closer to Sv than Pv
224  if (w > 0.5 && isv.value() <= ipv.value() && dR < dRCut && timeSig < maxTimeSignificance) {
225  selTracks.push_back(tt);
226 #ifdef VTXDEBUG
227  std::cout << " = ";
228 #endif
229  }
230  if (w > 0.5 && isv.value() <= ipv.value() && dR >= dRCut) {
231 #ifdef VTXDEBUG
232  std::cout << " - ";
233 #endif
234  }
235  }
236 #ifdef VTXDEBUG
237 
238  std::cout << "t : " << itrack << " ref " << ref.key() << " w: " << w << " svip: " << isv.significance() << " "
239  << isv.value() << " pvip: " << ipv.significance() << " " << ipv.value() << " res "
240  << tt.track().residualX(0) << "," << tt.track().residualY(0)
241 // << " tpvip: " << itpv.second.significance() << " " << itpv.second.value() << " dr: " << dR << std::endl;
242 #endif
243 
244  } else {
245 #ifdef VTXDEBUG
246 
247  std::cout << " . t : " << itrack << " ref " << ref.key() << " w: " << w
248  << " svip: " << isv.second.significance() << " " << isv.second.value()
249  << " pvip: " << ipv.significance() << " " << ipv.value() << " dr: " << dR << std::endl;
250 #endif
251  }
252  }
253 
254  if (selTracks.size() >= 2) {
255  TransientVertex singleFitVertex;
256  singleFitVertex = theAdaptiveFitter.vertex(selTracks, ssv);
257 
258  if (singleFitVertex.isValid()) {
259  if (pv.covariance(3, 3) > 0.)
260  svhelper::updateVertexTime(singleFitVertex);
261  recoVertices.push_back(VTX(singleFitVertex));
262 
263 #ifdef VTXDEBUG
264  const VTX &extVertex = recoVertices.back();
265  GlobalVector vtxDir = GlobalVector(extVertex.p4().X(), extVertex.p4().Y(), extVertex.p4().Z());
266 
267  std::cout << " pointing : " << Geom::deltaR(extVertex.position() - pv.position(), vtxDir) << std::endl;
268 #endif
269  }
270  }
271  }
272  return recoVertices;
273 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
const double w
Definition: UKUtility.cc:23
void setBeamSpot(const reco::BeamSpot &beamSpot)
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
std::pair< bool, Measurement1D > absoluteImpactParameter3D(const reco::TransientTrack &transientTrack, const reco::Vertex &vertex)
Definition: IPTools.cc:38
bool trackFilterArbitrator(const reco::TransientTrack &track) const
double timeExt() const
const MagneticField * field() const
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:148
constexpr bool isFinite(T x)
const Point & position() const
position
Definition: Vertex.h:127
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
double significance() const
Definition: Measurement1D.h:29
const Track & track() const
double dtErrorExt() const
double value() const
Definition: Measurement1D.h:25
void updateVertexTime(TransientVertex &vtx)
Definition: SVTimeHelpers.h:11
double cov33(const reco::Vertex &sv)
tuple cout
Definition: gather_cfg.py:144
TrajectoryStateOnSurface impactPointState() const
tuple secondaryVertices
bool isValid() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
float trackWeight(const reco::Vertex &sv, const reco::TransientTrack &track)
Global3DVector GlobalVector
Definition: GlobalVector.h:10

Member Data Documentation

template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::beamSpotCollection
private

Definition at line 69 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::distCut
private

Definition at line 71 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::dLenFraction
private

Definition at line 73 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::dRCut
private
template<class VTX >
double TrackVertexArbitration< VTX >::fitterRatio
private

Definition at line 76 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::fitterSigmacut
private

Definition at line 74 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::fitterTini
private

Definition at line 75 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::maxTimeSignificance
private

Definition at line 80 of file TrackVertexArbitration.h.

template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::primaryVertexCollection
private

Definition at line 66 of file TrackVertexArbitration.h.

template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::secondaryVertexCollection
private

Definition at line 67 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::sigCut
private

Definition at line 72 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::sign
private
template<class VTX >
edm::InputTag TrackVertexArbitration< VTX >::trackCollection
private

Definition at line 68 of file TrackVertexArbitration.h.

template<class VTX >
int TrackVertexArbitration< VTX >::trackMinLayers
private

Definition at line 77 of file TrackVertexArbitration.h.

template<class VTX >
int TrackVertexArbitration< VTX >::trackMinPixels
private

Definition at line 79 of file TrackVertexArbitration.h.

template<class VTX >
double TrackVertexArbitration< VTX >::trackMinPt
private

Definition at line 78 of file TrackVertexArbitration.h.