CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AdaptiveVertexReconstructor.h
Go to the documentation of this file.
1 #ifndef _AdaptiveVertexReconstructor_H_
2 #define _AdaptiveVertexReconstructor_H_
3 
7 #include <set>
8 
10 public:
11 
12  /***
13  *
14  * \paramname primcut sigma_cut for the first iteration
15  * (primary vertex)
16  * \paramname seccut sigma_cut for all subsequent vertex fits.
17  * \paramname minweight the minimum weight for a track to
18  * stay in a fitted vertex
19  * \paramname smoothing perform track smoothing?
20  */
21  AdaptiveVertexReconstructor( float primcut = 2.0, float seccut = 6.0,
22  float minweight = 0.5, bool smoothing = false );
23 
25 
34 
35  std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack> & v ) const;
36 
37  std::vector<TransientVertex>
38  vertices(const std::vector<reco::TransientTrack> &, const reco::BeamSpot & ) const;
39 
40  std::vector<TransientVertex>
41  vertices(const std::vector<reco::TransientTrack> & primaries,
42  const std::vector<reco::TransientTrack> & tracks,
43  const reco::BeamSpot & ) const;
44 
45  virtual AdaptiveVertexReconstructor * clone() const {
46  return new AdaptiveVertexReconstructor( * this );
47  }
48 
49 private:
53  std::vector<TransientVertex>
54  vertices( const std::vector<reco::TransientTrack> & primaries,
55  const std::vector<reco::TransientTrack> & trks,
56  const reco::BeamSpot &, bool has_primaries, bool usespot ) const;
57 
62  void erase ( const TransientVertex & newvtx,
63  std::set < reco::TransientTrack > & remainingtrks, float w ) const;
64 
69  std::vector<TransientVertex> cleanUpVertices (
70  const std::vector < TransientVertex > & ) const;
71 
74  void setupFitters ( float primcut, float primT, float primr,
75  float seccut, float secT, float secr, bool smoothing );
76 
77  TransientVertex cleanUp ( const TransientVertex & old ) const;
78 
79 private:
80  AdaptiveVertexFitter * thePrimaryFitter; // one fitter for the primaries ...
81  AdaptiveVertexFitter * theSecondaryFitter; // ... and one for the rest.
82 
83  // the minimum weight for a track to stay in a vertex.
84  float theMinWeight;
85  // the minimum weight for a track to be considered "significant".
87 };
88 
89 #endif
AdaptiveVertexReconstructor(float primcut=2.0, float seccut=6.0, float minweight=0.5, bool smoothing=false)
const double w
Definition: UKUtility.cc:23
TransientVertex cleanUp(const TransientVertex &old) const
virtual AdaptiveVertexReconstructor * clone() const
void setupFitters(float primcut, float primT, float primr, float seccut, float secT, float secr, bool smoothing)
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &v) const
std::vector< TransientVertex > cleanUpVertices(const std::vector< TransientVertex > &) const
AdaptiveVertexFitter * theSecondaryFitter
tuple tracks
Definition: testEve_cfg.py:39
void erase(const TransientVertex &newvtx, std::set< reco::TransientTrack > &remainingtrks, float w) const