CMS 3D CMS Logo

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