CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFDisplacedVertexFinder.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFTracking_PFDisplacedVertexFinder_h
2 #define RecoParticleFlow_PFTracking_PFDisplacedVertexFinder_h
3 
6 
13 
16 
18 
20 
21 
23 
28 class TrackingGeometry;
29 class TrackerGeometry;
30 class MagneticField;
31 
33 
34  public:
35 
37 
39 
41 
42  typedef std::set< reco::TrackBaseRef >::iterator IEset;
43  typedef reco::PFDisplacedVertexCandidateCollection::iterator IDVC;
44  typedef reco::PFDisplacedVertexSeedCollection::iterator IDVS;
45  typedef reco::PFDisplacedVertexCollection::iterator IDV;
46 
47  typedef std::pair <unsigned int, unsigned int> PFTrackHitInfo;
48  typedef std::pair <PFTrackHitInfo, PFTrackHitInfo> PFTrackHitFullInfo;
49 
51  enum FitterType {
56  };
57 
58 
60 
62  void setParameters(double transvSize, double longSize,
63  double primaryVertexCut, double tobCut,
64  double tecCut, double minAdaptWeight) {
65  transvSize_ = transvSize;
66  longSize_ = longSize;
67  primaryVertexCut_ = primaryVertexCut;
68  tobCut_ = tobCut;
69  tecCut_ = tecCut;
70  minAdaptWeight_ = minAdaptWeight;
71  }
72 
74  void setDebug( bool debug ) {debug_ = debug;}
75 
77  void setEdmParameters( const MagneticField* magField,
79  edm::ESHandle<TrackerGeometry> tkerGeomHandle){
80  magField_ = magField;
81  globTkGeomHandle_ = globTkGeomHandle;
82  tkerGeomHandle_ = tkerGeomHandle;
83  }
84 
87  }
88 
91  }
92 
94  edm::Handle< reco::BeamSpot > beamSpotHandle){
95  helper_.setPrimaryVertex(mainVertexHandle, beamSpotHandle);
96  }
97 
99  sigmacut_ = ps.getParameter<double>("sigmacut");
100  t_ini_ = ps.getParameter<double>("Tini");
101  ratio_ = ps.getParameter<double>("ratio");
102  }
103 
106 
107 
109  std::auto_ptr< reco::PFDisplacedVertexCollection > transferDisplacedVertices() {return displacedVertices_;}
110 
111  const std::auto_ptr< reco::PFDisplacedVertexCollection >& displacedVertices() const {return displacedVertices_;}
112 
113 
114 
116 
117  void findDisplacedVertices();
118 
119 
120  private:
121 
123 
126 
128  void mergeSeeds(reco::PFDisplacedVertexSeedCollection&, std::vector<bool>& bLocked);
129 
132 
134  void selectAndLabelVertices(reco::PFDisplacedVertexCollection&, std::vector <bool>&);
135 
137 
139 
141 
142  double getTransvDiff(const GlobalPoint&, const GlobalPoint&) const;
143  double getLongDiff(const GlobalPoint&, const GlobalPoint&) const;
144  double getLongProj(const GlobalPoint&, const GlobalVector&) const;
145 
147 
148  unsigned commonTracks(const reco::PFDisplacedVertex&, const reco::PFDisplacedVertex&) const;
149 
150  friend std::ostream& operator<<(std::ostream&, const PFDisplacedVertexFinder&);
151 
152 
154 
155  std::auto_ptr< reco::PFDisplacedVertexCandidateCollection > displacedVertexCandidates_;
156  std::auto_ptr< reco::PFDisplacedVertexCollection > displacedVertices_;
157 
159 
161 
162  double transvSize_;
163  double longSize_;
165  double tobCut_;
166  double tecCut_;
168 
170 
171  double sigmacut_; //= 6;
172  double t_ini_; //= 256.;
173  double ratio_; //= 0.25;
174 
175 
177  bool debug_;
178 
181 
184 
187 
188 
190 
192 
193 };
194 
195 #endif
196 
197 
double getTransvDiff(const GlobalPoint &, const GlobalPoint &) const
T getParameter(std::string const &) const
void setVertexIdentifier(const edm::ParameterSet &ps)
Set Vertex identifier parameters.
std::vector< PFDisplacedVertex > PFDisplacedVertexCollection
collection of PFDisplacedVertex objects
void setTracksSelector(const edm::ParameterSet &ps)
Set Tracks selector parameters.
std::pair< unsigned int, unsigned int > PFTrackHitInfo
A block of tracks linked together.
reco::PFDisplacedVertexSeedCollection::iterator IDVS
void setEdmParameters(const MagneticField *magField, edm::ESHandle< GlobalTrackingGeometry > globTkGeomHandle, edm::ESHandle< TrackerGeometry > tkerGeomHandle)
Sets parameters for track extrapolation and hits study.
std::auto_ptr< reco::PFDisplacedVertexCollection > displacedVertices_
edm::ESHandle< GlobalTrackingGeometry > globTkGeomHandle_
Tracker geometry for discerning hit positions.
bool rejectAndLabelVertex(reco::PFDisplacedVertex &dv)
void setPrimaryVertex(edm::Handle< reco::VertexCollection > mainVertexHandle, edm::Handle< reco::BeamSpot > beamSpotHandle)
Update the primary vertex information.
void findDisplacedVertices()
-----— Main function which find vertices -----— ///
bool debug_
If true, debug printouts activated.
std::pair< PFTrackHitInfo, PFTrackHitInfo > PFTrackHitFullInfo
bool fitVertexFromSeed(reco::PFDisplacedVertexSeed &, reco::PFDisplacedVertex &)
Fit one by one the vertex points with associated tracks to get displaced vertices.
unsigned commonTracks(const reco::PFDisplacedVertex &, const reco::PFDisplacedVertex &) const
void mergeSeeds(reco::PFDisplacedVertexSeedCollection &, std::vector< bool > &bLocked)
Sometimes two vertex candidates can be quite close and coming from the same vertex.
std::set< reco::TrackBaseRef >::iterator IEset
-----— Useful Types -----— ///
edm::ESHandle< TrackerGeometry > tkerGeomHandle_
doc?
std::auto_ptr< reco::PFDisplacedVertexCandidateCollection > displacedVertexCandidates_
-----— Members -----— ///
std::vector< PFDisplacedVertexSeed > PFDisplacedVertexSeedCollection
collection of PFDisplacedVertexSeed objects
void setDebug(bool debug)
Sets debug printout flag.
Displaced Vertex Finder Algorithm.
void selectAndLabelVertices(reco::PFDisplacedVertexCollection &, std::vector< bool > &)
Remove potentially fakes displaced vertices.
void setVertexIdentifier(edm::ParameterSet ps)
std::auto_ptr< reco::PFDisplacedVertexCollection > transferDisplacedVertices()
void setAVFParameters(edm::ParameterSet ps)
const MagneticField * magField_
to be able to extrapolate tracks f
void setTracksSelector(edm::ParameterSet ps)
friend std::ostream & operator<<(std::ostream &, const PFDisplacedVertexFinder &)
double getLongProj(const GlobalPoint &, const GlobalVector &) const
reco::PFDisplacedVertexCandidateCollection::iterator IDVC
void setInput(const edm::Handle< reco::PFDisplacedVertexCandidateCollection > &)
Set input collections of tracks.
const std::auto_ptr< reco::PFDisplacedVertexCollection > & displacedVertices() const
void findSeedsFromCandidate(reco::PFDisplacedVertexCandidate &, reco::PFDisplacedVertexSeedCollection &)
--—— Different steps of the finder algorithm --—— ///
double transvSize_
--—— Parameters --—— ///
Block of elements.
void setPrimaryVertex(edm::Handle< reco::VertexCollection > mainVertexHandle, edm::Handle< reco::BeamSpot > beamSpotHandle)
double sigmacut_
Adaptive Vertex Fitter parameters.
reco::PFDisplacedVertexCollection::iterator IDV
double getLongDiff(const GlobalPoint &, const GlobalPoint &) const
bool isCloseTo(const reco::PFDisplacedVertexSeed &, const reco::PFDisplacedVertexSeed &) const
-----— Tools -----— ///
#define debug
Definition: MEtoEDMFormat.h:34
void setParameters(double transvSize, double longSize, double primaryVertexCut, double tobCut, double tecCut, double minAdaptWeight)
--—— Set different algo parameters --—— ///
reco::PFDisplacedVertex::VertexTrackType getVertexTrackType(PFTrackHitFullInfo &) const
PFDisplacedVertexHelper helper_
PFCheckHitPatter.