CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrimmedVertexFitter.cc
Go to the documentation of this file.
5 
6 
7 
9 {
11  setPtCut(0.);
12 }
13 
15 {
17  setPtCut(pSet.getParameter<double>("minPt"));
18  setTrackCompatibilityCut ( pSet.getParameter<double>("trackCompatibilityCut") );
19  setVertexFitProbabilityCut ( pSet.getParameter<double>("vtxFitProbCut") );
20 }
21 
23 TrimmedVertexFitter::vertex(const std::vector<reco::TransientTrack> & tracks) const
24 {
25  std::vector<TransientVertex> vtces = theRector.vertices ( tracks );
26  if (vtces.size() )
27  {
28  const TransientVertex & rv = *(vtces.begin());
31  VertexState state ( rv.position(), rv.positionError() );
32  std::vector < RefCountedVertexTrack > vtrks;
33  std::vector<reco::TransientTrack> mytrks = rv.originalTracks();
34  for ( std::vector<reco::TransientTrack>::const_iterator rt=mytrks.begin();
35  rt!=mytrks.end() ; ++rt )
36  {
38  ( rv.position(), *rt );
39 
40  RefCountedVertexTrack vtrk = vfac.vertexTrack ( lstate, state, 1.0 );
41  vtrks.push_back ( vtrk );
42  };
43  return CachingVertex<5> ( rv.position(), rv.positionError(), vtrks, rv.totalChiSquared() );
44  };
45  return CachingVertex<5>();
46 }
47 
49  const std::vector<RefCountedVertexTrack> & tracks) const
50 {
51  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
52  throw VertexException("not implemented");
53 }
54 
56  const std::vector<RefCountedVertexTrack> & tracks,
57  const reco::BeamSpot & spot ) const
58 {
59  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
60  throw VertexException("not implemented");
61 }
62 
63 
65  const std::vector<reco::TransientTrack> & tracks, const GlobalPoint& linPoint) const
66 {
67  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
68  throw VertexException("not implemented");
69 }
70 
72  const std::vector<reco::TransientTrack> & tracks, const GlobalPoint& priorPos,
73  const GlobalError& priorError) const
74 {
75  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
76  throw VertexException("not implemented");
77 }
78 
80  const std::vector<RefCountedVertexTrack> & tracks,
81  const GlobalPoint& priorPos,
82  const GlobalError& priorError) const
83 {
84  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
85  throw VertexException("not implemented");
86 }
87 
89 TrimmedVertexFitter::vertex(const std::vector<reco::TransientTrack> & tracks,
90  const reco::BeamSpot& beamSpot) const
91 {
92  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
93  throw VertexException("not implemented");
94 }
95 
96 
98 {
99  return new TrimmedVertexFitter( * this );
100 }
101 
103 {
104  ptcut = cut;
105  theRector.setPtCut ( cut );
106 }
107 
109 {
111 }
112 
114 {
116 }
117 
GlobalError positionError() const
T getParameter(std::string const &) const
KalmanTrimmedVertexFinder theRector
RefCountedVertexTrack vertexTrack(const RefCountedLinearizedTrackState lt, const VertexState vs, float weight=1.0) const
Common base class.
float totalChiSquared() const
virtual std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
std::vector< reco::TransientTrack > originalTracks() const
GlobalPoint position() const
RefCountedLinearizedTrackState linearizedTrackState(const GlobalPoint &linP, const reco::TransientTrack &track) const
CachingVertex< 5 >::RefCountedVertexTrack RefCountedVertexTrack
TrimmedVertexFitter * clone() const
tuple tracks
Definition: testEve_cfg.py:39
void setTrackCompatibilityCut(float cut)
char state
Definition: procUtils.cc:75
tuple cout
Definition: gather_cfg.py:121
void setVertexFitProbabilityCut(float cut)