CMS 3D CMS Logo

GhostTrackVertexFinder.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <iterator>
3 #include <vector>
4 #include <map>
5 #include <set>
6 
7 #include <Math/SVector.h>
8 #include <Math/SMatrix.h>
9 #include <Math/MatrixFunctions.h>
10 
14 
16 
23 
37 
44 
46 
47 // #define DEBUG
48 
49 #ifdef DEBUG
51 #endif
52 
53 using namespace reco;
54 
55 namespace {
56  using namespace ROOT::Math;
57 
58  typedef SVector<double, 3> Vector3;
59 
60  typedef SMatrix<double, 3, 3, MatRepSym<double, 3> > Matrix3S;
61  typedef SMatrix<double, 2, 2, MatRepSym<double, 2> > Matrix2S;
62  typedef SMatrix<double, 2, 3> Matrix23;
63 
65  RefCountedVertexTrack;
67  RefCountedLinearizedTrackState;
68 
69  struct VtxTrackIs {
70  VtxTrackIs(const TransientTrack &track) : track(track) {}
71  bool operator()(const RefCountedVertexTrack &vtxTrack) const
72  { return vtxTrack->linearizedTrack()->track() == track; }
73 
74  const TransientTrack &track;
75  };
76 
77  inline Vector3 conv(const GlobalVector &vec)
78  {
79  Vector3 result;
80  result[0] = vec.x();
81  result[1] = vec.y();
82  result[2] = vec.z();
83  return result;
84  }
85 
86  inline double sqr(double arg) { return arg * arg; }
87 }
88 
91  const GhostTrack &ghostTrack, const BeamSpot &beamSpot,
92  bool hasBeamSpot, bool hasPrimaries) :
93  primaryVertex(primaryVertex), pred(ghostTrack.prediction()),
94  prior(ghostTrack.prior()), states(ghostTrack.states()),
95  beamSpot(beamSpot), hasBeamSpot(hasBeamSpot),
96  hasPrimaries(hasPrimaries), field(nullptr) {}
97 
101  std::vector<GhostTrackState> states;
107 };
108 
110  const MagneticField *field)
111 { return TransientTrackFromFTSFactory().build(pred.fts(field)); }
112 
113 static double vtxErrorLong(const GlobalError &error, const GlobalVector &dir)
114 {
115  return ROOT::Math::Similarity(conv(dir), error.matrix());
116 }
117 
118 static GlobalPoint vtxMean(const GlobalPoint &p1, const GlobalError &e1,
119  const GlobalPoint &p2, const GlobalError &e2)
120 {
121  GlobalVector diff = p2 - p1;
122 
123  double err1 = vtxErrorLong(e1, diff);
124  double err2 = vtxErrorLong(e2, diff);
125 
126  double weight = err1 / (err1 + err2);
127 
128  return p1 + weight * diff;
129 }
130 
131 static VertexState stateMean(const VertexState &v1, const VertexState &v2)
132 {
133  return VertexState(vtxMean(v1.position(), v1.error(),
134  v2.position(), v2.error()),
135  v1.error() + v2.error());
136 }
137 
138 static bool covarianceUpdate(Matrix3S &cov, const Vector3 &residual,
139  const Matrix3S &error, double &chi2,
140  double theta, double phi)
141 {
142  using namespace ROOT::Math;
143 
144  Matrix23 jacobian;
145  jacobian(0, 0) = std::cos(phi) * std::cos(theta);
146  jacobian(0, 1) = std::sin(phi) * std::cos(theta);
147  jacobian(0, 2) = -std::sin(theta);
148  jacobian(1, 0) = -std::sin(phi);
149  jacobian(1, 1) = std::cos(phi);
150 
151  Matrix2S measErr = Similarity(jacobian, error);
152  Matrix2S combErr = Similarity(jacobian, cov) + measErr;
153  if (!measErr.Invert() || !combErr.Invert())
154  return false;
155 
156  cov -= SimilarityT(jacobian * cov, combErr);
157  chi2 += Similarity(jacobian * residual, measErr);
158 
159  return true;
160 }
161 
163  const GhostTrackPrediction &pred,
164  const GhostTrackState &state)
165 {
166  LinearizedTrackStateFactory linTrackFactory;
167  VertexTrackFactory<5> vertexTrackFactory;
168 
169  if (!state.isValid())
170  return CachingVertex<5>();
171 
172  GlobalPoint pca1 = pred.position(state.lambda());
173  GlobalError err1 = pred.positionError(state.lambda());
174 
175  GlobalPoint pca2 = state.globalPosition();
176  GlobalError err2 = state.cartesianError();
177 
178  GlobalPoint point = vtxMean(pca1, err1, pca2, err2);
179 
180  const TransientTrack& recTrack = state.track();
181 
182  RefCountedLinearizedTrackState linState[2] = {
183  linTrackFactory.linearizedTrackState(point, ghostTrack),
184  linTrackFactory.linearizedTrackState(point, recTrack)
185  };
186  if ( !linState[0]->isValid() || !linState[1]->isValid() )
187  return CachingVertex<5>();
188 
189  Matrix3S cov = SMatrixIdentity();
190  cov *= 10000;
191  double chi2 = 0.;
192  if (!covarianceUpdate(cov, conv(pca1 - point), err1.matrix(), chi2,
193  linState[0]->predictedStateParameters()[1],
194  linState[0]->predictedStateParameters()[2]) ||
195  !covarianceUpdate(cov, conv(pca2 - point), err2.matrix(), chi2,
196  linState[1]->predictedStateParameters()[1],
197  linState[1]->predictedStateParameters()[2]))
198  return CachingVertex<5>();
199 
200  GlobalError error(cov);
201  VertexState vtxState(point, error);
202 
203  std::vector<RefCountedVertexTrack> linTracks(2);
204  linTracks[0] = vertexTrackFactory.vertexTrack(linState[0], vtxState);
205  linTracks[1] = vertexTrackFactory.vertexTrack(linState[1], vtxState);
206 
207  return CachingVertex<5>(point, error, linTracks, chi2);
208 }
209 
210 static RefCountedVertexTrack relinearizeTrack(
211  const RefCountedVertexTrack &track,
212  const VertexState &state,
213  const VertexTrackFactory<5> factory)
214 {
215  RefCountedLinearizedTrackState linTrack = track->linearizedTrack();
216  linTrack = linTrack->stateWithNewLinearizationPoint(state.position());
217  return factory.vertexTrack(linTrack, state);
218 }
219 
220 static std::vector<RefCountedVertexTrack> relinearizeTracks(
221  const std::vector<RefCountedVertexTrack> &tracks,
222  const VertexState &state)
223 {
224  VertexTrackFactory<5> vertexTrackFactory;
225 
226  std::vector<RefCountedVertexTrack> finalTracks;
227  finalTracks.reserve(tracks.size());
228 
229  for(std::vector<RefCountedVertexTrack>::const_iterator iter =
230  tracks.begin(); iter != tracks.end(); ++iter)
231  finalTracks.push_back(
232  relinearizeTrack(*iter, state, vertexTrackFactory));
233 
234  return finalTracks;
235 }
236 
238  const CachingVertex<5> &vtx2,
239  const FinderInfo &info,
240  double scale1, double scale2)
241 {
242  Vector3 diff = conv(vtx2.position() - vtx1.position());
243  Matrix3S cov = scale1 * vtx1.error().matrix() +
244  scale2 * vtx2.error().matrix();
245 
246  return sqr(ROOT::Math::Dot(diff, diff)) / ROOT::Math::Similarity(cov, diff);
247 }
248 
250  const RefCountedVertexTrack &vertexTrack)
251 {
252  using namespace ROOT::Math;
253 
254  TransientTrack track = vertexTrack->linearizedTrack()->track();
255  GlobalPoint point = vtx.position();
258  extrap.extrapolate(track.impactPointState(), point);
259 
260  if (!tsos.isValid())
261  return 1.0e6;
262 
263  GlobalPoint point1 = vtx.position();
264  GlobalPoint point2 = tsos.globalPosition();
265  Vector3 dir = conv(point2 - point1);
266  Matrix3S error = vtx.error().matrix() +
267  tsos.cartesianError().matrix().Sub<Matrix3S>(0, 0);
268  if (!error.Invert())
269  return 1.0e6;
270 
271  return ROOT::Math::Similarity(error, dir);
272 }
273 
275  maxFitChi2_(10.0),
276  mergeThreshold_(3.0),
277  primcut_(2.0),
278  seccut_(5.0),
279  fitType_(kAlwaysWithGhostTrack)
280 {
281 }
282 
284  double maxFitChi2, double mergeThreshold, double primcut,
285  double seccut, FitType fitType) :
286  maxFitChi2_(maxFitChi2),
287  mergeThreshold_(mergeThreshold),
288  primcut_(primcut),
289  seccut_(seccut),
290  fitType_(fitType)
291 {
292 }
293 
295 {
296 }
297 
299 {
300  if (!ghostTrackFitter_.get())
302 
303  return *ghostTrackFitter_;
304 }
305 
307 {
308  std::unique_ptr<VertexFitter<5> > *ptr =
309  primary ? &primVertexFitter_ : &secVertexFitter_;
310  if (!ptr->get())
311  ptr->reset(new AdaptiveVertexFitter(
312  GeometricAnnealing(primary
313  ? primcut_ : seccut_)));
314 
315  return **ptr;
316 }
317 
318 #ifdef DEBUG
319 static void debugTrack(const TransientTrack &track, const TransientVertex *vtx = 0)
320 {
321  const FreeTrajectoryState &fts = track.initialFreeState();
322  GlobalVector mom = fts.momentum();
323  std::cout << "\tpt = " << mom.perp() << ", eta = " << mom.eta() << ", phi = " << mom.phi() << ", charge = " << fts.charge();
324  if (vtx && vtx->hasTrackWeight())
325  std::cout << ", w = " << vtx->trackWeight(track) << std::endl;
326  else
327  std::cout << std::endl;
328 }
329 
330 static void debugTracks(const std::vector<TransientTrack> &tracks, const TransientVertex *vtx = 0)
331 {
332  TrackKinematics kin;
333  for(std::vector<TransientTrack>::const_iterator iter = tracks.begin();
334  iter != tracks.end(); ++iter) {
335  debugTrack(*iter, vtx);
336  try {
337  TrackBaseRef track = iter->trackBaseRef();
338  kin.add(*track);
339  } catch ( exception const & e ) {
340  // ignore, yeah this is debugging, so shut up ^^
341  }
342  }
343 
344  std::cout << "mass = " << kin.vectorSum().M() << std::endl;
345 }
346 
347 static void debugTracks(const std::vector<RefCountedVertexTrack> &tracks)
348 {
349  std::vector<TransientTrack> tracks_;
350  for(std::vector<RefCountedVertexTrack>::const_iterator iter = tracks.begin();
351  iter != tracks.end(); ++iter) {
352  std::cout << "+ " << (*iter)->linearizedTrack()->linearizationPoint() << std::endl;
353  tracks_.push_back((*iter)->linearizedTrack()->track());
354  }
355  debugTracks(tracks_);
356 }
357 
358 static void debugVertex(const TransientVertex &vtx, const GhostTrackPrediction &pred)
359 {
360  double origin = 0.;
361  std::cout << vtx.position() << "-> " << vtx.totalChiSquared() << " with " << vtx.degreesOfFreedom() << std::endl;
362 
363  double err = std::sqrt(vtxErrorLong(
364  vtx.positionError(),
365  pred.direction())
366  / pred.rho2());
367  std::cout << "* " << (pred.lambda(vtx.position()) * pred.rho() - origin)
368  << " +-" << err
369  << " => " << ((pred.lambda(vtx.position()) * pred.rho() - origin) / err)
370  << std::endl;
371 
372  std::vector<TransientTrack> tracks = vtx.originalTracks();
373  debugTracks(tracks, &vtx);
374 }
375 #endif
376 
377 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
378  const Vertex &primaryVertex,
379  const GlobalVector &direction,
380  double coneRadius,
381  const std::vector<TransientTrack> &tracks) const
382 {
383  return vertices(RecoVertex::convertPos(primaryVertex.position()),
384  RecoVertex::convertError(primaryVertex.error()),
385  direction, coneRadius, tracks);
386 }
387 
388 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
389  const Vertex &primaryVertex,
390  const GlobalVector &direction,
391  double coneRadius,
392  const reco::BeamSpot &beamSpot,
393  const std::vector<TransientTrack> &tracks) const
394 {
395  return vertices(RecoVertex::convertPos(primaryVertex.position()),
396  RecoVertex::convertError(primaryVertex.error()),
397  direction, coneRadius, beamSpot, tracks);
398 }
399 
400 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
401  const Vertex &primaryVertex,
402  const GlobalVector &direction,
403  double coneRadius,
404  const reco::BeamSpot &beamSpot,
405  const std::vector<TransientTrack> &primaries,
406  const std::vector<TransientTrack> &tracks) const
407 {
408  return vertices(RecoVertex::convertPos(primaryVertex.position()),
409  RecoVertex::convertError(primaryVertex.error()),
410  direction, coneRadius, beamSpot, primaries, tracks);
411 }
412 
413 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
414  const GlobalPoint &primaryPosition,
415  const GlobalError &primaryError,
416  const GlobalVector &direction,
417  double coneRadius,
418  const std::vector<TransientTrack> &tracks) const
419 {
420  GhostTrack ghostTrack =
421  ghostTrackFitter().fit(primaryPosition, primaryError,
422  direction, coneRadius, tracks);
423 
424  CachingVertex<5> primary(primaryPosition, primaryError,
425  std::vector<RefCountedVertexTrack>(), 0.);
426 
427  std::vector<TransientVertex> result = vertices(ghostTrack, primary);
428 
429 #ifdef DEBUG
430  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
431  iter != result.end(); ++iter)
432  debugVertex(*iter, ghostTrack.prediction());
433 #endif
434 
435  return result;
436 }
437 
438 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
439  const GlobalPoint &primaryPosition,
440  const GlobalError &primaryError,
441  const GlobalVector &direction,
442  double coneRadius,
443  const BeamSpot &beamSpot,
444  const std::vector<TransientTrack> &tracks) const
445 {
446  GhostTrack ghostTrack =
447  ghostTrackFitter().fit(primaryPosition, primaryError,
448  direction, coneRadius, tracks);
449 
450  CachingVertex<5> primary(primaryPosition, primaryError,
451  std::vector<RefCountedVertexTrack>(), 0.);
452 
453  std::vector<TransientVertex> result =
454  vertices(ghostTrack, primary, beamSpot, true);
455 
456 #ifdef DEBUG
457  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
458  iter != result.end(); ++iter)
459  debugVertex(*iter, ghostTrack.prediction());
460 #endif
461 
462  return result;
463 }
464 
465 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
466  const GlobalPoint &primaryPosition,
467  const GlobalError &primaryError,
468  const GlobalVector &direction,
469  double coneRadius,
470  const BeamSpot &beamSpot,
471  const std::vector<TransientTrack> &primaries,
472  const std::vector<TransientTrack> &tracks) const
473 {
474  GhostTrack ghostTrack =
475  ghostTrackFitter().fit(primaryPosition, primaryError,
476  direction, coneRadius, tracks);
477 
478  std::vector<RefCountedVertexTrack> primaryVertexTracks;
479  if (!primaries.empty()) {
480  LinearizedTrackStateFactory linTrackFactory;
481  VertexTrackFactory<5> vertexTrackFactory;
482 
483  VertexState state(primaryPosition, primaryError);
484 
485  for(std::vector<TransientTrack>::const_iterator iter =
486  primaries.begin(); iter != primaries.end(); ++iter) {
487 
488  RefCountedLinearizedTrackState linState =
489  linTrackFactory.linearizedTrackState(
490  primaryPosition, *iter);
491 
492  primaryVertexTracks.push_back(
493  vertexTrackFactory.vertexTrack(
494  linState, state));
495  }
496  }
497 
498  CachingVertex<5> primary(primaryPosition, primaryError,
499  primaryVertexTracks, 0.);
500 
501  std::vector<TransientVertex> result =
502  vertices(ghostTrack, primary, beamSpot, true, true);
503 
504 #ifdef DEBUG
505  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
506  iter != result.end(); ++iter)
507  debugVertex(*iter, ghostTrack.prediction());
508 #endif
509 
510  return result;
511 }
512 
513 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
514  const GlobalPoint &primaryPosition,
515  const GlobalError &primaryError,
516  const GhostTrack &ghostTrack) const
517 {
518  CachingVertex<5> primary(primaryPosition, primaryError,
519  std::vector<RefCountedVertexTrack>(), 0.);
520 
521  std::vector<TransientVertex> result = vertices(ghostTrack, primary);
522 
523 #ifdef DEBUG
524  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
525  iter != result.end(); ++iter)
526  debugVertex(*iter, ghostTrack.prediction());
527 #endif
528 
529  return result;
530 }
531 
532 
533 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
534  const GlobalPoint &primaryPosition,
535  const GlobalError &primaryError,
536  const BeamSpot &beamSpot,
537  const GhostTrack &ghostTrack) const
538 {
539  CachingVertex<5> primary(primaryPosition, primaryError,
540  std::vector<RefCountedVertexTrack>(), 0.);
541 
542  std::vector<TransientVertex> result =
543  vertices(ghostTrack, primary, beamSpot, true);
544 
545 #ifdef DEBUG
546  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
547  iter != result.end(); ++iter)
548  debugVertex(*iter, ghostTrack.prediction());
549 #endif
550 
551  return result;
552 }
553 
554 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
555  const GlobalPoint &primaryPosition,
556  const GlobalError &primaryError,
557  const BeamSpot &beamSpot,
558  const std::vector<TransientTrack> &primaries,
559  const GhostTrack &ghostTrack) const
560 {
561  std::vector<RefCountedVertexTrack> primaryVertexTracks;
562  if (!primaries.empty()) {
563  LinearizedTrackStateFactory linTrackFactory;
564  VertexTrackFactory<5> vertexTrackFactory;
565 
566  VertexState state(primaryPosition, primaryError);
567 
568  for(std::vector<TransientTrack>::const_iterator iter =
569  primaries.begin(); iter != primaries.end(); ++iter) {
570 
571  RefCountedLinearizedTrackState linState =
572  linTrackFactory.linearizedTrackState(
573  primaryPosition, *iter);
574 
575  primaryVertexTracks.push_back(
576  vertexTrackFactory.vertexTrack(
577  linState, state));
578  }
579  }
580 
581  CachingVertex<5> primary(primaryPosition, primaryError,
582  primaryVertexTracks, 0.);
583 
584  std::vector<TransientVertex> result =
585  vertices(ghostTrack, primary, beamSpot, true, true);
586 
587 #ifdef DEBUG
588  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
589  iter != result.end(); ++iter)
590  debugVertex(*iter, ghostTrack.prediction());
591 #endif
592 
593  return result;
594 }
595 
596 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
597  const Vertex &primaryVertex,
598  const GhostTrack &ghostTrack) const
599 {
600  return vertices(RecoVertex::convertPos(primaryVertex.position()),
601  RecoVertex::convertError(primaryVertex.error()),
602  ghostTrack);
603 }
604 
605 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
606  const Vertex &primaryVertex,
607  const BeamSpot &beamSpot,
608  const GhostTrack &ghostTrack) const
609 {
610  return vertices(RecoVertex::convertPos(primaryVertex.position()),
611  RecoVertex::convertError(primaryVertex.error()),
612  beamSpot, ghostTrack);
613 }
614 
615 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
616  const Vertex &primaryVertex,
617  const BeamSpot &beamSpot,
618  const std::vector<TransientTrack> &primaries,
619  const GhostTrack &ghostTrack) const
620 {
621  return vertices(RecoVertex::convertPos(primaryVertex.position()),
622  RecoVertex::convertError(primaryVertex.error()),
623  beamSpot, primaries, ghostTrack);
624 }
625 
627  const Vertex &primaryVertex, const Track &ghostTrack)
628 {
629  return GhostTrackPrediction(
630  RecoVertex::convertPos(primaryVertex.position()),
631  RecoVertex::convertError(primaryVertex.error()),
632  GlobalVector(ghostTrack.px(), ghostTrack.py(),
633  ghostTrack.pz()), 0.05);
634 }
635 
636 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
637  const Vertex &primaryVertex,
638  const Track &ghostTrack,
639  const std::vector<TransientTrack> &tracks,
640  const std::vector<float> &weights) const
641 {
642  GhostTrack ghostTrack_(ghostTrack, tracks, weights,
643  dummyPrediction(primaryVertex, ghostTrack),
645  primaryVertex.position()),
646  true);
647 
648  CachingVertex<5> primary(
649  RecoVertex::convertPos(primaryVertex.position()),
650  RecoVertex::convertError(primaryVertex.error()),
651  std::vector<RefCountedVertexTrack>(), 0.);
652 
653  std::vector<TransientVertex> result = vertices(ghostTrack_, primary);
654 
655 #ifdef DEBUG
656  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
657  iter != result.end(); ++iter)
658  debugVertex(*iter, ghostTrack_.prediction());
659 #endif
660 
661  return result;
662 }
663 
664 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
665  const Vertex &primaryVertex,
666  const Track &ghostTrack,
667  const BeamSpot &beamSpot,
668  const std::vector<TransientTrack> &tracks,
669  const std::vector<float> &weights) const
670 {
671  GhostTrack ghostTrack_(ghostTrack, tracks, weights,
672  dummyPrediction(primaryVertex, ghostTrack),
673  RecoVertex::convertPos(primaryVertex.position()),
674  true);
675 
676  CachingVertex<5> primary(
677  RecoVertex::convertPos(primaryVertex.position()),
678  RecoVertex::convertError(primaryVertex.error()),
679  std::vector<RefCountedVertexTrack>(), 0.);
680 
681  std::vector<TransientVertex> result =
682  vertices(ghostTrack_, primary, beamSpot, true);
683 
684 #ifdef DEBUG
685  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
686  iter != result.end(); ++iter)
687  debugVertex(*iter, ghostTrack_.prediction());
688 #endif
689 
690  return result;
691 }
692 
693 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
694  const Vertex &primaryVertex,
695  const Track &ghostTrack,
696  const BeamSpot &beamSpot,
697  const std::vector<TransientTrack> &primaries,
698  const std::vector<TransientTrack> &tracks,
699  const std::vector<float> &weights) const
700 {
701  GhostTrack ghostTrack_(ghostTrack, tracks, weights,
702  dummyPrediction(primaryVertex, ghostTrack),
703  RecoVertex::convertPos(primaryVertex.position()),
704  true);
705 
706  std::vector<RefCountedVertexTrack> primaryVertexTracks;
707  if (!primaries.empty()) {
708  LinearizedTrackStateFactory linTrackFactory;
709  VertexTrackFactory<5> vertexTrackFactory;
710 
711  VertexState state(
712  RecoVertex::convertPos(primaryVertex.position()),
713  RecoVertex::convertError(primaryVertex.error()));
714 
715  for(std::vector<TransientTrack>::const_iterator iter =
716  primaries.begin(); iter != primaries.end(); ++iter) {
717 
718  RefCountedLinearizedTrackState linState =
719  linTrackFactory.linearizedTrackState(
720  state.position(), *iter);
721 
722  primaryVertexTracks.push_back(
723  vertexTrackFactory.vertexTrack(
724  linState, state));
725  }
726  }
727 
728  CachingVertex<5> primary(
729  RecoVertex::convertPos(primaryVertex.position()),
730  RecoVertex::convertError(primaryVertex.error()),
731  primaryVertexTracks, 0.);
732 
733  std::vector<TransientVertex> result =
734  vertices(ghostTrack_, primary, beamSpot, true, true);
735 
736 #ifdef DEBUG
737  for(std::vector<TransientVertex>::const_iterator iter = result.begin();
738  iter != result.end(); ++iter)
739  debugVertex(*iter, ghostTrack_.prediction());
740 #endif
741 
742  return result;
743 }
744 
745 static double fitChi2(const CachingVertex<5> &vtx)
746 {
747  return vtx.totalChiSquared() / vtx.degreesOfFreedom();
748 }
749 
750 std::vector<CachingVertex<5> > GhostTrackVertexFinder::initialVertices(
751  const FinderInfo &info) const
752 {
753  std::vector<CachingVertex<5> > vertices;
754  for(std::vector<GhostTrackState>::const_iterator iter =
755  info.states.begin(); iter != info. states.end(); ++iter) {
756 
757  if (!iter->isValid())
758  continue;
759 
760  GhostTrackState state(*iter);
761  state.linearize(info.pred);
762 
763  if (!state.isValid())
764  continue;
765 
767  info.pred, state);
768 
769  if (vtx.isValid()) // && fitChi2(vtx) < maxFitChi2_)
770  vertices.push_back(vtx);
771  }
772 
773  return vertices;
774 }
775 
776 static void mergeTrackHelper(const std::vector<RefCountedVertexTrack> &tracks,
777  std::vector<RefCountedVertexTrack> &newTracks,
778  const VertexState &state,
779  const VtxTrackIs &ghostTrackFinder,
780  RefCountedVertexTrack &ghostTrack,
781  const VertexTrackFactory<5> &factory)
782 {
783  for(std::vector<RefCountedVertexTrack>::const_iterator iter =
784  tracks.begin(); iter != tracks.end(); ++iter) {
785  bool gt = ghostTrackFinder(*iter);
786  if (gt && ghostTrack)
787  continue;
788 
789  RefCountedVertexTrack track =
790  relinearizeTrack(*iter, state, factory);
791 
792  if (gt)
793  ghostTrack = *iter;
794  else
795  newTracks.push_back(*iter);
796  }
797 }
798 
800  const CachingVertex<5> &vertex1,
801  const CachingVertex<5> &vertex2,
802  const FinderInfo &info,
803  bool isPrimary) const
804 {
805  VertexTrackFactory<5> vertexTrackFactory;
806 
807  VertexState state = stateMean(vertex1.vertexState(),
808  vertex2.vertexState());
809  std::vector<RefCountedVertexTrack> linTracks;
810  VtxTrackIs isGhostTrack(info.ghostTrack);
811  RefCountedVertexTrack vtxGhostTrack;
812 
813  mergeTrackHelper(vertex1.tracks(), linTracks, state,
814  isGhostTrack, vtxGhostTrack, vertexTrackFactory);
815  mergeTrackHelper(vertex2.tracks(), linTracks, state,
816  isGhostTrack, vtxGhostTrack, vertexTrackFactory);
817 
818  if (vtxGhostTrack &&
819  (fitType_ == kAlwaysWithGhostTrack || linTracks.size() < 2))
820  linTracks.push_back(
821  vertexTrackFactory.vertexTrack(
822  vtxGhostTrack->linearizedTrack(),
823  vtxGhostTrack->vertexState()));
824 
825  try {
827  if (info.hasBeamSpot && isPrimary)
828  vtx = vertexFitter(true).vertex(linTracks,
829  info.beamSpot.position(),
830  info.beamSpot.error());
831  else
832  vtx = vertexFitter(isPrimary).vertex(linTracks);
833  if (vtx.isValid())
834  return vtx;
835  } catch(const VertexException &e) {
836  // fit failed
837  }
838 
839  return CachingVertex<5>();
840 }
841 
843  std::vector<CachingVertex<5> > &vertices,
844  const FinderInfo &info) const
845 {
846  typedef std::pair<unsigned int, unsigned int> Indices;
847 
848  std::multimap<float, Indices> compatMap;
849  unsigned int n = vertices.size();
850  for(unsigned int i = 0; i < n; i++) {
851  const CachingVertex<5> &v1 = vertices[i];
852  for(unsigned int j = i + 1; j < n; j++) {
853  const CachingVertex<5> &v2 = vertices[j];
854 
855  float compat = vertexCompat(v1, v2, info,
856  i == 0 ? primcut_ : seccut_, seccut_);
857 
858  if (compat > mergeThreshold_)
859  continue;
860 
861  compatMap.insert(
862  std::make_pair(compat, Indices(i, j)));
863  }
864  }
865 
866  bool changesMade = false;
867  bool repeat = true;
868  while(repeat) {
869  repeat = false;
870  for(std::multimap<float, Indices>::const_iterator iter =
871  compatMap.begin(); iter != compatMap.end(); ++iter) {
872  unsigned int v1 = iter->second.first;
873  unsigned int v2 = iter->second.second;
874 
875  CachingVertex<5> newVtx =
877  info, v1 == 0);
878  if (!newVtx.isValid() ||
879  (v1 != 0 && fitChi2(newVtx) > maxFitChi2_))
880  continue;
881 
882  std::swap(vertices[v1], newVtx);
883  vertices.erase(vertices.begin() + v2);
884  n--;
885 
886  std::multimap<float, Indices> newCompatMap;
887  for(++iter; iter != compatMap.end(); ++iter) {
888  if (iter->second.first == v1 ||
889  iter->second.first == v2 ||
890  iter->second.second == v1 ||
891  iter->second.second == v2)
892  continue;
893 
894  Indices indices = iter->second;
895  indices.first -= indices.first > v2;
896  indices.second -= indices.second > v2;
897 
898  newCompatMap.insert(std::make_pair(
899  iter->first, indices));
900  }
901  std::swap(compatMap, newCompatMap);
902 
903  for(unsigned int i = 0; i < n; i++) {
904  if (i == v1)
905  continue;
906 
907  const CachingVertex<5> &other = vertices[i];
908  float compat = vertexCompat(
909  vertices[v1], other, info,
910  v1 == 0 ? primcut_ : seccut_,
911  i == 0 ? primcut_ : seccut_);
912 
913  if (compat > mergeThreshold_)
914  continue;
915 
916  compatMap.insert(
917  std::make_pair(
918  compat,
919  Indices(std::min(i, v1),
920  std::max(i, v1))));
921  }
922 
923  changesMade = true;
924  repeat = true;
925  break;
926  }
927  }
928 
929  return changesMade;
930 }
931 
933  std::vector<CachingVertex<5> > &vertices_,
934  const FinderInfo &info) const
935 {
936  std::vector<std::pair<RefCountedVertexTrack,
937  std::vector<RefCountedVertexTrack> > >
938  trackBundles(vertices_.size());
939 
940  VtxTrackIs isGhostTrack(info.ghostTrack);
941 
942  bool assignmentChanged = false;
943  for(std::vector<CachingVertex<5> >::const_iterator iter =
944  vertices_.begin(); iter != vertices_.end(); ++iter) {
945  std::vector<RefCountedVertexTrack> vtxTracks = iter->tracks();
946 
947  if (vtxTracks.empty()) {
948  LinearizedTrackStateFactory linTrackFactory;
949  VertexTrackFactory<5> vertexTrackFactory;
950 
951  RefCountedLinearizedTrackState linState =
952  linTrackFactory.linearizedTrackState(
953  iter->position(), info.ghostTrack);
954 
955  trackBundles[iter - vertices_.begin()].first =
956  vertexTrackFactory.vertexTrack(
957  linState, iter->vertexState());
958  }
959 
960  for(std::vector<RefCountedVertexTrack>::const_iterator track =
961  vtxTracks.begin(); track != vtxTracks.end(); ++track) {
962 
963  if (isGhostTrack(*track)) {
964  trackBundles[iter - vertices_.begin()]
965  .first = *track;
966  continue;
967  }
968 
969  if ((*track)->weight() < 1e-3) {
970  trackBundles[iter - vertices_.begin()]
971  .second.push_back(*track);
972  continue;
973  }
974 
975  unsigned int idx = iter - vertices_.begin();
976  double best = 1.0e9;
977  for(std::vector<CachingVertex<5> >::const_iterator
978  vtx = vertices_.begin();
979  vtx != vertices_.end(); ++vtx) {
980  if (info.pred.lambda(vtx->position()) <
981  info.pred.lambda(vertices_[0].position()))
982  continue;
983 
984  double compat =
985  sqr(trackVertexCompat(*vtx, *track));
986 
987  compat /= (vtx == vertices_.begin()) ?
988  primcut_ : seccut_;
989 
990  if (compat < best) {
991  best = compat;
992  idx = vtx - vertices_.begin();
993  }
994  }
995 
996  if ((int)idx != iter - vertices_.begin())
997  assignmentChanged = true;
998 
999  trackBundles[idx].second.push_back(*track);
1000  }
1001  }
1002 
1003  if (!assignmentChanged)
1004  return false;
1005 
1006  VertexTrackFactory<5> vertexTrackFactory;
1007  std::vector<CachingVertex<5> > vertices;
1008  vertices.reserve(vertices_.size());
1009 
1010  for(std::vector<CachingVertex<5> >::const_iterator iter =
1011  vertices_.begin(); iter != vertices_.end(); ++iter) {
1012  const std::vector<RefCountedVertexTrack> &tracks =
1013  trackBundles[iter - vertices_.begin()].second;
1014  if (tracks.empty())
1015  continue;
1016 
1018 
1019  if (tracks.size() == 1) {
1020  const TransientTrack &track =
1021  tracks[0]->linearizedTrack()->track();
1022 
1023  int idx = -1;
1024  for(std::vector<GhostTrackState>::const_iterator iter =
1025  info.states.begin();
1026  iter != info.states.end(); ++iter) {
1027  if (iter->isTrack() && iter->track() == track) {
1028  idx = iter - info.states.begin();
1029  break;
1030  }
1031  }
1032  if (idx < 0)
1033  continue;
1034 
1035  vtx = vertexAtState(info.ghostTrack, info.pred,
1036  info.states[idx]);
1037  if (!vtx.isValid())
1038  continue;
1039  } else {
1040  std::vector<RefCountedVertexTrack> linTracks =
1041  relinearizeTracks(tracks, iter->vertexState());
1042 
1044  linTracks.push_back(relinearizeTrack(
1045  trackBundles[iter - vertices_.begin()].first,
1046  iter->vertexState(), vertexTrackFactory));
1047 
1048  bool primary = iter == vertices_.begin();
1049  try {
1050  if (primary && info.hasBeamSpot)
1051  vtx = vertexFitter(true).vertex(
1052  linTracks,
1053  info.beamSpot.position(),
1054  info.beamSpot.error());
1055  else
1056  vtx = vertexFitter(primary).vertex(
1057  linTracks);
1058  } catch(const VertexException &e) {
1059  // fit failed;
1060  }
1061  if (!vtx.isValid())
1062  return false;
1063  }
1064 
1065  vertices.push_back(vtx);
1066  };
1067 
1068  std::swap(vertices_, vertices);
1069  return true;
1070 }
1071 
1073  std::vector<CachingVertex<5> > &vertices,
1074  FinderInfo &info) const
1075 {
1076  VtxTrackIs isGhostTrack(info.ghostTrack);
1077 
1078  std::vector<GhostTrackState> states;
1079  std::vector<unsigned int> oldStates;
1080  oldStates.reserve(info.states.size());
1081 
1082  for(std::vector<CachingVertex<5> >::const_iterator iter =
1083  vertices.begin(); iter != vertices.end(); ++iter) {
1084  std::vector<RefCountedVertexTrack> vtxTracks = iter->tracks();
1085 
1086  oldStates.clear();
1087  for(std::vector<RefCountedVertexTrack>::const_iterator track =
1088  vtxTracks.begin(); track != vtxTracks.end(); ++track) {
1089 
1090  if (isGhostTrack(*track) || (*track)->weight() < 1e-3)
1091  continue;
1092 
1093  const TransientTrack &tt =
1094  (*track)->linearizedTrack()->track();
1095 
1096  int idx = -1;
1097  for(std::vector<GhostTrackState>::const_iterator iter =
1098  info.states.begin();
1099  iter != info.states.end(); ++iter) {
1100  if (iter->isTrack() && iter->track() == tt) {
1101  idx = iter - info.states.begin();
1102  break;
1103  }
1104  }
1105 
1106  if (idx >= 0)
1107  oldStates.push_back(idx);
1108  }
1109 
1110  if (oldStates.size() == 1)
1111  states.push_back(info.states[oldStates[0]]);
1112  else
1113  states.push_back(GhostTrackState(iter->vertexState()));
1114  }
1115 
1116  KalmanGhostTrackUpdater updater;
1118  double ndof, chi2;
1119  info.pred = fitter.fit(updater, info.prior, states, ndof, chi2);
1120  TransientTrack ghostTrack = transientGhostTrack(info.pred, info.field);
1121 
1122  std::swap(info.states, states);
1123  states.clear();
1124 
1125  std::vector<CachingVertex<5> > newVertices;
1126  for(std::vector<CachingVertex<5> >::const_iterator iter =
1127  vertices.begin(); iter != vertices.end(); ++iter) {
1128  std::vector<RefCountedVertexTrack> vtxTracks = iter->tracks();
1129 
1130  int idx = -1;
1131  bool redo = false;
1132  for(std::vector<RefCountedVertexTrack>::iterator track =
1133  vtxTracks.begin(); track != vtxTracks.end(); ++track) {
1134 
1135  if (isGhostTrack(*track)) {
1136  LinearizedTrackStateFactory linTrackFactory;
1137  VertexTrackFactory<5> vertexTrackFactory;
1138 
1139  *track = vertexTrackFactory.vertexTrack(
1140  linTrackFactory.linearizedTrackState(
1141  iter->position(), ghostTrack),
1142  iter->vertexState());
1143  redo = true;
1144  continue;
1145  }
1146 
1147  const TransientTrack &tt =
1148  (*track)->linearizedTrack()->track();
1149 
1150  if (idx >= 0) {
1151  idx = -1;
1152  break;
1153  }
1154 
1155  for(std::vector<GhostTrackState>::const_iterator it =
1156  info.states.begin();
1157  it != info.states.end(); ++it) {
1158  if (!it->isTrack())
1159  continue;
1160 
1161  if (it->track() == tt) {
1162  idx = it - info.states.begin();
1163  break;
1164  }
1165  }
1166  }
1167 
1168  if (idx >= 0) {
1169  CachingVertex<5> vtx =
1170  vertexAtState(ghostTrack, info.pred,
1171  info.states[idx]);
1172  if (vtx.isValid())
1173  newVertices.push_back(vtx);
1174  } else if (redo) {
1175  bool primary = iter == vertices.begin();
1177  if (primary && info.hasBeamSpot)
1178  vtx = vertexFitter(true).vertex(
1179  vtxTracks,
1180  info.beamSpot.position(),
1181  info.beamSpot.error());
1182  else
1183  vtx = vertexFitter(primary).vertex(vtxTracks);
1184  if (vtx.isValid())
1185  newVertices.push_back(vtx);
1186  } else
1187  newVertices.push_back(*iter);
1188  }
1189 
1190  std::swap(newVertices, vertices);
1191  info.ghostTrack = ghostTrack;
1192 }
1193 
1194 // implementation
1195 
1196 std::vector<TransientVertex> GhostTrackVertexFinder::vertices(
1197  const GhostTrack &ghostTrack,
1198  const CachingVertex<5> &primary,
1199  const BeamSpot &beamSpot,
1200  bool hasBeamSpot, bool hasPrimaries) const
1201 {
1202  FinderInfo info(primary, ghostTrack, beamSpot,
1203  hasBeamSpot, hasPrimaries);
1204 
1205  std::vector<TransientVertex> result;
1206  if (info.states.empty())
1207  return result;
1208 
1209  info.field = info.states[0].track().field();
1210  info.ghostTrack = transientGhostTrack(info.pred, info.field);
1211 
1212  std::vector<CachingVertex<5> > vertices = initialVertices(info);
1213  if (primary.isValid()) {
1214  vertices.push_back(primary);
1215  if (vertices.size() > 1)
1216  std::swap(vertices.front(), vertices.back());
1217  }
1218 
1219  unsigned int reassigned = 0;
1220  while(reassigned < 3) {
1221  if (vertices.size() < 2)
1222  break;
1223 
1224 #ifdef DEBUG
1225  for(std::vector<CachingVertex<5> >::const_iterator iter =
1226  vertices.begin(); iter != vertices.end(); ++iter)
1227 
1228  debugVertex(*iter, ghostTrack.prediction());
1229 
1230  std::cout << "----- recursive merging: ---------" << std::endl;
1231 #endif
1232 
1233  bool changed = recursiveMerge(vertices, info);
1234  if ((!changed && !reassigned) || vertices.size() < 2)
1235  break;
1236  if (changed)
1237  reassigned = 0;
1238 
1240  refitGhostTrack(vertices, info);
1241  changed = true;
1242  }
1243 
1244  try {
1245 #ifdef DEBUG
1246  for(std::vector<CachingVertex<5> >::const_iterator
1247  iter = vertices.begin();
1248  iter != vertices.end(); ++iter)
1249  debugVertex(*iter, ghostTrack.prediction());
1250  std::cout << "----- reassignment: ---------" << std::endl;
1251 #endif
1252  if (reassignTracks(vertices, info)) {
1253  reassigned++;
1254  changed = true;
1255  } else
1256  reassigned = 0;
1257  } catch(const VertexException &e) {
1258  // just keep vertices as they are
1259  }
1260 
1261  if (!changed)
1262  break;
1263  }
1264 
1265  for(std::vector<CachingVertex<5> >::const_iterator iter =
1266  vertices.begin(); iter != vertices.end(); ++iter) {
1267  std::vector<RefCountedVertexTrack> tracks = iter->tracks();
1268  std::vector<RefCountedVertexTrack> newTracks;
1269  newTracks.reserve(tracks.size());
1270 
1271  std::remove_copy_if(tracks.begin(), tracks.end(),
1272  std::back_inserter(newTracks),
1273  VtxTrackIs(info.ghostTrack));
1274 
1275  if (newTracks.empty())
1276  continue;
1277 
1278  CachingVertex<5> vtx(iter->vertexState(), newTracks,
1279  iter->totalChiSquared());
1280  result.push_back(vtx);
1281  }
1282 
1283  return result;
1284 }
1285 
GlobalError positionError() const
reco::Vertex::Point convertPos(const GlobalPoint &p)
std::vector< RefCountedVertexTrack > tracks() const
static const TGPicture * info(bool iBackgroundIsBlack)
std::unique_ptr< VertexFitter< 5 > > primVertexFitter_
static bool covarianceUpdate(Matrix3S &cov, const Vector3 &residual, const Matrix3S &error, double &chi2, double theta, double phi)
double lambda(const GlobalPoint &point) const
const TransientTrack & track() const
static HepMC::IO_HEPEVT conv
void add(const reco::Track &track, double weight=1.0)
TrackBaseRef trackBaseRef() const
RefCountedVertexTrack vertexTrack(const RefCountedLinearizedTrackState lt, const VertexState vs, float weight=1.0) const
Common base class.
VertexState const & vertexState() const
const AlgebraicSymMatrix33 matrix() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
float totalChiSquared() const
static RefCountedVertexTrack relinearizeTrack(const RefCountedVertexTrack &track, const VertexState &state, const VertexTrackFactory< 5 > factory)
#define nullptr
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:63
FreeTrajectoryState fts(const MagneticField *fieldProvider) const
Definition: weight.py:1
static double vertexCompat(const CachingVertex< 5 > &vtx1, const CachingVertex< 5 > &vtx2, const FinderInfo &info, double scale1=1.0, double scale2=1.0)
const MagneticField * field() const
bool reassignTracks(std::vector< CachingVertex< 5 > > &vertices, const FinderInfo &info) const
bool recursiveMerge(std::vector< CachingVertex< 5 > > &vertices, const FinderInfo &info) const
GlobalError cartesianError() const
TrackCharge charge() const
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:660
static TransientTrack transientGhostTrack(const GhostTrackPrediction &pred, const MagneticField *field)
static double fitChi2(const CachingVertex< 5 > &vtx)
GlobalPoint position() const
Definition: VertexState.h:69
const Point & position() const
position
Definition: Vertex.h:109
A arg
Definition: Factorize.h:38
static VertexState stateMean(const VertexState &v1, const VertexState &v2)
GhostTrack fit(const GlobalPoint &priorPosition, const GlobalError &priorError, const GlobalVector &direction, double coneRadius, const std::vector< TransientTrack > &tracks) const
GhostTrackPrediction fit(const GhostTrackFitter::PredictionUpdater &updater, const GhostTrackPrediction &prior, std::vector< GhostTrackState > &states, double &ndof, double &chi2) override
std::vector< TransientVertex > vertices(const reco::Vertex &primaryVertex, const GlobalVector &direction, double coneRadius, const std::vector< TransientTrack > &tracks) const
Indices
Definition: EdmEventSize.cc:29
static std::vector< RefCountedVertexTrack > relinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const VertexState &state)
std::vector< reco::TransientTrack > const & originalTracks() const
float degreesOfFreedom() const
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:18
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
float totalChiSquared() const
static GhostTrackPrediction dummyPrediction(const Vertex &primaryVertex, const Track &ghostTrack)
static GlobalPoint vtxMean(const GlobalPoint &p1, const GlobalError &e1, const GlobalPoint &p2, const GlobalError &e2)
float degreesOfFreedom() const
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const =0
T min(T a, T b)
Definition: MathUtil.h:58
double lambda() const
double p2[4]
Definition: TauolaWrapper.h:90
GlobalVector momentum() const
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:672
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:27
FreeTrajectoryState initialFreeState() const
RefCountedLinearizedTrackState linearizedTrackState(const GlobalPoint &linP, const reco::TransientTrack &track) const override
const Track & track() const
CachingVertex< 5 > mergeVertices(const CachingVertex< 5 > &vertex1, const CachingVertex< 5 > &vertex2, const FinderInfo &info, bool isPrimary) const
bool linearize(const GhostTrackPrediction &pred, bool initial=false, double lambda=0.)
reco::TransientTrack build(const FreeTrajectoryState &fts) const
std::unique_ptr< VertexFitter< 5 > > secVertexFitter_
std::vector< CachingVertex< 5 > > initialVertices(const FinderInfo &info) const
Error error() const
return SMatrix
Definition: Vertex.h:139
FinderInfo(const CachingVertex< 5 > &primaryVertex, const GhostTrack &ghostTrack, const BeamSpot &beamSpot, bool hasBeamSpot, bool hasPrimaries)
VertexFitter< 5 > & vertexFitter(bool primary) const
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:17
GlobalPoint position() const
static double trackVertexCompat(const CachingVertex< 5 > &vtx, const RefCountedVertexTrack &vertexTrack)
GhostTrackFitter & ghostTrackFitter() const
fixed size matrix
const std::vector< reco::CandidatePtr > & tracks_
double p1[4]
Definition: TauolaWrapper.h:89
const GhostTrackPrediction & prediction() const
Definition: GhostTrack.h:41
const math::XYZTLorentzVector & vectorSum() const
Square< F >::type sqr(const F &f)
Definition: Square.h:13
GlobalPoint position(double lambda=0.) const
void refitGhostTrack(std::vector< CachingVertex< 5 > > &vertices, FinderInfo &info) const
bool isValid() const
static double vtxErrorLong(const GlobalError &error, const GlobalVector &dir)
static void mergeTrackHelper(const std::vector< RefCountedVertexTrack > &tracks, std::vector< RefCountedVertexTrack > &newTracks, const VertexState &state, const VtxTrackIs &ghostTrackFinder, RefCountedVertexTrack &ghostTrack, const VertexTrackFactory< 5 > &factory)
GlobalError error() const
Definition: VertexState.h:74
TrajectoryStateOnSurface impactPointState() const
std::unique_ptr< GhostTrackFitter > ghostTrackFitter_
static CachingVertex< 5 > vertexAtState(const TransientTrack &ghostTrack, const GhostTrackPrediction &pred, const GhostTrackState &state)
dbl *** dir
Definition: mlp_gen.cc:35
T x() const
Definition: PV3DBase.h:62
GlobalPoint globalPosition() const
GlobalError error() const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:666
Global3DVector GlobalVector
Definition: GlobalVector.h:10
GlobalError positionError(double lambda=0.) const
const GlobalVector direction() const