CMS 3D CMS Logo

TrimmedVertexFitter.cc
Go to the documentation of this file.
5 
8  setPtCut(0.);
9 }
10 
13  setPtCut(pSet.getParameter<double>("minPt"));
14  setTrackCompatibilityCut(pSet.getParameter<double>("trackCompatibilityCut"));
15  setVertexFitProbabilityCut(pSet.getParameter<double>("vtxFitProbCut"));
16 }
17 
18 CachingVertex<5> TrimmedVertexFitter::vertex(const std::vector<reco::TransientTrack>& tracks) const {
19  std::vector<TransientVertex> vtces = theRector.vertices(tracks);
20  if (!vtces.empty()) {
21  const TransientVertex& rv = *(vtces.begin());
24  VertexState state(rv.position(), rv.positionError());
25  std::vector<RefCountedVertexTrack> vtrks;
26  std::vector<reco::TransientTrack> mytrks = rv.originalTracks();
27  for (std::vector<reco::TransientTrack>::const_iterator rt = mytrks.begin(); rt != mytrks.end(); ++rt) {
29 
30  RefCountedVertexTrack vtrk = vfac.vertexTrack(lstate, state, 1.0);
31  vtrks.push_back(vtrk);
32  };
33  return CachingVertex<5>(rv.position(), rv.positionError(), vtrks, rv.totalChiSquared());
34  };
35  return CachingVertex<5>();
36 }
37 
38 CachingVertex<5> TrimmedVertexFitter::vertex(const std::vector<RefCountedVertexTrack>& tracks) const {
39  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
40  throw VertexException("not implemented");
41 }
42 
43 CachingVertex<5> TrimmedVertexFitter::vertex(const std::vector<RefCountedVertexTrack>& tracks,
44  const reco::BeamSpot& spot) const {
45  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
46  throw VertexException("not implemented");
47 }
48 
49 CachingVertex<5> TrimmedVertexFitter::vertex(const std::vector<reco::TransientTrack>& tracks,
50  const GlobalPoint& linPoint) const {
51  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
52  throw VertexException("not implemented");
53 }
54 
55 CachingVertex<5> TrimmedVertexFitter::vertex(const std::vector<reco::TransientTrack>& tracks,
56  const GlobalPoint& priorPos,
57  const GlobalError& priorError) const {
58  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
59  throw VertexException("not implemented");
60 }
61 
62 CachingVertex<5> TrimmedVertexFitter::vertex(const std::vector<RefCountedVertexTrack>& tracks,
63  const GlobalPoint& priorPos,
64  const GlobalError& priorError) const {
65  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
66  throw VertexException("not implemented");
67 }
68 
69 CachingVertex<5> TrimmedVertexFitter::vertex(const std::vector<reco::TransientTrack>& tracks,
70  const reco::BeamSpot& beamSpot) const {
71  std::cout << "[TrimmedVertexFitter] method not implemented" << std::endl;
72  throw VertexException("not implemented");
73 }
74 
76 
78  ptcut = cut;
80 }
81 
83 
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
pwdgSkimBPark_cfi.beamSpot
beamSpot
Definition: pwdgSkimBPark_cfi.py:5
MessageLogger.h
TkAlMuonSelectors_cfi.cut
cut
Definition: TkAlMuonSelectors_cfi.py:5
VertexException
Common base class.
Definition: VertexException.h:12
LinearizedTrackStateFactory.h
CachingVertex< 5 >
KalmanTrimmedVertexFinder::setTrackCompatibilityCut
void setTrackCompatibilityCut(float cut)
Definition: KalmanTrimmedVertexFinder.h:58
gather_cfg.cout
cout
Definition: gather_cfg.py:144
hcal_runs.rt
rt
Definition: hcal_runs.py:76
TransientVertex::position
GlobalPoint position() const
Definition: TransientVertex.h:169
ReferenceCountingPointer< LinearizedTrackState< 5 > >
TrimmedVertexFitter::setVertexFitProbabilityCut
void setVertexFitProbabilityCut(float cut)
Definition: TrimmedVertexFitter.cc:84
TrimmedVertexFitter
Definition: TrimmedVertexFitter.h:12
LinearizedTrackStateFactory::linearizedTrackState
RefCountedLinearizedTrackState linearizedTrackState(const GlobalPoint &linP, const reco::TransientTrack &track) const override
Definition: LinearizedTrackStateFactory.cc:9
KalmanTrimmedVertexFinder::vertices
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
Definition: KalmanTrimmedVertexFinder.h:23
VertexTrackFactory.h
TrimmedVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition: TrimmedVertexFitter.cc:18
TrimmedVertexFitter::setTrackCompatibilityCut
void setTrackCompatibilityCut(float cut)
Definition: TrimmedVertexFitter.cc:82
LinearizedTrackStateFactory
Definition: LinearizedTrackStateFactory.h:14
TrimmedVertexFitter.h
TrimmedVertexFitter::RefCountedVertexTrack
CachingVertex< 5 >::RefCountedVertexTrack RefCountedVertexTrack
Definition: TrimmedVertexFitter.h:14
reco::BeamSpot
Definition: BeamSpot.h:21
TrimmedVertexFitter::theRector
KalmanTrimmedVertexFinder theRector
Definition: TrimmedVertexFitter.h:49
Point3DBase< float, GlobalTag >
edm::ParameterSet
Definition: ParameterSet.h:36
VertexTrackFactory< 5 >
GlobalErrorBase< double, ErrorMatrixTag >
KalmanTrimmedVertexFinder::setVertexFitProbabilityCut
void setVertexFitProbabilityCut(float cut)
Definition: KalmanTrimmedVertexFinder.h:60
TransientVertex
Definition: TransientVertex.h:18
KalmanTrimmedVertexFinder::setMaxNbOfVertices
void setMaxNbOfVertices(int max)
Definition: KalmanTrimmedVertexFinder.h:61
TrimmedVertexFitter::TrimmedVertexFitter
TrimmedVertexFitter()
Definition: TrimmedVertexFitter.cc:6
TrimmedVertexFitter::clone
TrimmedVertexFitter * clone() const override
Definition: TrimmedVertexFitter.cc:75
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
VertexTrackFactory::vertexTrack
RefCountedVertexTrack vertexTrack(const RefCountedLinearizedTrackState lt, const VertexState vs, float weight=1.0) const
Definition: VertexTrackFactory.h:27
VertexState
Definition: VertexState.h:13
KalmanTrimmedVertexFinder::setPtCut
void setPtCut(float cut)
Definition: KalmanTrimmedVertexFinder.h:57
TransientVertex::originalTracks
std::vector< reco::TransientTrack > const & originalTracks() const
Definition: TransientVertex.h:200
TrimmedVertexFitter::setPtCut
void setPtCut(float cut)
Definition: TrimmedVertexFitter.cc:77
TrimmedVertexFitter::ptcut
double ptcut
Definition: TrimmedVertexFitter.h:50
TransientVertex::positionError
GlobalError positionError() const
Definition: TransientVertex.h:170