CMS 3D CMS Logo

PrimaryVertexProducerAlgorithm.cc
Go to the documentation of this file.
3 
9 
13 
18 
20  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
21  trackLabel = conf.getParameter<edm::InputTag>("TrackLabel");
22  beamSpotLabel = conf.getParameter<edm::InputTag>("beamSpotLabel");
23 
24  // select and configure the track selection
25  std::string trackSelectionAlgorithm =
26  conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
27  if (trackSelectionAlgorithm == "filter") {
28  theTrackFilter = new TrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
29  } else if (trackSelectionAlgorithm == "filterWithThreshold") {
31  } else {
32  throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " +
33  trackSelectionAlgorithm);
34  }
35 
36  // select and configure the track clusterizer
37  std::string clusteringAlgorithm =
38  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
39  if (clusteringAlgorithm == "gap") {
41  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
42  }
43  // provide the vectorized version of the clusterizer, if supported by the build
44  else if (clusteringAlgorithm == "DA_vect") {
46  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
47  }
48 
49  else {
50  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
51  }
52 
53  // select and configure the vertex fitters
54  std::vector<edm::ParameterSet> vertexCollections =
55  conf.getParameter<std::vector<edm::ParameterSet> >("vertexCollections");
56 
57  for (std::vector<edm::ParameterSet>::const_iterator algoconf = vertexCollections.begin();
58  algoconf != vertexCollections.end();
59  algoconf++) {
61  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
62  if (fitterAlgorithm == "KalmanVertexFitter") {
63  algorithm.fitter = new KalmanVertexFitter();
64  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
65  algorithm.fitter = new AdaptiveVertexFitter();
66  } else {
67  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
68  }
69  algorithm.label = algoconf->getParameter<std::string>("label");
70  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
71  algorithm.useBeamConstraint = algoconf->getParameter<bool>("useBeamConstraint");
72  algorithm.vertexSelector =
73  new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
74  algorithms.push_back(algorithm);
75  }
76 }
77 
79  if (theTrackFilter)
80  delete theTrackFilter;
82  delete theTrackClusterizer;
83  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
84  if (algorithm->fitter)
85  delete algorithm->fitter;
86  if (algorithm->vertexSelector)
87  delete algorithm->vertexSelector;
88  }
89 }
90 
91 //
92 // member functions
93 //
94 
95 // obsolete method, unfortunately required through inheritance from VertexReconstructor
96 std::vector<TransientVertex> PrimaryVertexProducerAlgorithm::vertices(
97  const std::vector<reco::TransientTrack>& tracks) const {
98  throw VertexException("PrimaryVertexProducerAlgorithm: cannot make a Primary Vertex without a beam spot");
99 
100  return std::vector<TransientVertex>();
101 }
102 
103 std::vector<TransientVertex> PrimaryVertexProducerAlgorithm::vertices(const std::vector<reco::TransientTrack>& t_tks,
104  const reco::BeamSpot& beamSpot,
105  const std::string& label) const {
106  bool validBS = true;
107  VertexState beamVertexState(beamSpot);
108  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
109  (beamVertexState.error().czz() <= 0.)) {
110  validBS = false;
111  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
112  }
113 
114  // // get RECO tracks from the event
115  // // `tks` can be used as a ptr to a reco::TrackCollection
116  // edm::Handle<reco::TrackCollection> tks;
117  // iEvent.getByLabel(trackLabel, tks);
118 
119  // select tracks
120  std::vector<reco::TransientTrack> seltks = theTrackFilter->select(t_tks);
121 
122  // clusterize tracks in Z
123  std::vector<std::vector<reco::TransientTrack> > clusters = theTrackClusterizer->clusterize(seltks);
124  if (fVerbose) {
125  std::cout << " clustering returned " << clusters.size() << " clusters from " << seltks.size()
126  << " selected tracks" << std::endl;
127  }
128 
129  // vertex fits
130  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
131  if (!(algorithm->label == label))
132  continue;
133 
134  //std::auto_ptr<reco::VertexCollection> result(new reco::VertexCollection);
135  // reco::VertexCollection vColl;
136 
137  std::vector<TransientVertex> pvs;
138  for (std::vector<std::vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin();
139  iclus != clusters.end();
140  iclus++) {
142  if (algorithm->useBeamConstraint && validBS && ((*iclus).size() > 1)) {
143  v = algorithm->fitter->vertex(*iclus, beamSpot);
144 
145  } else if (!(algorithm->useBeamConstraint) && ((*iclus).size() > 1)) {
146  v = algorithm->fitter->vertex(*iclus);
147 
148  } // else: no fit ==> v.isValid()=False
149 
150  if (fVerbose) {
151  if (v.isValid())
152  std::cout << "x,y,z=" << v.position().x() << " " << v.position().y() << " " << v.position().z() << std::endl;
153  else
154  std::cout << "Invalid fitted vertex\n";
155  }
156 
157  if (v.isValid() && (v.degreesOfFreedom() >= algorithm->minNdof) &&
158  (!validBS || (*(algorithm->vertexSelector))(v, beamVertexState)))
159  pvs.push_back(v);
160  } // end of cluster loop
161 
162  if (fVerbose) {
163  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
164  }
165 
166  // sort vertices by pt**2 vertex (aka signal vertex tagging)
167  if (pvs.size() > 1) {
168  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
169  }
170 
171  return pvs;
172  }
173 
174  std::vector<TransientVertex> dummy;
175  return dummy; //avoid compiler warning, should never be here
176 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
Common base class.
Log< level::Error, false > LogError
T getUntrackedParameter(std::string const &, T const &) const
char const * label
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
GlobalError error() const
Definition: VertexState.h:64
const AlgebraicSymMatrix33 matrix() const
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
PrimaryVertexProducerAlgorithm(const edm::ParameterSet &)
TrackFilterForPVFindingBase * theTrackFilter