CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PrimaryVertexProducerAlgorithm Class Reference

#include <RecoVertex/PrimaryVertexProducerAlgorithm/src/PrimaryVertexProducerAlgorithm.cc>

Inheritance diagram for PrimaryVertexProducerAlgorithm:
VertexReconstructor

Public Member Functions

virtual
PrimaryVertexProducerAlgorithm
clone () const
 
 PrimaryVertexProducerAlgorithm (const edm::ParameterSet &)
 
virtual std::vector
< TransientVertex
vertices (const std::vector< reco::TransientTrack > &tracks) const
 
virtual std::vector
< TransientVertex
vertices (const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &beamSpot) const
 
 ~PrimaryVertexProducerAlgorithm ()
 
- Public Member Functions inherited from VertexReconstructor
 VertexReconstructor ()
 
virtual std::vector
< TransientVertex
vertices (const std::vector< reco::TransientTrack > &primaries, const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &spot) const
 
virtual ~VertexReconstructor ()
 

Private Attributes

bool fapply_finder
 
bool fFailsafe
 
double fMinNdof
 
bool fUseBeamConstraint
 
bool fVerbose
 
edm::ParameterSet theConfig
 
KalmanTrimmedVertexFinder theFinder
 
VertexFitter< 5 > * theFitter
 
TrackClusterizerInZtheTrackClusterizer
 
TrackFilterForPVFindingBasetheTrackFilter
 
VertexCompatibleWithBeam theVertexSelector
 

Detailed Description

Description: finds primary vertices, compatible with the beam line

Implementation: <Notes on="" implementation>="">

Definition at line 35 of file PrimaryVertexProducerAlgorithm.h.

Constructor & Destructor Documentation

PrimaryVertexProducerAlgorithm::PrimaryVertexProducerAlgorithm ( const edm::ParameterSet conf)
explicit

Definition at line 21 of file PrimaryVertexProducerAlgorithm.cc.

References ExpressReco_HICollisions_FallBack::algorithm, fapply_finder, fFailsafe, fMinNdof, fUseBeamConstraint, fVerbose, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), KalmanTrimmedVertexFinder::setParameters(), theFinder, theFitter, theTrackClusterizer, and theTrackFilter.

Referenced by clone().

23  : theConfig(conf),
25  conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"))
26 {
27  edm::LogInfo("PVDebugInfo")
28  << "PVSelParameters::maxDistanceToBeam = "
29  << conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam") << "\n";
30 
31 
32  fUseBeamConstraint = conf.getParameter<bool>("useBeamConstraint");
33  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
34  fMinNdof = conf.getParameter<double>("minNdof");
35  fFailsafe = true; //conf.getUntrackedParameter<bool>("failsafe",true);
36 
37 
38  // select and configure the track selection
39  std::string trackSelectionAlgorithm=conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
40  if(trackSelectionAlgorithm=="filter"){
41  theTrackFilter= new TrackFilterForPVFinding( conf.getParameter<edm::ParameterSet>("TkFilterParameters") );
42  }else if (trackSelectionAlgorithm=="filterWithThreshold"){
44  }else{
45  throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " + trackSelectionAlgorithm);
46  }
47 
48 
49  // select and configure the track clusterizer
50  std::string clusteringAlgorithm=conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
51  if (clusteringAlgorithm=="gap"){
52  theTrackClusterizer = new GapClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
53  }else if(clusteringAlgorithm=="DA"){
54  theTrackClusterizer = new DAClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
55  }else{
56  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
57  }
58 
59  // select and configure the vertex fitter
60  std::string algorithm = conf.getParameter<std::string>("algorithm");
61  fapply_finder = false;
62  if (algorithm == "TrimmedKalmanFinder") {
63  fapply_finder = true;
64  theFinder.setParameters(conf.getParameter<edm::ParameterSet>("VtxFinderParameters"));
65  } else if (algorithm=="KalmanVertexFitter") {
67  } else if( algorithm=="AdaptiveVertexFitter") {
69  } else {
70  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + algorithm);
71  }
72 
73  edm::LogInfo("PVDebugInfo")
74  << "Using " << algorithm << "\n";
75  edm::LogInfo("PVDebugInfo")
76  << "beam-constraint " << fUseBeamConstraint << "\n";
77 
78  edm::LogInfo("PVDebugInfo")
79  << "PV producer algorithm initialization: done" << "\n";
80 
81 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
Common base class.
void setParameters(const edm::ParameterSet &)
TrackFilterForPVFindingBase * theTrackFilter
PrimaryVertexProducerAlgorithm::~PrimaryVertexProducerAlgorithm ( )

