Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef SequentialVertexFitter_H
0002 #define SequentialVertexFitter_H
0003 
0004 #include "RecoVertex/VertexPrimitives/interface/VertexFitter.h"
0005 #include "RecoVertex/VertexTools/interface/LinearizationPointFinder.h"
0006 #include "RecoVertex/VertexPrimitives/interface/VertexUpdator.h"
0007 #include "RecoVertex/VertexPrimitives/interface/VertexSmoother.h"
0008 #include "RecoVertex/VertexTools/interface/LinearizedTrackStateFactory.h"
0009 #include "RecoVertex/VertexTools/interface/VertexTrackFactory.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Utilities/interface/isFinite.h"
0012 #include <cmath>
0013 // #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0014 // #include "Vertex/VertexPrimitives/interface/VertexSeedFactory.h"
0015 
0016 /**
0017  * Sequential vertex fitter, where the vertex is updated one track
0018  * at the time.
0019  * The fitter will iterate over the set of tracks until the transverse
0020  * distance between vertices computed in the previous and the current
0021  * iterations is less than the specified convergence criteria, or until 
0022  * the maximum number of iterations is reached. 
0023  * The transverse distance determines the linearization error. 
0024  * The default convergence criterion is 1 mm. The default maximum 
0025  * number of steps is 10. 
0026  * These parameters can be configured in .orcarc (
0027  * SequentialVertexFitter:maximumDistance and 
0028  * SequentialVertexFitter:maximumNumberOfIterations). 
0029  * After the vertex fit, the tracks can be refit with the additional
0030  * constraint of the vertex position.
0031  */
0032 
0033 template <unsigned int N>
0034 class SequentialVertexFitter : public VertexFitter<N> {
0035 public:
0036   typedef ReferenceCountingPointer<RefittedTrackState<N> > RefCountedRefittedTrackState;
0037   typedef ReferenceCountingPointer<VertexTrack<N> > RefCountedVertexTrack;
0038   typedef ReferenceCountingPointer<LinearizedTrackState<N> > RefCountedLinearizedTrackState;
0039 
0040   /**
0041    *   Reimplemented constructors to use any kind of
0042    *   linearisation point finder, vertex updator and smoother.
0043    *   If no smoother is to be used, do not specify an instance for it.
0044    */
0045 
0046   SequentialVertexFitter(const LinearizationPointFinder& linP,
0047                          const VertexUpdator<N>& updator,
0048                          const VertexSmoother<N>& smoother,
0049                          const AbstractLTSFactory<N>& ltsf);
0050 
0051   /**
0052    *   Same as above, using a ParameterSet to set the convergence criteria
0053    */
0054 
0055   SequentialVertexFitter(const edm::ParameterSet& pSet,
0056                          const LinearizationPointFinder& linP,
0057                          const VertexUpdator<N>& updator,
0058                          const VertexSmoother<N>& smoother,
0059                          const AbstractLTSFactory<N>& ltsf);
0060 
0061   /**
0062    * Copy constructor
0063    */
0064 
0065   SequentialVertexFitter(const SequentialVertexFitter& original);
0066 
0067   ~SequentialVertexFitter() override;
0068 
0069   /**
0070    *  Method to set the convergence criterion 
0071    *  (the maximum distance between the vertex computed in the previous
0072    *   and the current iterations to consider the fit to have converged)
0073    */
0074 
0075   void setMaximumDistance(float maxShift) { theMaxShift = maxShift; }
0076 
0077   /**
0078    *   Method to set the maximum number of iterations to perform
0079    */
0080 
0081   void setMaximumNumberOfIterations(int maxIterations) { theMaxStep = maxIterations; }
0082 
0083   /**
0084   * Method returning the fitted vertex, from a container of reco::TransientTracks.
0085   * The linearization point will be searched with the given LP finder.
0086   * No prior vertex position will be used in the vertex fit.
0087   * \param tracks The container of RecTracks to fit.
0088   * \return The fitted vertex
0089   */
0090   CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks) const override;
0091 
0092   /**
0093   * Method returning the fitted vertex, from a container of VertexTracks.
0094   * For the first loop, the LinearizedTrackState contained in the VertexTracks
0095   * will be used. If subsequent loops are needed, the new VertexTracks will
0096   * be created with the last estimate of the vertex as linearization point.
0097   * No prior vertex position will be used in the vertex fit.
0098   * \param tracks The container of VertexTracks to fit.
0099   * \return The fitted vertex
0100   */
0101   CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack>& tracks) const override;
0102 
0103   /**
0104    * Same as above, only now also with BeamSpot!
0105    */
0106   CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack>& tracks, const reco::BeamSpot& spot) const override;
0107 
0108   /** Fit vertex out of a set of RecTracks. Uses the specified linearization point.
0109    */
0110   CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks, const GlobalPoint& linPoint) const override;
0111 
0112   /** Fit vertex out of a set of TransientTracks. 
0113    *  The specified BeamSpot will be used as priot, but NOT for the linearization.
0114    * The specified LinearizationPointFinder will be used to find the linearization point.
0115    */
0116   CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks,
0117                           const reco::BeamSpot& beamSpot) const override;
0118 
0119   /** Fit vertex out of a set of RecTracks. 
0120    *   Uses the position as both the linearization point AND as prior
0121    *   estimate of the vertex position. The error is used for the 
0122    *   weight of the prior estimate.
0123    */
0124   CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks,
0125                           const GlobalPoint& priorPos,
0126                           const GlobalError& priorError) const override;
0127 
0128   /** Fit vertex out of a set of VertexTracks
0129    *   Uses the position and error for the prior estimate of the vertex.
0130    *   This position is not used to relinearize the tracks.
0131    */
0132   CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack>& tracks,
0133                           const GlobalPoint& priorPos,
0134                           const GlobalError& priorError) const override;
0135 
0136   /**
0137   * Method returning the fitted vertex, from a VertexSeed.
0138   * For the first loop, the position of the VertexSeed will be used as
0139   * linearization point. If subsequent loops are needed, the new VertexTracks
0140   * will be created with the last estimate of the vertex as linearization point.
0141   * In case a non-sero error is given, the position and error of the
0142   * VertexSeed will be used as prior estimate in the vertex fit.
0143   * \param seed The VertexSeed to fit.
0144   * \return The fitted vertex
0145   */
0146   //  virtual CachingVertex<N> vertex(const RefCountedVertexSeed seed) const;
0147 
0148   /**
0149    * Access methods
0150    */
0151   const LinearizationPointFinder* linearizationPointFinder() const { return theLinP; }
0152 
0153   const VertexUpdator<N>* vertexUpdator() const { return theUpdator; }
0154 
0155   const VertexSmoother<N>* vertexSmoother() const { return theSmoother; }
0156 
0157   const float maxShift() const { return theMaxShift; }
0158 
0159   const int maxStep() const { return theMaxStep; }
0160 
0161   const edm::ParameterSet parameterSet() const { return thePSet; }
0162 
0163   SequentialVertexFitter* clone() const override { return new SequentialVertexFitter(*this); }
0164 
0165   const AbstractLTSFactory<N>* linearizedTrackStateFactory() const { return theLTrackFactory; }
0166 
0167   /**
0168   * Method checking whether a point is within the "tracker" bounds.
0169   * The default values are set to the CMS inner tracker and vertices
0170   * outsides these bounds will be rejected.
0171   * To reconstruct vertices within the full detector including the 
0172   * muon system, set the tracker bounds to larger values.
0173   */
0174   const bool insideTrackerBounds(const GlobalPoint& point) const {
0175     return ((point.transverse() < trackerBoundsRadius) && (abs(point.z()) < trackerBoundsHalfLength));
0176   }
0177 
0178   void setTrackerBounds(float radius, float halfLength) {
0179     trackerBoundsRadius = radius;
0180     trackerBoundsHalfLength = halfLength;
0181   }
0182 
0183 protected:
0184   /**
0185    *   Default constructor. Is here, as we do not want anybody to use it.
0186    */
0187 
0188   SequentialVertexFitter() {}
0189 
0190 private:
0191   /**
0192    * The methode where the vrte fit is actually done. The seed is used as the
0193    * prior estimate in the vertex fit (in case its error is large, it will have
0194    * little influence on the fit.
0195    * The tracks will be relinearized in case further loops are needed.
0196    *   \parameter tracks The tracks to use in the fit.
0197    *   \paraemter priorVertex The prior estimate of the vertex
0198    *   \return The fitted vertex
0199    */
0200   CachingVertex<N> fit(const std::vector<RefCountedVertexTrack>& tracks,
0201                        const VertexState priorVertex,
0202                        bool withPrior) const;
0203 
0204   /**
0205    * Construct a container of VertexTrack from a set of RecTracks.
0206    * \param tracks The container of RecTracks.
0207    * \param seed The seed to use for the VertexTracks. This position will
0208    *    also be used as the new linearization point.
0209    * \return The container of VertexTracks which are to be used in the next fit.
0210    */
0211   std::vector<RefCountedVertexTrack> linearizeTracks(const std::vector<reco::TransientTrack>& tracks,
0212                                                      const VertexState state) const;
0213 
0214   /**
0215    * Construct new a container of VertexTrack with a new linearization point
0216    * and vertex seed, from an existing set of VertexTrack, from which only the
0217    * recTracks will be used.
0218    * \param tracks The original container of VertexTracks, from which the RecTracks
0219    *    will be extracted.
0220    * \param seed The seed to use for the VertexTracks. This position will
0221    *    also be used as the new linearization point.
0222    * \return The container of VertexTracks which are to be used in the next fit.
0223    */
0224   std::vector<RefCountedVertexTrack> reLinearizeTracks(const std::vector<RefCountedVertexTrack>& tracks,
0225                                                        const VertexState state) const;
0226 
0227   /**
0228    *   Reads the configurable parameters.
0229    */
0230 
0231   void readParameters();
0232   void setDefaultParameters();
0233 
0234   /**
0235    * Checks whether any of the three coordinates is a Nan
0236    */
0237   inline bool hasNan(const GlobalPoint& point) const {
0238     using namespace std;
0239     return (edm::isNotFinite(point.x()) || edm::isNotFinite(point.y()) || edm::isNotFinite(point.z()));
0240   }
0241 
0242   float theMaxShift;
0243   int theMaxStep;
0244 
0245   // FIXME using hard-coded tracker bounds as default instead of taking them from geometry service
0246   float trackerBoundsRadius{112.};
0247   float trackerBoundsHalfLength{273.5};
0248 
0249   edm::ParameterSet thePSet;
0250   LinearizationPointFinder* theLinP;
0251   VertexUpdator<N>* theUpdator;
0252   VertexSmoother<N>* theSmoother;
0253   const AbstractLTSFactory<N>* theLTrackFactory;
0254   //   LinearizedTrackStateFactoryAbstractLTSFactory theLTrackFactory;
0255   VertexTrackFactory<N> theVTrackFactory;
0256 
0257   // VertexSeedFactory theVSeedFactory;
0258 };
0259 
0260 #endif