CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
GsfVertexFitter Class Reference

#include <GsfVertexFitter.h>

Inheritance diagram for GsfVertexFitter:
VertexFitter< 5 >

Public Types

typedef CachingVertex< 5 >::RefCountedVertexTrack RefCountedVertexTrack
 

Public Member Functions

GsfVertexFitterclone () const override
 
 GsfVertexFitter (const edm::ParameterSet &pSet, const LinearizationPointFinder &linP=DefaultLinearizationPointFinder())
 
 GsfVertexFitter (const GsfVertexFitter &original)
 
CachingVertex< 5 > vertex (const std::vector< reco::TransientTrack > &tracks) const override
 
CachingVertex< 5 > vertex (const std::vector< RefCountedVertexTrack > &tracks) const override
 
CachingVertex< 5 > vertex (const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &linPoint) const override
 
CachingVertex< 5 > vertex (const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const override
 
CachingVertex< 5 > vertex (const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &beamSpot) const override
 
CachingVertex< 5 > vertex (const std::vector< RefCountedVertexTrack > &tracks, const reco::BeamSpot &spot) const override
 
CachingVertex< 5 > vertex (const std::vector< RefCountedVertexTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const override
 
 ~GsfVertexFitter () override
 
- Public Member Functions inherited from VertexFitter< 5 >
virtual CachingVertex< Nvertex (const std::vector< typename CachingVertex< N >::RefCountedVertexTrack > &tracks) const=0
 
virtual CachingVertex< Nvertex (const std::vector< typename CachingVertex< N >::RefCountedVertexTrack > &tracks, const reco::BeamSpot &spot) const=0
 
virtual CachingVertex< Nvertex (const std::vector< typename CachingVertex< N >::RefCountedVertexTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const=0
 
 VertexFitter ()
 
virtual ~VertexFitter ()
 

Private Attributes

SequentialVertexFitter< 5 > * theSequentialFitter
 

Detailed Description

Sequential vertex fitter, to be used with the Gaussian Sum Vertex Filter After the vertes fit, the tracks can be refit with the additional constraint of the vertex position.

Definition at line 20 of file GsfVertexFitter.h.

Member Typedef Documentation

◆ RefCountedVertexTrack

Definition at line 22 of file GsfVertexFitter.h.

Constructor & Destructor Documentation

◆ GsfVertexFitter() [1/2]

GsfVertexFitter::GsfVertexFitter ( const edm::ParameterSet pSet,
const LinearizationPointFinder linP = DefaultLinearizationPointFinder() 
)

Default constructor, using the given linearization point finder.

Parameters
linPThe LinearizationPointFinder to use
useSmoothingSpecifies whether the tracks sould be refit.

Definition at line 6 of file GsfVertexFitter.cc.

References edm::ParameterSet::getParameter(), SequentialVertexFitter< N >::setMaximumDistance(), SequentialVertexFitter< N >::setMaximumNumberOfIterations(), and theSequentialFitter.

Referenced by clone().

6  {
7  float theMaxShift = pSet.getParameter<double>("maxDistance"); //0.01
8  int theMaxStep = pSet.getParameter<int>("maxNbrOfIterations"); //10
9  bool limitComponents_ = pSet.getParameter<bool>("limitComponents");
10  bool useSmoothing = pSet.getParameter<bool>("smoothTracks");
11 
12  VertexSmoother<5>* theSmoother;
14 
15  if (limitComponents_) {
16  edm::ParameterSet mergerPSet = pSet.getParameter<edm::ParameterSet>("GsfMergerParameters");
17  theMerger = new GsfVertexMerger(mergerPSet);
18  }
19 
20  if (useSmoothing)
21  theSmoother = new GsfVertexSmoother(limitComponents_, &*theMerger);
22  else
23  theSmoother = new DummyVertexSmoother<5>();
24 
26  linP, GsfVertexUpdator(limitComponents_, &*theMerger), *theSmoother, MultiPerigeeLTSFactory());
29 
30  delete theSmoother;
31 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void setMaximumNumberOfIterations(int maxIterations)
SequentialVertexFitter< 5 > * theSequentialFitter
void setMaximumDistance(float maxShift)

◆ ~GsfVertexFitter()

GsfVertexFitter::~GsfVertexFitter ( )
override

Definition at line 37 of file GsfVertexFitter.cc.

References theSequentialFitter.

37 { delete theSequentialFitter; }
SequentialVertexFitter< 5 > * theSequentialFitter

◆ GsfVertexFitter() [2/2]

GsfVertexFitter::GsfVertexFitter ( const GsfVertexFitter original)

Copy constructor

Definition at line 33 of file GsfVertexFitter.cc.

References definitions::original, and theSequentialFitter.

33  {
34  theSequentialFitter = original.theSequentialFitter->clone();
35 }
SequentialVertexFitter< 5 > * theSequentialFitter

Member Function Documentation

◆ clone()

GsfVertexFitter* GsfVertexFitter::clone ( void  ) const
inlineoverridevirtual

Fit vertex out of a VertexSeed

Implements VertexFitter< 5 >.

Definition at line 40 of file GsfVertexFitter.h.

References GsfVertexFitter().

40 { return new GsfVertexFitter(*this); }
GsfVertexFitter(const edm::ParameterSet &pSet, const LinearizationPointFinder &linP=DefaultLinearizationPointFinder())

◆ vertex() [1/7]

CachingVertex<5> GsfVertexFitter::vertex ( const std::vector< reco::TransientTrack > &  tracks) const
inlineoverridevirtual

Fit vertex out of a set of RecTracks

Implements VertexFitter< 5 >.

Definition at line 45 of file GsfVertexFitter.h.

References theSequentialFitter, tracks, and SequentialVertexFitter< N >::vertex().

Referenced by Tau.Tau::dxy().

45  {
47  }
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
SequentialVertexFitter< 5 > * theSequentialFitter
auto const & tracks
cannot be loose

◆ vertex() [2/7]

CachingVertex<5> GsfVertexFitter::vertex ( const std::vector< RefCountedVertexTrack > &  tracks) const
inlineoverride

Fit vertex out of a set of VertexTracks

Definition at line 51 of file GsfVertexFitter.h.

References theSequentialFitter, tracks, and SequentialVertexFitter< N >::vertex().

Referenced by Tau.Tau::dxy().

51  {
53  }
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
SequentialVertexFitter< 5 > * theSequentialFitter
auto const & tracks
cannot be loose

◆ vertex() [3/7]

CachingVertex<5> GsfVertexFitter::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const GlobalPoint linPoint 
) const
inlineoverridevirtual

Fit vertex out of a set of RecTracks. Uses the specified linearization point.

Implements VertexFitter< 5 >.

Definition at line 58 of file GsfVertexFitter.h.

References theSequentialFitter, tracks, and SequentialVertexFitter< N >::vertex().

Referenced by Tau.Tau::dxy().

59  {
60  return theSequentialFitter->vertex(tracks, linPoint);
61  }
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
SequentialVertexFitter< 5 > * theSequentialFitter
auto const & tracks
cannot be loose

◆ vertex() [4/7]

CachingVertex<5> GsfVertexFitter::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const GlobalPoint priorPos,
const GlobalError priorError 
) const
inlineoverridevirtual

Fit vertex out of a set of RecTracks. Uses the specified point as both the linearization point AND as prior estimate of the vertex position. The error is used for the weight of the prior estimate.

Implements VertexFitter< 5 >.

Definition at line 68 of file GsfVertexFitter.h.

References theSequentialFitter, tracks, and SequentialVertexFitter< N >::vertex().

Referenced by Tau.Tau::dxy().

70  {
71  return theSequentialFitter->vertex(tracks, priorPos, priorError);
72  }
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
SequentialVertexFitter< 5 > * theSequentialFitter
auto const & tracks
cannot be loose

◆ vertex() [5/7]

CachingVertex<5> GsfVertexFitter::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const reco::BeamSpot beamSpot 
) const
inlineoverridevirtual

