Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef _ConfigurableTrimmedVertexFinder_H_
0002 #define _ConfigurableTrimmedVertexFinder_H_
0003 
0004 #include "RecoVertex/VertexPrimitives/interface/VertexReconstructor.h"
0005 #include "RecoVertex/TrimmedKalmanVertexFinder/interface/TrimmedVertexFinder.h"
0006 #include "RecoVertex/TrimmedKalmanVertexFinder/interface/TrimmedTrackFilter.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include <vector>
0009 
0010 /** Algorithm to find a series of distinct vertices among the given set 
0011  *  of tracks. The technique is: <BR>
0012  *    1) use `TrimmedTrackFilter` to select tracks 
0013  *    above a certain pT; <BR>
0014  *    2) use `TrimmedVertexFinder` to split the set of tracks into 
0015  *    a cluster of compatible tracks and a set of remaining tracks; <BR>
0016  *    3) repeat 2) with the remaining set, and repeat as long as 
0017  *    (a) a cluster of compatible tracks can be found or (b) 
0018  *    the maximum number of clusters asked for is reached; <BR>
0019  *    4) reject vertices with a low fit probability. <BR>
0020  *
0021  *  This algorithm has 5 parameters that can be set at runtime 
0022  *  via the corresponding set() methods, or a ParameterSet: <BR>
0023  *   - "ptCut" (default: 1.5 GeV/c) 
0024  *  which defines the minimum pT of the tracks used to make vertices. 
0025  *  This value overrides the corresponding configurable parameter 
0026  *  of the TrimmedTrackFilter used internally. <BR>
0027  *   - "trackCompatibilityCut" (default: 0.05) 
0028  *  which defines the probability below which a track is considered 
0029  *  incompatible with the 1st vertex candidate formed. 
0030  *  This value overrides the corresponding configurable parameter 
0031  *  of the TrimmedVertexFinder used internally. <BR>
0032  *   - "trackCompatibilityToSV" (default: 0.01) 
0033  *  which defines the probability below which a track is considered 
0034  *  incompatible with the next vertex candidates formed. 
0035  *  This value overrides the corresponding configurable parameter 
0036  *  of the TrimmedVertexFinder used internally. <BR>
0037  *   - "vertexFitProbabilityCut" (default: 0.01) 
0038  *  which defines the probability below which a vertex is rejected. <BR>
0039  *   - "maxNbOfVertices" (default: 0) 
0040  *  which defines the maximum number of vertices searched for. 
0041  *  0 means all vertices which can be found. 
0042  */
0043 
0044 class ConfigurableTrimmedVertexFinder : public VertexReconstructor {
0045 public:
0046   ConfigurableTrimmedVertexFinder(const VertexFitter<5>* vf,
0047                                   const VertexUpdator<5>* vu,
0048                                   const VertexTrackCompatibilityEstimator<5>* ve);
0049 
0050   ~ConfigurableTrimmedVertexFinder() override {}
0051 
0052   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack>& tracks) const override;
0053 
0054   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack>& tracks,
0055                                         const reco::BeamSpot& spot) const override;
0056 
0057   std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack>& tracks,
0058                                         std::vector<reco::TransientTrack>& unused,
0059                                         const reco::BeamSpot& spot,
0060                                         bool use_spot) const;
0061 
0062   /** Access to parameters
0063    */
0064   float ptCut() const { return theFilter.ptCut(); }
0065   const TrimmedTrackFilter& trackFilter() const { return theFilter; }
0066   float trackCompatibilityCut() const { return theTrackCompatibilityToPV; }
0067   float trackCompatibilityToSV() const { return theTrackCompatibilityToSV; }
0068   float vertexFitProbabilityCut() const { return theVtxFitProbCut; }
0069   int maxNbOfVertices() const { return theMaxNbOfVertices; }
0070 
0071   /** Set parameters
0072    */
0073   void setParameters(const edm::ParameterSet&);
0074 
0075   void setPtCut(float cut) { theFilter.setPtCut(cut); }
0076   void setTrackCompatibilityCut(float cut) { theTrackCompatibilityToPV = cut; }
0077   void setTrackCompatibilityToSV(float cut) { theTrackCompatibilityToSV = cut; }
0078   void setVertexFitProbabilityCut(float cut) { theVtxFitProbCut = cut; }
0079   void setMaxNbOfVertices(int max) { theMaxNbOfVertices = max; }
0080 
0081   /** Clone method
0082    */
0083   ConfigurableTrimmedVertexFinder* clone() const override { return new ConfigurableTrimmedVertexFinder(*this); }
0084 
0085 protected:
0086   virtual void resetEvent(const std::vector<reco::TransientTrack>& tracks) const {}
0087 
0088   virtual void analyseInputTracks(const std::vector<reco::TransientTrack>& tracks) const {}
0089 
0090   virtual void analyseClusterFinder(const std::vector<TransientVertex>& vts,
0091                                     const std::vector<reco::TransientTrack>& remain) const {}
0092 
0093   virtual void analyseVertexCandidates(const std::vector<TransientVertex>& vts) const {}
0094 
0095   virtual void analyseFoundVertices(const std::vector<TransientVertex>& vts) const {}
0096 
0097 private:
0098   using VertexReconstructor::vertices;
0099 
0100   // separate tracks passing the filter
0101   //  void separateTracks(std::vector<TransientTrack>& filtered,
0102   //              std::vector<TransientTrack>& unused) const;
0103 
0104   // find vertex candidates
0105   std::vector<TransientVertex> vertexCandidates(const std::vector<reco::TransientTrack>& tracks,
0106                                                 std::vector<reco::TransientTrack>& unused,
0107                                                 const reco::BeamSpot& spot,
0108                                                 bool use_spot) const;
0109 
0110   // remove bad candidates
0111   std::vector<TransientVertex> clean(const std::vector<TransientVertex>& candidates) const;
0112 
0113   //
0114   mutable TrimmedVertexFinder theClusterFinder;
0115   float theVtxFitProbCut;
0116   float theTrackCompatibilityToPV;
0117   float theTrackCompatibilityToSV;
0118   int theMaxNbOfVertices;
0119   TrimmedTrackFilter theFilter;
0120 };
0121 
0122 #endif