Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef AdaptiveVertexFitter_H
0002 #define AdaptiveVertexFitter_H
0003 
0004 #include "RecoVertex/VertexPrimitives/interface/VertexFitter.h"
0005 #include "RecoVertex/LinearizationPointFinders/interface/DefaultLinearizationPointFinder.h"
0006 #include "RecoVertex/VertexTools/interface/GeometricAnnealing.h"
0007 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexUpdator.h"
0008 #include "RecoVertex/VertexTools/interface/DummyVertexSmoother.h"
0009 // #include "RecoVertex/AdaptiveVertexFit/interface/KalmanVertexSmoother.h"
0010 #include "RecoVertex/VertexPrimitives/interface/VertexState.h"
0011 #include "RecoVertex/VertexPrimitives/interface/VertexTrackCompatibilityEstimator.h"
0012 #include "RecoVertex/VertexTools/interface/LinearizedTrackStateFactory.h"
0013 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexTrackCompatibilityEstimator.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 
0016 /**
0017  * \class AdaptiveVertexFitter
0018  * An iterative reweighted fitter.
0019  * Very robust, very adaptive.
0020  *
0021  * See CMS Note 2007/008.
0022  *
0023  * Exceptions
0024  * VertexException( "Supplied fewer than two tracks" )
0025  * VertexException( "fewer than 2 significant tracks (w>threshold)" )
0026  *
0027  */
0028 
0029 class AdaptiveVertexFitter : public VertexFitter<5> {
0030 public:
0031   typedef ReferenceCountingPointer<VertexTrack<5> > RefCountedVertexTrack;
0032   typedef ReferenceCountingPointer<LinearizedTrackState<5> > RefCountedLinearizedTrackState;
0033 
0034   /**
0035    *   Reimplemented constructors to use any kind of
0036    *   linearisation point finder, vertex updator and smoother.
0037    *   If no smoother is to be used, do not specify an instance for it.
0038    */
0039   AdaptiveVertexFitter(const AnnealingSchedule &ann = GeometricAnnealing(),
0040                        const LinearizationPointFinder &linP = DefaultLinearizationPointFinder(),
0041                        const VertexUpdator<5> &updator = KalmanVertexUpdator<5>(),
0042                        const VertexTrackCompatibilityEstimator<5> &estor = KalmanVertexTrackCompatibilityEstimator<5>(),
0043                        const VertexSmoother<5> &smoother = DummyVertexSmoother<5>(),
0044                        const AbstractLTSFactory<5> &ltsf = LinearizedTrackStateFactory());
0045 
0046   AdaptiveVertexFitter(const AdaptiveVertexFitter &original);
0047 
0048   ~AdaptiveVertexFitter() override;
0049 
0050   /**
0051   * Method returning the fitted vertex, from a container of reco::TransientTracks.
0052   * The linearization point will be searched with the given LP finder.
0053   * No prior vertex position will be used in the vertex fit.
0054   * \return The fitted vertex
0055   */
0056   CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &) const override;
0057 
0058   /**
0059   * Method returning the fitted vertex, from a container of VertexTracks.
0060   * For the first loop, the LinearizedTrack contained in the VertexTracks
0061   * will be used. If subsequent loops are needed, the new VertexTracks will
0062   * be created with the last estimate of the vertex as linearization point.
0063   * No prior vertex position will be used in the vertex fit.
0064   * \return The fitted vertex
0065   */
0066   CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> &) const override;
0067 
0068   /**
0069    *  Same as above, only now also with BeamSpot constraint.
0070    */
0071   CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> &, const reco::BeamSpot &spot) const override;
0072 
0073   /** Fit vertex out of a std::vector of reco::TransientTracks. Uses the specified
0074    * linearization point.
0075    */
0076   CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &, const GlobalPoint &linPoint) const override;
0077 
0078   /** Fit vertex out of a set of reco::TransientTracks.
0079    *   Uses the position as both the linearization point AND as prior
0080    *   estimate of the vertex position. The error is used for the
0081    *   weight of the prior estimate.
0082    */
0083   CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &,
0084                           const GlobalPoint &priorPos,
0085                           const GlobalError &priorError) const override;
0086 
0087   /** Fit vertex out of a set of TransientTracks. 
0088    *  The specified BeamSpot will be used as priot, but NOT for the linearization.
0089    * The specified LinearizationPointFinder will be used to find the linearization point.
0090    */
0091   CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &tracks,
0092                           const reco::BeamSpot &beamSpot) const override;
0093 
0094   /**  Fit vertex out of a set of VertexTracks
0095    *   Uses the position and error for the prior estimate of the vertex.
0096    *   This position is not used to relinearize the tracks.
0097    */
0098   CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> &,
0099                           const GlobalPoint &priorPos,
0100                           const GlobalError &priorError) const override;
0101 
0102   AdaptiveVertexFitter *clone() const override;
0103 
0104   /**
0105    *  Set the weight threshold
0106    *  should be used only to find (once)
0107    *  a good value
0108    *  FIXME this should disappear in the final version
0109    */
0110   void setWeightThreshold(float w);
0111 
0112   /**
0113    *   Reads the configurable parameters.
0114    *   \param maxshift if the vertex moves further than this (in cm),
0115    *   then we re-iterate.
0116    *   \param maxlpshift if the vertex moves further than this,
0117    *   then we re-linearize the tracks.
0118    *   \param maxstep that's the maximum of iterations that we allow for.
0119    *   \param weightthreshold that's the minimum track weight
0120    *   for a track to be considered "significant".
0121    *   If fewer than two tracks are significant, an exception is thrown.
0122    */
0123   void setParameters(double maxshift = 0.0001,
0124                      double maxlpshift = 0.1,
0125                      unsigned maxstep = 30,
0126                      double weightthreshold = .001);
0127 
0128   /**
0129    *  Sets parameters.
0130    *  The following parameters are expected:
0131    *  maxshift,  maxlpshift,  maxstep,  weightthreshold
0132    */
0133   void setParameters(const edm::ParameterSet &);
0134 
0135   void gsfIntermediarySmoothing(bool sm) { gsfIntermediarySmoothing_ = sm; }
0136 
0137   bool gsfIntermediarySmoothing() const { return gsfIntermediarySmoothing_; }
0138 
0139 private:
0140   /**
0141    * Construct new a container of VertexTrack with a new linearization point
0142    * and vertex seed, from an existing set of VertexTrack, from which only the
0143    * recTracks will be used.
0144    * \param tracks The original container of VertexTracks, from which the reco::TransientTracks
0145    *     will be extracted.
0146    * \param vertex The seed to use for the VertexTracks. This position will
0147    *    also be used as the new linearization point.
0148    * \return The container of VertexTracks which are to be used in the next fit.
0149    */
0150   std::vector<RefCountedVertexTrack> reLinearizeTracks(const std::vector<RefCountedVertexTrack> &tracks,
0151                                                        const CachingVertex<5> &vertex) const;
0152 
0153   /**
0154    * Construct a new container of VertexTracks with new weights 
0155    * accounting for vertex error, from an existing set of LinearizedTracks. 
0156    */
0157   std::vector<RefCountedVertexTrack> reWeightTracks(const std::vector<RefCountedLinearizedTrackState> &,
0158                                                     const CachingVertex<5> &seed) const;
0159 
0160   /**
0161    * Construct new a container of VertexTracks with new weights 
0162    * accounting for vertex error, from an existing set of VertexTracks. 
0163    * From these the LinearizedTracks will be reused.
0164    */
0165   std::vector<RefCountedVertexTrack> reWeightTracks(const std::vector<RefCountedVertexTrack> &,
0166                                                     const CachingVertex<5> &seed) const;
0167 
0168   /**
0169    *  Weight the tracks, for the first time, using
0170    *  KalmanChiSquare.
0171    */
0172 
0173   std::vector<RefCountedVertexTrack> weightTracks(const std::vector<RefCountedLinearizedTrackState> &,
0174                                                   const VertexState &seed) const;
0175 
0176   /**
0177    *  Linearize tracks, for the first time in the iteration.
0178    */
0179   std::vector<RefCountedVertexTrack> linearizeTracks(const std::vector<reco::TransientTrack> &,
0180                                                      const VertexState &) const;
0181   /**
0182    *  perform the fit
0183    */
0184   CachingVertex<5> fit(const std::vector<RefCountedVertexTrack> &tracks,
0185                        const VertexState &priorSeed,
0186                        bool withPrior) const;
0187 
0188   double getWeight(float chi2) const;
0189 
0190 private:
0191   double theMaxShift;
0192   double theMaxLPShift;
0193   int theMaxStep;
0194   double theWeightThreshold;
0195   mutable int theNr;
0196 
0197   LinearizationPointFinder *theLinP;
0198   VertexUpdator<5> *theUpdator;
0199   VertexSmoother<5> *theSmoother;
0200   AnnealingSchedule *theAssProbComputer;
0201   VertexTrackCompatibilityEstimator<5> *theComp;
0202   const AbstractLTSFactory<5> *theLinTrkFactory;
0203   bool gsfIntermediarySmoothing_;
0204 };
0205 
0206 #endif