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