Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef _AdaptiveVertexReconstructor_H_
0002 #define _AdaptiveVertexReconstructor_H_
0003 
0004 #include "RecoVertex/VertexPrimitives/interface/VertexReconstructor.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
0007 #include <set>
0008 
0009 class AdaptiveVertexReconstructor : public VertexReconstructor {
0010 public:
0011   /***
0012    *
0013    * \paramname primcut sigma_cut for the first iteration
0014    *   (primary vertex)
0015    * \paramname seccut sigma_cut for all subsequent vertex fits.
0016    * \paramname minweight the minimum weight for a track to
0017    * stay in a fitted vertex
0018    * \paramname smoothing perform track smoothing?
0019    */
0020   AdaptiveVertexReconstructor(float primcut = 2.0, float seccut = 6.0, float minweight = 0.5, bool smoothing = false);
0021 
0022   ~AdaptiveVertexReconstructor() override;
0023 
0024   /**
0025    *  The ParameterSet should have the following defined:
0026    *  double primcut
0027    *  double seccut
0028    *  double minweight
0029    *  for descriptions see 
0030    */
0031   AdaptiveVertexReconstructor(const edm::ParameterSet &s);
0032 
0033   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack> &v) const override;
0034 
0035   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack> &,
0036                                         const reco::BeamSpot &) const override;
0037 
0038   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack> &primaries,
0039                                         const std::vector<reco::TransientTrack> &tracks,
0040                                         const reco::BeamSpot &) const override;
0041 
0042   AdaptiveVertexReconstructor *clone() const override { return new AdaptiveVertexReconstructor(*this); }
0043 
0044 private:
0045   /**
0046    *  the actual fit to avoid code duplication
0047    */
0048   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack> &primaries,
0049                                         const std::vector<reco::TransientTrack> &trks,
0050                                         const reco::BeamSpot &,
0051                                         bool has_primaries,
0052                                         bool usespot) const;
0053 
0054   /**
0055    *  contrary to what its name has you believe, ::erase removes all
0056    *  newvtx.originalTracks() above theMinWeight from remainingtrks.
0057    */
0058   void erase(const TransientVertex &newvtx, std::set<reco::TransientTrack> &remainingtrks, float w) const;
0059 
0060   /**
0061    *  cleanup reconstructed vertices. discard all with too few significant
0062    *  tracks.
0063    */
0064   std::vector<TransientVertex> cleanUpVertices(const std::vector<TransientVertex> &) const;
0065 
0066   /** setup the vertex fitters.
0067    */
0068   void setupFitters(float primcut, float primT, float primr, float seccut, float secT, float secr, bool smoothing);
0069 
0070   TransientVertex cleanUp(const TransientVertex &old) const;
0071 
0072 private:
0073   AdaptiveVertexFitter *thePrimaryFitter;    // one fitter for the primaries ...
0074   AdaptiveVertexFitter *theSecondaryFitter;  // ... and one for the rest.
0075 
0076   // the minimum weight for a track to stay in a vertex.
0077   float theMinWeight;
0078   // the minimum weight for a track to be considered "significant".
0079   float theWeightThreshold;
0080 };
0081 
0082 #endif