Fit vertex out of a set of TransientTracks. The specified BeamSpot will be used as priot, but NOT for the linearization. The specified LinearizationPointFinder will be used to find the linearization point.

Implements VertexFitter< 5 >.

Definition at line 78 of file GsfVertexFitter.h.

References pwdgSkimBPark_cfi::beamSpot, theSequentialFitter, tracks, and SequentialVertexFitter< N >::vertex().

Referenced by Tau.Tau::dxy().

79  {
81  }
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
SequentialVertexFitter< 5 > * theSequentialFitter
auto const & tracks
cannot be loose

◆ vertex() [6/7]

CachingVertex<5> GsfVertexFitter::vertex ( const std::vector< RefCountedVertexTrack > &  tracks,
const reco::BeamSpot spot 
) const
inlineoverride

Definition at line 83 of file GsfVertexFitter.h.

References theSequentialFitter, tracks, and SequentialVertexFitter< N >::vertex().

Referenced by Tau.Tau::dxy().

84  {
85  return theSequentialFitter->vertex(tracks, spot);
86  }
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
SequentialVertexFitter< 5 > * theSequentialFitter
auto const & tracks
cannot be loose

◆ vertex() [7/7]

CachingVertex<5> GsfVertexFitter::vertex ( const std::vector< RefCountedVertexTrack > &  tracks,
const GlobalPoint priorPos,
const GlobalError priorError 
) const
inlineoverride

Fit vertex out of a set of VertexTracks. Uses the specified point and error as the prior estimate of the vertex. This position is not used to relinearize the tracks.

Definition at line 92 of file GsfVertexFitter.h.

References theSequentialFitter, tracks, and SequentialVertexFitter< N >::vertex().

Referenced by Tau.Tau::dxy().

94  {
95  return theSequentialFitter->vertex(tracks, priorPos, priorError);
96  }
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
SequentialVertexFitter< 5 > * theSequentialFitter
auto const & tracks
cannot be loose

Member Data Documentation

◆ theSequentialFitter

SequentialVertexFitter<5>* GsfVertexFitter::theSequentialFitter
private

Definition at line 99 of file GsfVertexFitter.h.

Referenced by GsfVertexFitter(), vertex(), and ~GsfVertexFitter().