Definition at line 84 of file PrimaryVertexProducerAlgorithm.cc.

References theFitter, theTrackClusterizer, and theTrackFilter.

85 {
86  if (theFitter) delete theFitter;
87  if (theTrackFilter) delete theTrackFilter;
89 }
TrackFilterForPVFindingBase * theTrackFilter

Member Function Documentation

virtual PrimaryVertexProducerAlgorithm* PrimaryVertexProducerAlgorithm::clone ( void  ) const
inlinevirtual

Clone method

Implements VertexReconstructor.

Definition at line 54 of file PrimaryVertexProducerAlgorithm.h.

References PrimaryVertexProducerAlgorithm().

54  {
55  return new PrimaryVertexProducerAlgorithm(*this);
56  }
PrimaryVertexProducerAlgorithm(const edm::ParameterSet &)
std::vector< TransientVertex > PrimaryVertexProducerAlgorithm::vertices ( const std::vector< reco::TransientTrack > &  tracks) const
virtual

Find primary vertices

Implements VertexReconstructor.

Definition at line 98 of file PrimaryVertexProducerAlgorithm.cc.

Referenced by PrimaryVertexProducer::produce().

99 {
100 
101  throw VertexException("PrimaryVertexProducerAlgorithm: cannot make a Primary Vertex without a beam spot constraint " );
102 
103  /* std::cout<< "PrimaryVertexProducer::vertices> Obsolete function, using dummy beamspot " << std::endl;
104  reco::BeamSpot dummyBeamSpot;
105  dummyBeamSpot.dummy();
106  return vertices(tracks,dummyBeamSpot); */
107  return std::vector<TransientVertex>();
108 }
Common base class.
std::vector< TransientVertex > PrimaryVertexProducerAlgorithm::vertices ( const std::vector< reco::TransientTrack > &  t,
const reco::BeamSpot  
) const
virtual

Reconstruct vertices, exploiting the beamspot constraint for the primary vertex

Reimplemented from VertexReconstructor.

Definition at line 112 of file PrimaryVertexProducerAlgorithm.cc.

