CMS 3D CMS Logo

Public Types | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes

SequentialVertexFitter< N > Class Template Reference

#include <SequentialVertexFitter.h>

Inheritance diagram for SequentialVertexFitter< N >:
VertexFitter< N >

List of all members.

Public Types

typedef
ReferenceCountingPointer
< LinearizedTrackState< N > > 
RefCountedLinearizedTrackState
typedef
ReferenceCountingPointer
< RefittedTrackState< N > > 
RefCountedRefittedTrackState
typedef
ReferenceCountingPointer
< VertexTrack< N > > 
RefCountedVertexTrack

Public Member Functions

SequentialVertexFitterclone () const
const LinearizationPointFinderlinearizationPointFinder () const
const AbstractLTSFactory< N > * linearizedTrackStateFactory () const
const float maxShift () const
const int maxStep () const
const edm::ParameterSet parameterSet () const
 SequentialVertexFitter (const edm::ParameterSet &pSet, const LinearizationPointFinder &linP, const VertexUpdator< N > &updator, const VertexSmoother< N > &smoother, const AbstractLTSFactory< N > &ltsf)
 SequentialVertexFitter (const LinearizationPointFinder &linP, const VertexUpdator< N > &updator, const VertexSmoother< N > &smoother, const AbstractLTSFactory< N > &ltsf)
 SequentialVertexFitter (const SequentialVertexFitter &original)
