CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
AdaptiveVertexFitter.cc
Go to the documentation of this file.
8 
9 #include <algorithm>
10 
11 using namespace edm;
12 
13 // #define STORE_WEIGHTS
14 #ifdef STORE_WEIGHTS
15 #include <dataharvester/Writer.h>
16 #endif
17 
18 using namespace std;
19 
20 namespace {
21  void sortTracksByPt(std::vector<reco::TransientTrack>& cont) {
22  auto s = cont.size();
23  float pt2[s];
24  int ind[s];
25  int i = 0;
26  for (auto const& tk : cont) {
27  ind[i] = i;
28  pt2[i++] = tk.impactPointState().globalMomentum().perp2();
29  }
30  //clang can not handle lambdas with variable length arrays
31  auto* p_pt2 = pt2;
32  std::sort(ind, ind + s, [p_pt2](int i, int j) { return p_pt2[i] > p_pt2[j]; });
33  std::vector<reco::TransientTrack> tmp;
34  tmp.reserve(s);
35  for (auto i = 0U; i < s; ++i)
36  tmp.emplace_back(std::move(cont[ind[i]]));
37  cont.swap(tmp);
38  }
39 
40  // typedef ReferenceCountingPointer<VertexTrack<5> > RefCountedVertexTrack;
41  typedef AdaptiveVertexFitter::RefCountedVertexTrack RefCountedVertexTrack;
42 
43  AlgebraicSymMatrix33 initFitError() {
44  // that's how we model the lin pt error for the initial seed!
45  const float initialError = 10000;
47  ret(0, 0) = initialError;
48  ret(1, 1) = initialError;
49  ret(2, 2) = initialError;
50  return ret;
51  }
52 
53  GlobalError const fitError = initFitError();
54 
55  AlgebraicSymMatrix33 initLinePointError() {
56  // that's how we model the error of the linearization point.
57  // for track weighting!
59  ret(0, 0) = .3;
60  ret(1, 1) = .3;
61  ret(2, 2) = 3.;
62  // ret(0,0)=1e-7; ret(1,1)=1e-7; ret(2,2)=1e-7;
63  return ret;
64  }
65 
66  GlobalError const linPointError = initLinePointError();
67 
68  void sortByDistanceToRefPoint(std::vector<RefCountedVertexTrack>& cont, const GlobalPoint ref) {
69  auto s = cont.size();
70  float d2[s];
71  int ind[s];
72  int i = 0;
73  for (auto const& tk : cont) {
74  ind[i] = i;
75  d2[i++] = (tk->linearizedTrack()->track().initialFreeState().position() - ref).mag2();
76  }
77  //clang can not handle lambdas with variable length arrays
78  auto* p_d2 = d2;
79  std::sort(ind, ind + s, [p_d2](int i, int j) { return p_d2[i] < p_d2[j]; });
80  std::vector<RefCountedVertexTrack> tmp;
81  tmp.reserve(s);
82  for (auto i = 0U; i < s; ++i)
83  tmp.emplace_back(std::move(cont[ind[i]]));
84  cont.swap(tmp);
85  }
86 
87 #ifdef STORE_WEIGHTS
88  //NOTE: This is not thread safe
89  map<RefCountedLinearizedTrackState, int> ids;
90  int iter = 0;
91 
92  int getId(const RefCountedLinearizedTrackState& r) {
93  static int ctr = 1;
94  if (ids.count(r) == 0) {
95  ids[r] = ctr++;
96  }
97  return ids[r];
98  }
99 #endif
100 } // namespace
101 
103  const LinearizationPointFinder& linP,
104  const VertexUpdator<5>& updator,
106  const VertexSmoother<5>& smoother,
107  const AbstractLTSFactory<5>& ltsf)
108  : theNr(0),
109  theLinP(linP.clone()),
110  theUpdator(updator.clone()),
111  theSmoother(smoother.clone()),
112  theAssProbComputer(ann.clone()),
113  theComp(crit.clone()),
114  theLinTrkFactory(ltsf.clone()),
115  gsfIntermediarySmoothing_(false) {
116  setParameters();
117 }
118 
120 
122  : theMaxShift(o.theMaxShift),
123  theMaxLPShift(o.theMaxLPShift),
124  theMaxStep(o.theMaxStep),
125  theWeightThreshold(o.theWeightThreshold),
126  theNr(o.theNr),
127  theLinP(o.theLinP->clone()),
128  theUpdator(o.theUpdator->clone()),
129  theSmoother(o.theSmoother->clone()),
130  theAssProbComputer(o.theAssProbComputer->clone()),
131  theComp(o.theComp->clone()),
132  theLinTrkFactory(o.theLinTrkFactory->clone()),
133  gsfIntermediarySmoothing_(o.gsfIntermediarySmoothing_) {}
134 
136  delete theLinP;
137  delete theUpdator;
138  delete theSmoother;
139  delete theAssProbComputer;
140  delete theComp;
141  delete theLinTrkFactory;
142 }
143 
144 void AdaptiveVertexFitter::setParameters(double maxshift, double maxlpshift, unsigned maxstep, double weightthreshold) {
145  theMaxShift = maxshift;
146  theMaxLPShift = maxlpshift;
147  theMaxStep = maxstep;
148  theWeightThreshold = weightthreshold;
149 }
150 
152  setParameters(s.getParameter<double>("maxshift"),
153  s.getParameter<double>("maxlpshift"),
154  s.getParameter<int>("maxstep"),
155  s.getParameter<double>("weightthreshold"));
156 }
157 
158 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack>& unstracks) const {
159  if (unstracks.size() < 2) {
160  LogError("RecoVertex|AdaptiveVertexFitter") << "Supplied fewer than two tracks. Vertex is invalid.";
161  return CachingVertex<5>(); // return invalid vertex
162  };
163  vector<reco::TransientTrack> tracks = unstracks;
164  sortTracksByPt(tracks);
165  // Linearization Point
166  GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
167  // Initial vertex seed, with a very large error matrix
168  VertexState lseed(linP, linPointError);
169  vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, lseed);
170 
171  VertexState seed(linP, fitError);
172  return fit(vtContainer, seed, false);
173 }
174 
175 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<RefCountedVertexTrack>& tracks) const {
176  if (tracks.size() < 2) {
177  LogError("RecoVertex|AdaptiveVertexFitter") << "Supplied fewer than two tracks. Vertex is invalid.";
178  return CachingVertex<5>(); // return invalid vertex
179  };
180  // Initial vertex seed, with a very small weight matrix
181  GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
182  VertexState seed(linP, fitError);
183  return fit(tracks, seed, false);
184 }
185 
186 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<RefCountedVertexTrack>& tracks,
187  const reco::BeamSpot& spot) const {
188  if (tracks.empty()) {
189  LogError("RecoVertex|AdaptiveVertexFitter") << "Supplied no tracks. Vertex is invalid.";
190  return CachingVertex<5>(); // return invalid vertex
191  };
192  VertexState beamSpotState(spot);
193  return fit(tracks, beamSpotState, true);
194 }
195 
199 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack>& tracks,
200  const GlobalPoint& linPoint) const {
201  if (tracks.size() < 2) {
202  LogError("RecoVertex|AdaptiveVertexFitter") << "Supplied fewer than two tracks. Vertex is invalid.";
203  return CachingVertex<5>(); // return invalid vertex
204  };
205  // Initial vertex seed, with a very large error matrix
206  VertexState seed(linPoint, linPointError);
207  vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, seed);
208  VertexState fitseed(linPoint, fitError);
209  return fit(vtContainer, fitseed, false);
210 }
211 
216 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack>& unstracks,
217  const reco::BeamSpot& beamSpot) const {
218  if (unstracks.empty()) {
219  LogError("RecoVertex|AdaptiveVertexFitter") << "Supplied no tracks. Vertex is invalid.";
220  return CachingVertex<5>(); // return invalid vertex
221  };
222 
223  VertexState beamSpotState(beamSpot);
224  vector<RefCountedVertexTrack> vtContainer;
225 
226  vector<reco::TransientTrack> tracks = unstracks;
227  sortTracksByPt(tracks);
228 
229  if (tracks.size() > 1) {
230  // Linearization Point search if there are more than 1 track
231  GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
232  VertexState lpState(linP, linPointError);
233  vtContainer = linearizeTracks(tracks, lpState);
234  } else {
235  // otherwise take the beamspot position.
236  vtContainer = linearizeTracks(tracks, beamSpotState);
237  }
238 
239  return fit(vtContainer, beamSpotState, true);
240 }
241 
247 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<reco::TransientTrack>& tracks,
248  const GlobalPoint& priorPos,
249  const GlobalError& priorError) const
250 
251 {
252  if (tracks.empty()) {
253  LogError("RecoVertex|AdaptiveVertexFitter") << "Supplied no tracks. Vertex is invalid.";
254  return CachingVertex<5>(); // return invalid vertex
255  };
256  VertexState seed(priorPos, priorError);
257  vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, seed);
258  return fit(vtContainer, seed, true);
259 }
260 
265 CachingVertex<5> AdaptiveVertexFitter::vertex(const vector<RefCountedVertexTrack>& tracks,
266  const GlobalPoint& priorPos,
267  const GlobalError& priorError) const {
268  if (tracks.empty()) {
269  LogError("RecoVertex|AdaptiveVertexFitter") << "Supplied no tracks. Vertex is invalid.";
270  return CachingVertex<5>(); // return invalid vertex
271  };
272  VertexState seed(priorPos, priorError);
273  return fit(tracks, seed, true);
274 }
275 
283 vector<AdaptiveVertexFitter::RefCountedVertexTrack> AdaptiveVertexFitter::linearizeTracks(
284  const vector<reco::TransientTrack>& tracks, const VertexState& seed) const {
285  const GlobalPoint& linP(seed.position());
286  vector<RefCountedLinearizedTrackState> lTracks;
287  for (vector<reco::TransientTrack>::const_iterator i = tracks.begin(); i != tracks.end(); ++i) {
288  try {
290  lTracks.push_back(lTrData);
291  } catch (exception& e) {
292  LogInfo("RecoVertex/AdaptiveVertexFitter") << "Exception " << e.what() << " in ::linearizeTracks."
293  << "Your future vertex has just lost a track.";
294  };
295  }
296  return weightTracks(lTracks, seed);
297 }
298 
304 vector<AdaptiveVertexFitter::RefCountedVertexTrack> AdaptiveVertexFitter::reLinearizeTracks(
305  const vector<RefCountedVertexTrack>& tracks, const CachingVertex<5>& vertex) const {
306  const VertexState& seed = vertex.vertexState();
307  GlobalPoint linP = seed.position();
308  vector<RefCountedLinearizedTrackState> lTracks;
309  for (vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
310  try {
312  theLinTrkFactory->linearizedTrackState(linP, (**i).linearizedTrack()->track());
313  /*
314  RefCountedLinearizedTrackState lTrData =
315  (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
316  */
317  lTracks.push_back(lTrData);
318  } catch (exception& e) {
319  LogInfo("RecoVertex/AdaptiveVertexFitter") << "Exception " << e.what() << " in ::relinearizeTracks. "
320  << "Will not relinearize this track.";
321  lTracks.push_back((**i).linearizedTrack());
322  };
323  };
324  return reWeightTracks(lTracks, vertex);
325 }
326 
328 
330  double weight = theAssProbComputer->weight(chi2);
331 
332  if (weight > 1.0) {
333  LogInfo("RecoVertex/AdaptiveVertexFitter") << "Weight " << weight << " > 1.0!";
334  weight = 1.0;
335  };
336 
337  if (weight < 1e-20) {
338  // LogInfo("RecoVertex/AdaptiveVertexFitter") << "Weight " << weight << " < 0.0!";
339  weight = 1e-20;
340  };
341  return weight;
342 }
343 
344 vector<AdaptiveVertexFitter::RefCountedVertexTrack> AdaptiveVertexFitter::reWeightTracks(
345  const vector<RefCountedLinearizedTrackState>& lTracks, const CachingVertex<5>& vertex) const {
346  const VertexState& seed = vertex.vertexState();
347  // cout << "[AdaptiveVertexFitter] now reweight around " << seed.position() << endl;
348  theNr++;
349  // GlobalPoint pos = seed.position();
350 
351  vector<RefCountedVertexTrack> finalTracks;
352  VertexTrackFactory<5> vTrackFactory;
353 #ifdef STORE_WEIGHTS
354  iter++;
355 #endif
356  for (vector<RefCountedLinearizedTrackState>::const_iterator i = lTracks.begin(); i != lTracks.end(); i++) {
357  double weight = 0.;
358  // cout << "[AdaptiveVertexFitter] estimate " << endl;
359  pair<bool, double> chi2Res(false, 0.);
360  try {
361  chi2Res = theComp->estimate(vertex, *i, std::distance(lTracks.begin(), i));
362  } catch (exception const& e) {
363  };
364  // cout << "[AdaptiveVertexFitter] /estimate " << endl;
365  if (!chi2Res.first) {
366  // cout << "[AdaptiveVertexFitter] aie... vertex candidate is at " << vertex.position() << endl;
367  LogInfo("AdaptiveVertexFitter") << "When reweighting, chi2<0. Will add this track with w=0.";
368  // edm::LogInfo("AdaptiveVertexFitter" ) << "pt=" << (**i).track().pt();
369  } else {
370  weight = getWeight(chi2Res.second);
371  }
372 
373  RefCountedVertexTrack vTrData = vTrackFactory.vertexTrack(*i, seed, weight);
374 
375 #ifdef STORE_WEIGHTS
376  map<string, dataharvester::MultiType> m;
377  m["chi2"] = chi2;
378  m["w"] = theAssProbComputer->weight(chi2);
379  m["T"] = theAssProbComputer->currentTemp();
380  m["n"] = iter;
381  m["pos"] = "reweight";
382  m["id"] = getId(*i);
383  dataharvester::Writer::file("w.txt").save(m);
384 #endif
385 
386  finalTracks.push_back(vTrData);
387  }
388  sortByDistanceToRefPoint(finalTracks, vertex.position());
389  // cout << "[AdaptiveVertexFitter] /now reweight" << endl;
390  return finalTracks;
391 }
392 
393 vector<AdaptiveVertexFitter::RefCountedVertexTrack> AdaptiveVertexFitter::weightTracks(
394  const vector<RefCountedLinearizedTrackState>& lTracks, const VertexState& seed) const {
395  theNr++;
396  CachingVertex<5> seedvtx(seed, vector<RefCountedVertexTrack>(), 0.);
400 
401  vector<RefCountedVertexTrack> finalTracks;
402  VertexTrackFactory<5> vTrackFactory;
403 #ifdef STORE_WEIGHTS
404  iter++;
405 #endif
406  for (vector<RefCountedLinearizedTrackState>::const_iterator i = lTracks.begin(); i != lTracks.end(); i++) {
407  double weight = 0.;
408  pair<bool, double> chi2Res = theComp->estimate(seedvtx, *i, std::distance(lTracks.begin(), i));
409  if (!chi2Res.first) {
410  // cout << "[AdaptiveVertexFitter] Aiee! " << endl;
411  LogInfo("AdaptiveVertexFitter") << "When weighting a track, chi2 calculation failed;"
412  << " will add with w=0.";
413  } else {
414  weight = getWeight(chi2Res.second);
415  }
416  RefCountedVertexTrack vTrData = vTrackFactory.vertexTrack(*i, seed, weight);
417 #ifdef STORE_WEIGHTS
418  map<string, dataharvester::MultiType> m;
419  m["chi2"] = chi2;
420  m["w"] = theAssProbComputer->weight(chi2);
421  m["T"] = theAssProbComputer->currentTemp();
422  m["n"] = iter;
423  m["id"] = getId(*i);
424  m["pos"] = "weight";
425  dataharvester::Writer::file("w.txt").save(m);
426 #endif
427  finalTracks.push_back(vTrData);
428  }
429  return finalTracks;
430 }
431 
437 vector<AdaptiveVertexFitter::RefCountedVertexTrack> AdaptiveVertexFitter::reWeightTracks(
438  const vector<RefCountedVertexTrack>& tracks, const CachingVertex<5>& seed) const {
439  vector<RefCountedLinearizedTrackState> lTracks;
440  for (vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
441  lTracks.push_back((**i).linearizedTrack());
442  }
443 
444  return reWeightTracks(lTracks, seed);
445 }
446 
447 /*
448  * The method where the vertex fit is actually done!
449  */
450 
451 CachingVertex<5> AdaptiveVertexFitter::fit(const vector<RefCountedVertexTrack>& tracks,
452  const VertexState& priorSeed,
453  bool withPrior) const {
454  // cout << "[AdaptiveVertexFit] fit with " << tracks.size() << endl;
456 
457  vector<RefCountedVertexTrack> initialTracks;
458  GlobalPoint priorVertexPosition = priorSeed.position();
459  GlobalError priorVertexError = priorSeed.error();
460 
461  CachingVertex<5> returnVertex(priorVertexPosition, priorVertexError, initialTracks, 0);
462  if (withPrior) {
463  returnVertex = CachingVertex<5>(
464  priorVertexPosition, priorVertexError, priorVertexPosition, priorVertexError, initialTracks, 0);
465  }
466 
467  std::vector<RefCountedVertexTrack> globalVTracks = tracks;
468  // sort the tracks, according to distance to seed!
469  sortByDistanceToRefPoint(globalVTracks, priorSeed.position());
470 
471  // main loop through all the VTracks
472  int lpStep = 0;
473  int step = 0;
474 
475  CachingVertex<5> initialVertex = returnVertex;
476 
477  GlobalPoint newPosition = priorVertexPosition;
478  GlobalPoint previousPosition = newPosition;
479 
480  int ns_trks = 0; // number of significant tracks.
481  // If we have only two significant tracks, we return an invalid vertex
482 
483  // cout << "[AdaptiveVertexFit] start " << tracks.size() << endl;
484  /*
485  for ( vector< RefCountedVertexTrack >::const_iterator
486  i=globalVTracks.begin(); i!=globalVTracks.end() ; ++i )
487  {
488  cout << " " << (**i).linearizedTrack()->track().initialFreeState().momentum() << endl;
489  }*/
490  do {
491  ns_trks = 0;
492  CachingVertex<5> fVertex = initialVertex;
493  // cout << "[AdaptiveVertexFit] step " << step << " at " << fVertex.position() << endl;
494  if ((previousPosition - newPosition).transverse() > theMaxLPShift) {
495  // relinearize and reweight.
496  // (reLinearizeTracks also reweights tracks)
497  // cout << "[AdaptiveVertexFit] relinearize at " << returnVertex.position() << endl;
499  returnVertex = theSmoother->smooth(returnVertex);
500  globalVTracks = reLinearizeTracks(globalVTracks, returnVertex);
501  lpStep++;
502  } else if (step) {
503  // reweight, if it is not the first step
504  // cout << "[AdaptiveVertexFit] reweight at " << returnVertex.position() << endl;
506  returnVertex = theSmoother->smooth(returnVertex);
507  globalVTracks = reWeightTracks(globalVTracks, returnVertex);
508  }
509  // cout << "[AdaptiveVertexFit] relinarized, reweighted" << endl;
510  // update sequentially the vertex estimate
511  CachingVertex<5> nVertex;
512  for (vector<RefCountedVertexTrack>::const_iterator i = globalVTracks.begin(); i != globalVTracks.end(); i++) {
513  if ((**i).weight() > 0.)
514  nVertex = theUpdator->add(fVertex, *i);
515  else
516  nVertex = fVertex;
517  if (nVertex.isValid()) {
518  if ((**i).weight() >= theWeightThreshold)
519  ns_trks++;
520 
521  if (fabs(nVertex.position().z()) > 10000. || nVertex.position().perp() > 120.) {
522  // were more than 100 m off!!
523  LogInfo("AdaptiveVertexFitter")
524  << "Vertex candidate just took off to " << nVertex.position() << "! Will discard this update!";
525  // //<< "track pt was " << (**i).linearizedTrack()->track().pt()
526  // << "track momentum was " << (**i).linearizedTrack()->track().initialFreeState().momentum()
527  // << "track position was " << (**i).linearizedTrack()->track().initialFreeState().position()
528  // << "track chi2 was " << (**i).linearizedTrack()->track().chi2()
529  // << "track ndof was " << (**i).linearizedTrack()->track().ndof()
530  // << "track w was " << (**i).weight()
531  // << "track schi2 was " << (**i).smoothedChi2();
532  } else {
533  fVertex = nVertex;
534  }
535  } else {
536  LogInfo("RecoVertex/AdaptiveVertexFitter")
537  << "The updator returned an invalid vertex when adding track " << i - globalVTracks.begin()
538  << ".\n Your vertex might just have lost one good track.";
539  }
540  }
541  previousPosition = newPosition;
542  newPosition = fVertex.position();
543  returnVertex = fVertex;
545  step++;
546  if (step >= theMaxStep)
547  break;
548 
549  } while (
550  // repeat as long as
551  // - vertex moved too much or
552  // - we're not yet annealed
553  (((previousPosition - newPosition).mag() > theMaxShift) || (!(theAssProbComputer->isAnnealed()))));
554 
555  if (theWeightThreshold > 0. && ns_trks < 2 && !withPrior) {
556  LogDebug("AdaptiveVertexFitter") << "fewer than two significant tracks (w>" << theWeightThreshold << ")."
557  << " Fitted vertex is invalid.";
558  return CachingVertex<5>(); // return invalid vertex
559  }
560 
561 #ifdef STORE_WEIGHTS
562  map<string, dataharvester::MultiType> m;
563  m["chi2"] = chi2;
564  m["w"] = theAssProbComputer->weight(chi2);
565  m["T"] = theAssProbComputer->currentTemp();
566  m["n"] = iter;
567  m["id"] = getId(*i);
568  m["pos"] = "final";
569  dataharvester::Writer::file("w.txt").save(m);
570 #endif
571  // cout << "[AdaptiveVertexFit] /fit" << endl;
572  return theSmoother->smooth(returnVertex);
573 }
virtual bool isAnnealed() const =0
AdaptiveVertexFitter * clone() const override
tuple ret
prodAgent to be discontinued
static unsigned int getId()
T perp() const
Definition: PV3DBase.h:69
std::vector< RefCountedVertexTrack > reWeightTracks(const std::vector< RefCountedLinearizedTrackState > &, const CachingVertex< 5 > &seed) const
tuple cont
load Luminosity info ##
Definition: generateEDF.py:628
const double w
Definition: UKUtility.cc:23
virtual RefCountedLinearizedTrackState linearizedTrackState(const GlobalPoint &linP, const reco::TransientTrack &track) const =0
RefCountedVertexTrack vertexTrack(const RefCountedLinearizedTrackState lt, const VertexState vs, float weight=1.0) const
LinearizationPointFinder * theLinP
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
VertexState const & vertexState() const
VertexUpdator< 5 > * theUpdator
auto const & tracks
cannot be loose
CachingVertex< 5 > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState &priorSeed, bool withPrior) const
Log< level::Error, false > LogError
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
GlobalPoint position() const
Definition: VertexState.h:62
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const CachingVertex< 5 > &vertex) const
virtual CachingVertex< N > add(const CachingVertex< N > &v, const typename CachingVertex< N >::RefCountedVertexTrack t) const =0
virtual BDpair estimate(const CachingVertex< N > &v, const RefCountedLinearizedTrackState track, unsigned int hint=UINT_MAX) const =0
T z() const
Definition: PV3DBase.h:61
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
def move
Definition: eostools.py:511
VertexSmoother< 5 > * theSmoother
AnnealingSchedule * theAssProbComputer
Log< level::Info, false > LogInfo
virtual double currentTemp() const =0
void setParameters(double maxshift=0.0001, double maxlpshift=0.1, unsigned maxstep=30, double weightthreshold=.001)
double getWeight(float chi2) const
AdaptiveVertexFitter(const AnnealingSchedule &ann=GeometricAnnealing(), const LinearizationPointFinder &linP=DefaultLinearizationPointFinder(), const VertexUpdator< 5 > &updator=KalmanVertexUpdator< 5 >(), const VertexTrackCompatibilityEstimator< 5 > &estor=KalmanVertexTrackCompatibilityEstimator< 5 >(), const VertexSmoother< 5 > &smoother=DummyVertexSmoother< 5 >(), const AbstractLTSFactory< 5 > &ltsf=LinearizedTrackStateFactory())
virtual CachingVertex< N > smooth(const CachingVertex< N > &vertex) const =0
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
GlobalPoint position() const
virtual GlobalPoint getLinearizationPoint(const std::vector< reco::TransientTrack > &) const =0
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
bool isValid() const
virtual void anneal()=0
virtual void resetAnnealing()=0
GlobalError error() const
Definition: VertexState.h:64
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &, const VertexState &) const
step
Definition: StallMonitor.cc:94
int weight
Definition: histoStyle.py:51
std::vector< RefCountedVertexTrack > weightTracks(const std::vector< RefCountedLinearizedTrackState > &, const VertexState &seed) const
tmp
align.sh
Definition: createJobs.py:716
virtual double weight(double chi2) const =0
const AbstractLTSFactory< 5 > * theLinTrkFactory
VertexTrackCompatibilityEstimator< 5 > * theComp
#define LogDebug(id)