Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:09

0001 #ifndef KalmanVertexFitter_H
0002 #define KalmanVertexFitter_H
0003 
0004 #include "RecoVertex/VertexTools/interface/SequentialVertexFitter.h"
0005 
0006 /** Least-squares vertex fitter implemented in the Kalman Filter formalism 
0007  *  Fits vertex position and, if smoothing is requested at construction, 
0008  *  constrains track parameters at the vertex. 
0009  *  A beam spot constraint can also be specified in form of a 3D point 
0010  *  and a 3D matrix describing a Gaussian beam spread. 
0011  *
0012  *  References: 
0013  *
0014  * <a href="http://cmsdoc.cern.ch/documents/03/in/in03_008.pdf"> `Vertex fitting with the Kalman filter formalism', T.speer, K.Prokofiev, R.Fruehwirth, CMS IN 2003/008</a><br>
0015  * <a href="http://cms.cern.ch/iCMS/jsp/openfile.jsp?type=NOTE&year=2006&files=NOTE2006_032.pdf"> T.Speer, K.Prokofiev, R.Fruhwirth, W.Waltenberger, P.Vanlaer, "Vertex Fitting in the CMS Tracker", CMS Note 2006/032</a><br>
0016  * P.Billoir, S.Qian, "Fast vertex fitting...", NIM A311 (1992) 139. <br>
0017  * R.Fruhwirth, R.Kubinec, W.Mitaroff, M.Regler, Comp. Physics Comm. 96, 189 (1996)<br>
0018  * <a href="http://www.phys.ufl.edu/~avery/fitting.html"> P.Avery, lectures 
0019 on track and vertex fitting" </a>
0020  */
0021 
0022 class KalmanVertexFitter : public VertexFitter<5> {
0023 public:
0024   /**
0025    * The constructor, setting everything up to have a VertexFitter using the 
0026    * Kalman algorithm.
0027    * \param useSmoothing Specifies whether the tracks should be refit or not.
0028    */
0029 
0030   KalmanVertexFitter(bool useSmoothing = false);
0031 
0032   /**
0033    * Same as above, using a ParameterSet to set the convergence criteria
0034    */
0035 
0036   KalmanVertexFitter(const edm::ParameterSet& pSet, bool useSmoothing = false);
0037 
0038   KalmanVertexFitter(const KalmanVertexFitter& other) : theSequentialFitter(other.theSequentialFitter->clone()) {}
0039 
0040   ~KalmanVertexFitter() override { delete theSequentialFitter; }
0041 
0042   KalmanVertexFitter* clone() const override { return new KalmanVertexFitter(*this); }
0043 
0044 public:
0045   typedef CachingVertex<5>::RefCountedVertexTrack RefCountedVertexTrack;
0046 
0047   /** Fit vertex out of a set of RecTracks
0048    */
0049   inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks) const override {
0050     return theSequentialFitter->vertex(tracks);
0051   }
0052 
0053   /** Fit vertex out of a set of VertexTracks
0054    */
0055   inline CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack>& tracks) const override {
0056     return theSequentialFitter->vertex(tracks);
0057   }
0058 
0059   inline CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack>& tracks,
0060                                  const reco::BeamSpot& spot) const override {
0061     return theSequentialFitter->vertex(tracks, spot);
0062   }
0063 
0064   /** Fit vertex out of a set of RecTracks. 
0065    *  Uses the specified linearization point.
0066    */
0067   inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks,
0068                                  const GlobalPoint& linPoint) const override {
0069     return theSequentialFitter->vertex(tracks, linPoint);
0070   }
0071 
0072   /** Fit vertex out of a set of RecTracks. 
0073    *  Uses the specified point as both the linearization point AND as prior
0074    *  estimate of the vertex position. The error is used for the 
0075    *  weight of the prior estimate.
0076    */
0077   inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks,
0078                                  const GlobalPoint& priorPos,
0079                                  const GlobalError& priorError) const override {
0080     return theSequentialFitter->vertex(tracks, priorPos, priorError);
0081   }
0082 
0083   /** Fit vertex out of a set of TransientTracks. 
0084    *  The specified BeamSpot will be used as priot, but NOT for the linearization.
0085    * The specified LinearizationPointFinder will be used to find the linearization point.
0086    */
0087   inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks,
0088                                  const reco::BeamSpot& beamSpot) const override {
0089     return theSequentialFitter->vertex(tracks, beamSpot);
0090   }
0091 
0092   /** Fit vertex out of a set of VertexTracks.
0093    *  Uses the specified point and error as the prior estimate of the vertex.
0094    *  This position is not used to relinearize the tracks.
0095    */
0096   inline CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack>& tracks,
0097                                  const GlobalPoint& priorPos,
0098                                  const GlobalError& priorError) const override {
0099     return theSequentialFitter->vertex(tracks, priorPos, priorError);
0100   }
0101 
0102   /** Default convergence criteria
0103    */
0104   //  edm::ParameterSet defaultParameters() const;
0105 
0106 private:
0107   void setup(const edm::ParameterSet& pSet, bool useSmoothing);
0108 
0109   edm::ParameterSet defaultParameters() const;
0110 
0111   const SequentialVertexFitter<5>* theSequentialFitter;
0112 };
0113 
0114 #endif