void setMaximumDistance (float maxShift)
void setMaximumNumberOfIterations (int maxIterations)
virtual CachingVertex< N > vertex (const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &beamSpot) const
virtual CachingVertex< N > vertex (const std::vector< RefCountedVertexTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const
virtual CachingVertex< N > vertex (const std::vector< RefCountedVertexTrack > &tracks) const
virtual CachingVertex< N > vertex (const std::vector< RefCountedVertexTrack > &tracks, const reco::BeamSpot &spot) const
virtual CachingVertex< N > vertex (const std::vector< reco::TransientTrack > &tracks) const
virtual CachingVertex< N > vertex (const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &linPoint) const
virtual CachingVertex< N > vertex (const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const
const VertexSmoother< N > * vertexSmoother () const
const VertexUpdator< N > * vertexUpdator () const
virtual ~SequentialVertexFitter ()

Protected Member Functions

 SequentialVertexFitter ()

Private Member Functions

CachingVertex< N > fit (const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
bool hasNan (const GlobalPoint &point) const
std::vector
< RefCountedVertexTrack
linearizeTracks (const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
void readParameters ()
std::vector
< RefCountedVertexTrack
reLinearizeTracks (const std::vector< RefCountedVertexTrack > &tracks, const VertexState state) const
void setDefaultParameters ()

Private Attributes

LinearizationPointFindertheLinP
const AbstractLTSFactory< N > * theLTrackFactory
float theMaxShift
int theMaxStep
edm::ParameterSet thePSet
VertexSmoother< N > * theSmoother
VertexUpdator< N > * theUpdator
VertexTrackFactory< N > theVTrackFactory

Detailed Description

template<unsigned int N>
class SequentialVertexFitter< N >

Sequential vertex fitter, where the vertex is updated one track at the time. The fitter will iterate over the set of tracks until the transverse distance between vertices computed in the previous and the current iterations is less than the specified convergence criteria, or until the maximum number of iterations is reached. The transverse distance determines the linearization error. The default convergence criterion is 1 mm. The default maximum number of steps is 10. These parameters can be configured in .orcarc ( SequentialVertexFitter:maximumDistance and SequentialVertexFitter:maximumNumberOfIterations). After the vertex fit, the tracks can be refit with the additional constraint of the vertex position.

Definition at line 34 of file SequentialVertexFitter.h.


Member Typedef Documentation

Definition at line 40 of file SequentialVertexFitter.h.

Definition at line 38 of file SequentialVertexFitter.h.

Definition at line 39 of file SequentialVertexFitter.h.


Constructor & Destructor Documentation

template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( const LinearizationPointFinder linP,
const VertexUpdator< N > &  updator,
const VertexSmoother< N > &  smoother,
const AbstractLTSFactory< N > &  ltsf 
)

Reimplemented constructors to use any kind of linearisation point finder, vertex updator and smoother. If no smoother is to be used, do not specify an instance for it.

Definition at line 23 of file SequentialVertexFitter.cc.

References SequentialVertexFitter< N >::setDefaultParameters().

                                       : 
  theLinP(linP.clone()), theUpdator(updator.clone()), 
  theSmoother(smoother.clone()), theLTrackFactory ( ltsf.clone() )
{
  setDefaultParameters();
}
template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( const edm::ParameterSet pSet,
const LinearizationPointFinder linP,
const VertexUpdator< N > &  updator,
const VertexSmoother< N > &  smoother,
const AbstractLTSFactory< N > &  ltsf 
)

Same as above, using a ParameterSet to set the convergence criteria

Definition at line 34 of file SequentialVertexFitter.cc.

References SequentialVertexFitter< N >::readParameters().

                                      : 
  thePSet(pSet), theLinP(linP.clone()), theUpdator(updator.clone()), 
  theSmoother(smoother.clone()), theLTrackFactory ( ltsf.clone() )
{
  readParameters();
}
template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( const SequentialVertexFitter< N > &  original)
template<unsigned int N>
SequentialVertexFitter< N >::~SequentialVertexFitter ( ) [virtual]

Definition at line 60 of file SequentialVertexFitter.cc.

{
  delete theLinP;
  delete theUpdator;
  delete theSmoother;
  delete theLTrackFactory;
}
template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( ) [inline, protected]

Default constructor. Is here, as we do not want anybody to use it.

Definition at line 192 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< 5 >::clone().

{}

Member Function Documentation

template<unsigned int N>
SequentialVertexFitter* SequentialVertexFitter< N >::clone ( ) const [inline, virtual]

Fit vertex out of a VertexSeed

Implements VertexFitter< N >.

Definition at line 179 of file SequentialVertexFitter.h.

Referenced by GsfVertexFitter::GsfVertexFitter().

                                         {
    return new SequentialVertexFitter(* this);
  }
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::fit ( const std::vector< RefCountedVertexTrack > &  tracks,
const VertexState  priorVertex,
bool  withPrior 
) const [private]

The methode where the vrte fit is actually done. The seed is used as the prior estimate in the vertex fit (in case its error is large, it will have little influence on the fit. The tracks will be relinearized in case further loops are needed. tracks The tracks to use in the fit. priorVertex The prior estimate of the vertex

Returns:
The fitted vertex

Definition at line 254 of file SequentialVertexFitter.cc.

References VertexState::error(), error, i, CachingVertex< N >::isValid(), LogDebug, VertexState::position(), CachingVertex< N >::position(), launcher::step, testEve_cfg::tracks, transverse(), and CachingVertex< N >::vertexState().

{
  std::vector<RefCountedVertexTrack> initialTracks;
  GlobalPoint priorVertexPosition = priorVertex.position();
  GlobalError priorVertexError = priorVertex.error();
  
  CachingVertex<N> returnVertex(priorVertexPosition,priorVertexError,initialTracks,0);
  if (withPrior) {
    returnVertex = CachingVertex<N>(priorVertexPosition,priorVertexError,
                priorVertexPosition,priorVertexError,initialTracks,0);
  }
  CachingVertex<N> initialVertex = returnVertex;
  std::vector<RefCountedVertexTrack> globalVTracks = tracks;

  // main loop through all the VTracks
  bool validVertex = true;
  int step = 0;
  GlobalPoint newPosition = priorVertexPosition;
  GlobalPoint previousPosition;
  do {
    CachingVertex<N> fVertex = initialVertex;
    // make new linearized and vertex tracks for the next iteration
    if(step != 0) globalVTracks = reLinearizeTracks(tracks, 
                                        returnVertex.vertexState());

    // update sequentially the vertex estimate
    for (typename std::vector<RefCountedVertexTrack>::const_iterator i 
           = globalVTracks.begin(); i != globalVTracks.end(); i++) {
      fVertex = theUpdator->add(fVertex,*i);
      if (!fVertex.isValid()) break;
    }

    validVertex = fVertex.isValid();
    // check tracker bounds and NaN in position
    if (validVertex && hasNan(fVertex.position())) {
      LogDebug("RecoVertex/SequentialVertexFitter") 
         << "Fitted position is NaN.\n";
      validVertex = false;
    }

    if (validVertex && !insideTrackerBounds(fVertex.position())) {
      LogDebug("RecoVertex/SequentialVertexFitter") 
         << "Fitted position is out of tracker bounds.\n";
      validVertex = false;
    }

    if (!validVertex) {
      // reset initial vertex position to (0,0,0) and force new iteration 
      // if number of steps not exceeded
      GlobalError error(AlgebraicSymMatrix(3,1)*10000);
      fVertex = CachingVertex<N>(GlobalPoint(0,0,0), error,
                              initialTracks, 0);
    }

    previousPosition = newPosition;
    newPosition = fVertex.position();

    returnVertex = fVertex;
    globalVTracks.clear();
    step++;
  } while ( (step != theMaxStep) &&
            (((previousPosition - newPosition).transverse() > theMaxShift) ||
                (!validVertex) ) );

  if (!validVertex) {
    LogDebug("RecoVertex/SequentialVertexFitter") 
       << "Fitted position is invalid (out of tracker bounds or has NaN). Returned vertex is invalid\n";
    return CachingVertex<N>(); // return invalid vertex
  }

  if (step >= theMaxStep) {
    LogDebug("RecoVertex/SequentialVertexFitter") 
       << "The maximum number of steps has been exceeded. Returned vertex is invalid\n";
    return CachingVertex<N>(); // return invalid vertex
  }

  // smoothing
  returnVertex = theSmoother->smooth(returnVertex);

  return returnVertex;
}
template<unsigned int N>
bool SequentialVertexFitter< N >::hasNan ( const GlobalPoint point) const [inline, private]

Checks whether any of the three coordinates is a Nan

Definition at line 243 of file SequentialVertexFitter.h.

                                                     {
    using namespace std;
    return (isnan(point.x())|| isnan(point.y()) || isnan(point.z()));
  }
template<unsigned int N>
const LinearizationPointFinder* SequentialVertexFitter< N >::linearizationPointFinder ( void  ) const [inline]

Method returning the fitted vertex, from a VertexSeed. For the first loop, the position of the VertexSeed will be used as linearization point. If subsequent loops are needed, the new VertexTracks will be created with the last estimate of the vertex as linearization point. In case a non-sero error is given, the position and error of the VertexSeed will be used as prior estimate in the vertex fit.

Parameters:
seedThe VertexSeed to fit.
Returns:
The fitted vertex Access methods

Definition at line 161 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

  {return theLinP;}
template<unsigned int N>
const AbstractLTSFactory<N>* SequentialVertexFitter< N >::linearizedTrackStateFactory ( ) const [inline]
template<unsigned int N>
std::vector<RefCountedVertexTrack> SequentialVertexFitter< N >::linearizeTracks ( const std::vector< reco::TransientTrack > &  tracks,
const VertexState  state 
) const [private]

Construct a container of VertexTrack from a set of RecTracks.

Parameters:
tracksThe container of RecTracks.
seedThe seed to use for the VertexTracks. This position will also be used as the new linearization point.
Returns:
The container of VertexTracks which are to be used in the next fit.
template<unsigned int N>
const float SequentialVertexFitter< N >::maxShift ( ) const [inline]
template<unsigned int N>
const int SequentialVertexFitter< N >::maxStep ( ) const [inline]
template<unsigned int N>
const edm::ParameterSet SequentialVertexFitter< N >::parameterSet ( ) const [inline]

Definition at line 176 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

  {return thePSet;}
template<unsigned int N>
void SequentialVertexFitter< N >::readParameters ( ) [private]

Reads the configurable parameters.

Definition at line 70 of file SequentialVertexFitter.cc.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

{
  theMaxShift = thePSet.getParameter<double>("maxDistance"); //0.01
  theMaxStep = thePSet.getParameter<int>("maxNbrOfIterations"); //10
}
template<unsigned int N>
vector< typename SequentialVertexFitter< N >::RefCountedVertexTrack > SequentialVertexFitter< N >::reLinearizeTracks ( const std::vector< RefCountedVertexTrack > &  tracks,
const VertexState  state 
) const [private]

Construct new a container of VertexTrack with a new linearization point and vertex seed, from an existing set of VertexTrack, from which only the recTracks will be used.

Parameters:
tracksThe original container of VertexTracks, from which the RecTracks will be extracted.
seedThe seed to use for the VertexTracks. This position will also be used as the new linearization point.
Returns:
The container of VertexTracks which are to be used in the next fit.

Definition at line 227 of file SequentialVertexFitter.cc.

References i, VertexState::position(), and VertexState::weight().

{

  GlobalPoint linP = state.position();
  std::vector<RefCountedVertexTrack> finalTracks;
  finalTracks.reserve(tracks.size());
  for(typename std::vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); 
      i != tracks.end(); i++) {
    RefCountedLinearizedTrackState lTrData = 
          (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
    //    RefCountedLinearizedTrackState lTrData = 
    //      theLTrackFactory->linearizedTrackState(linP, 
    //                              (**i).linearizedTrack()->track());
    RefCountedVertexTrack vTrData =
      theVTrackFactory.vertexTrack(lTrData,state, (**i).weight() );
    finalTracks.push_back(vTrData);
  }
  return finalTracks;
}
template<unsigned int N>
void SequentialVertexFitter< N >::setDefaultParameters ( ) [private]

Definition at line 77 of file SequentialVertexFitter.cc.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

{
  thePSet.addParameter<double>("maxDistance", 0.01);
  thePSet.addParameter<int>("maxNbrOfIterations", 10); //10
  readParameters();
}
template<unsigned int N>
void SequentialVertexFitter< N >::setMaximumDistance ( float  maxShift) [inline]

Method to set the convergence criterion (the maximum distance between the vertex computed in the previous and the current iterations to consider the fit to have converged)

Definition at line 78 of file SequentialVertexFitter.h.

template<unsigned int N>
void SequentialVertexFitter< N >::setMaximumNumberOfIterations ( int  maxIterations) [inline]

Method to set the maximum number of iterations to perform

Definition at line 85 of file SequentialVertexFitter.h.

        {theMaxStep = maxIterations;}
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< RefCountedVertexTrack > &  tracks,
const GlobalPoint priorPos,
const GlobalError priorError 
) const [virtual]

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

Definition at line 189 of file SequentialVertexFitter.cc.

References evf::utils::state.

{
  VertexState state(priorPos, priorError);
  return fit(tracks, state, true);
}
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< RefCountedVertexTrack > &  tracks,
const reco::BeamSpot spot 
) const [virtual]

Same as above, only now also with BeamSpot!

Definition at line 101 of file SequentialVertexFitter.cc.

References evf::utils::state.

{
  VertexState state(spot);
  return fit(tracks, state, true );
}
template<unsigned int N>
virtual CachingVertex<N> SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const GlobalPoint linPoint 
) const [virtual]

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

Implements VertexFitter< N >.

template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< RefCountedVertexTrack > &  tracks) const [virtual]

Method returning the fitted vertex, from a container of VertexTracks. For the first loop, the LinearizedTrackState contained in the VertexTracks will be used. If subsequent loops are needed, the new VertexTracks will be created with the last estimate of the vertex as linearization point. No prior vertex position will be used in the vertex fit.

Parameters:
tracksThe container of VertexTracks to fit.
Returns:
The fitted vertex

Definition at line 111 of file SequentialVertexFitter.cc.

References error, and evf::utils::state.

{
  // Initial vertex state, with a very small weight matrix
  GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
  AlgebraicSymMatrix we(3,1);
  GlobalError error(we*10000);
  VertexState state(linP, error);
  return fit(tracks, state, false);
}
template<unsigned int N>
virtual CachingVertex<N> SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks) const [virtual]