References TrackClusterizerInZ::clusterize(), gather_cfg::cout, GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::cyy(), GlobalErrorBase< T, ErrorWeightType >::czz(), TransientVertex::degreesOfFreedom(), VertexState::error(), fapply_finder, fFailsafe, fMinNdof, fUseBeamConstraint, fVerbose, TransientVertex::isValid(), GlobalErrorBase< T, ErrorWeightType >::matrix(), TransientVertex::position(), TrackFilterForPVFindingBase::select(), python.multivaluedict::sort(), theFinder, theFitter, theTrackClusterizer, theTrackFilter, theVertexSelector, v, VertexFitter< N >::vertex(), KalmanTrimmedVertexFinder::vertices(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

114 {
115  bool validBS = true;
116  VertexState beamVertexState(beamSpot);
117  if ( (beamVertexState.error().cxx() <= 0.) ||
118  (beamVertexState.error().cyy() <= 0.) ||
119  (beamVertexState.error().czz() <= 0.) ) {
120  validBS = false;
121  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors "<<beamVertexState.error().matrix();
122  }
123 
124  if ( fapply_finder) {
125  return theFinder.vertices( tracks );
126  }
127  std::vector<TransientVertex> pvs;
128 
129 
130  // select tracks
131  std::vector<TransientTrack> seltks = theTrackFilter->select( tracks );
132 
133 
134  // clusterize tracks in Z
135  std::vector< std::vector<reco::TransientTrack> > clusters = theTrackClusterizer->clusterize(seltks);
136  if (fVerbose){std::cout << " clustering returned "<< clusters.size() << " clusters from " << seltks.size() << " selected tracks" <<std::endl;}
137 
138 
139  // look for primary vertices in each cluster
140  std::vector<TransientVertex> pvCand;
141  int nclu=0;
142  for (std::vector< std::vector<reco::TransientTrack> >::const_iterator iclus
143  = clusters.begin(); iclus != clusters.end(); iclus++) {
144 
145 
147  if( fUseBeamConstraint && validBS &&((*iclus).size()>1) ){
148  if (fVerbose){std::cout << " constrained fit with "<< (*iclus).size() << " tracks" <<std::endl;}
149  v = theFitter->vertex(*iclus, beamSpot);
150  if (v.isValid() && (v.degreesOfFreedom()>=fMinNdof)) pvCand.push_back(v);
151 
152  if (fVerbose){
153  if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " << v.position().z() << std::endl;
154  else std::cout <<"Invalid fitted vertex\n";
155  }
156 
157  }else if((*iclus).size()>1){
158  if (fVerbose){std::cout << " unconstrained fit with "<< (*iclus).size() << " tracks" << std::endl;}
159 
160  v = theFitter->vertex(*iclus);
161  if (v.isValid() && (v.degreesOfFreedom()>=fMinNdof)) pvCand.push_back(v);
162 
163  if (fVerbose){
164  if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " << v.position().z() << std::endl;
165  else std::cout <<"Invalid fitted vertex\n";
166  }
167 
168  }
169 
170  nclu++;
171 
172  }// end of cluster loop
173 
174  if(fVerbose){
175  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvCand.size() << std::endl;
176  }
177 
178 
179 
180  // select vertices compatible with beam
181  int npv=0;
182  for (std::vector<TransientVertex>::const_iterator ipv = pvCand.begin();
183  ipv != pvCand.end(); ipv++) {
184  if(fVerbose){
185  std::cout << "PrimaryVertexProducerAlgorithm::vertices cand " << npv++ << " sel=" <<
186  (validBS && theVertexSelector(*ipv,beamVertexState)) << " z=" << ipv->position().z() << std::endl;
187  }
188  if (!validBS || theVertexSelector(*ipv,beamVertexState)) pvs.push_back(*ipv);
189  }
190 
191 
192  if(pvs.size()>0){
193 
194  // sort vertices by pt**2 vertex (aka signal vertex tagging)
195  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
196 
197  }else{
198 
199  if ( fFailsafe
200  && (seltks.size()>1)
201  && ( (clusters.size()!=1) || ( (clusters.size()==1) && (clusters.begin()->size()<seltks.size())) )
202  )
203  {
204  // if no vertex was found, try fitting all selected tracks, unless this has already been tried
205  // in low/no pile-up situations with low multiplicity vertices, this can recover vertices lost in clustering
206  // with very small zSep. Only makes sense when used with a robust fitter, like the AdaptiveVertexFitter
208  if( fUseBeamConstraint && validBS ){
209  v = theFitter->vertex(seltks, beamSpot);
210  }else{
211  v = theFitter->vertex(seltks);
212  }
213  if (fVerbose){ std::cout << "PrimaryVertexProducerAlgorithm: failsafe vertex "
214  <<" tracks=" << seltks.size()
215  <<" valid()=" << v.isValid() << " ndof()=" << v.degreesOfFreedom()
216  <<" selected="<< theVertexSelector(v,beamVertexState) << std::endl; }
217  if ( v.isValid() && (v.degreesOfFreedom()>=fMinNdof) && (theVertexSelector(v,beamVertexState)) ) pvs.push_back(v);
218  }
219 
220  }
221 
222 
223  return pvs;
224 
225 }
virtual std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const
T y() const
Definition: PV3DBase.h:57
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const =0
float degreesOfFreedom() const
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:58
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
TrackFilterForPVFindingBase * theTrackFilter
tuple tracks
Definition: testEve_cfg.py:39
tuple cout
Definition: gather_cfg.py:41
T x() const
Definition: PV3DBase.h:56
bool isValid() const
mathSSE::Vec4< T > v

Member Data Documentation

bool PrimaryVertexProducerAlgorithm::fapply_finder
private

Definition at line 71 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm(), and vertices().

bool PrimaryVertexProducerAlgorithm::fFailsafe
private

Definition at line 72 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm(), and vertices().

double PrimaryVertexProducerAlgorithm::fMinNdof
private

Definition at line 69 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm(), and vertices().

bool PrimaryVertexProducerAlgorithm::fUseBeamConstraint
private

Definition at line 68 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm(), and vertices().

bool PrimaryVertexProducerAlgorithm::fVerbose
private

Definition at line 67 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm(), and vertices().

edm::ParameterSet PrimaryVertexProducerAlgorithm::theConfig
private

Definition at line 61 of file PrimaryVertexProducerAlgorithm.h.

KalmanTrimmedVertexFinder PrimaryVertexProducerAlgorithm::theFinder
private

Definition at line 64 of file PrimaryVertexProducerAlgorithm.h.

Referenced by PrimaryVertexProducerAlgorithm(), and vertices().

VertexFitter<5>* PrimaryVertexProducerAlgorithm::theFitter
private
TrackClusterizerInZ* PrimaryVertexProducerAlgorithm::theTrackClusterizer
private
TrackFilterForPVFindingBase* PrimaryVertexProducerAlgorithm::theTrackFilter
private
VertexCompatibleWithBeam PrimaryVertexProducerAlgorithm::theVertexSelector
private

Definition at line 65 of file PrimaryVertexProducerAlgorithm.h.

Referenced by vertices().