Method returning the fitted vertex, from a container of reco::TransientTracks. The linearization point will be searched with the given LP finder. No prior vertex position will be used in the vertex fit.

Parameters:
tracksThe container of RecTracks to fit.
Returns:
The fitted vertex

Implements VertexFitter< N >.

Referenced by KalmanVertexFitter::vertex(), and GsfVertexFitter::vertex().

template<unsigned int N>
virtual CachingVertex<N> SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const reco::BeamSpot beamSpot 
) const [virtual]

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< N >.

template<unsigned int N>
virtual CachingVertex<N> SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const GlobalPoint priorPos,
const GlobalError priorError 
) const [virtual]

Fit vertex out of a set of RecTracks. Uses the position 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< N >.

template<unsigned int N>
const VertexSmoother<N>* SequentialVertexFitter< N >::vertexSmoother ( ) const [inline]
template<unsigned int N>
const VertexUpdator<N>* SequentialVertexFitter< N >::vertexUpdator ( ) const [inline]

Member Data Documentation

template<unsigned int N>
LinearizationPointFinder* SequentialVertexFitter< N >::theLinP [private]
template<unsigned int N>
const AbstractLTSFactory<N>* SequentialVertexFitter< N >::theLTrackFactory [private]
template<unsigned int N>
float SequentialVertexFitter< N >::theMaxShift [private]
template<unsigned int N>
int SequentialVertexFitter< N >::theMaxStep [private]
template<unsigned int N>
edm::ParameterSet SequentialVertexFitter< N >::thePSet [private]

Definition at line 252 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< 5 >::parameterSet().

template<unsigned int N>
VertexSmoother<N>* SequentialVertexFitter< N >::theSmoother [private]
template<unsigned int N>
VertexUpdator<N>* SequentialVertexFitter< N >::theUpdator [private]
template<unsigned int N>
VertexTrackFactory<N> SequentialVertexFitter< N >::theVTrackFactory [private]

Definition at line 258 of file SequentialVertexFitter.h.