CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Attributes
PrimaryVertexProducer Class Reference

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

Inheritance diagram for PrimaryVertexProducer:
edm::stream::EDProducer<>

Classes

struct  algo
 

Public Member Functions

edm::ParameterSet config () const
 
 PrimaryVertexProducer (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~PrimaryVertexProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

Private Attributes

std::vector< algoalgorithms
 
edm::EDGetTokenT< reco::BeamSpotbsToken
 
bool f4D
 
bool fVerbose
 
edm::ParameterSet theConfig
 
TrackClusterizerInZtheTrackClusterizer
 
TrackFilterForPVFindingBasetheTrackFilter
 
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
 
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
 
edm::EDGetTokenT< reco::TrackCollectiontrkToken
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Description: steers tracker primary vertex reconstruction and storage

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

Definition at line 54 of file PrimaryVertexProducer.h.

Constructor & Destructor Documentation

PrimaryVertexProducer::PrimaryVertexProducer ( const edm::ParameterSet conf)

Definition at line 21 of file PrimaryVertexProducer.cc.

References qcdUeDQM_cfi::algorithm, algorithms, bsToken, edm::ParameterSet::exists(), f4D, PrimaryVertexProducer::algo::fitter, fVerbose, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), PrimaryVertexProducer::algo::label, PrimaryVertexProducer::algo::minNdof, AlCaHLTBitMon_QueryRunRegistry::string, theTrackClusterizer, theTrackFilter, trkTimeResosToken, trkTimesToken, trkToken, PrimaryVertexProducer::algo::useBeamConstraint, HLT_2018_cff::vertexCollections, and PrimaryVertexProducer::algo::vertexSelector.

21  : theConfig(conf) {
22  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
23 
24  trkToken = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackLabel"));
25  bsToken = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotLabel"));
26  f4D = false;
27 
28  // select and configure the track selection
29  std::string trackSelectionAlgorithm =
30  conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
31  if (trackSelectionAlgorithm == "filter") {
32  theTrackFilter = new TrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
33  } else if (trackSelectionAlgorithm == "filterWithThreshold") {
35  } else {
36  throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " +
37  trackSelectionAlgorithm);
38  }
39 
40  // select and configure the track clusterizer
41  std::string clusteringAlgorithm =
42  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
43  if (clusteringAlgorithm == "gap") {
45  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
46  } else if (clusteringAlgorithm == "DA") {
48  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
49  }
50  // provide the vectorized version of the clusterizer, if supported by the build
51  else if (clusteringAlgorithm == "DA_vect") {
53  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
54  } else if (clusteringAlgorithm == "DA2D_vect") {
56  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
57  f4D = true;
58  }
59 
60  else {
61  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
62  }
63 
64  if (f4D) {
65  trkTimesToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimesLabel"));
66  trkTimeResosToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimeResosLabel"));
67  }
68 
69  // select and configure the vertex fitters
70  if (conf.exists("vertexCollections")) {
71  std::vector<edm::ParameterSet> vertexCollections =
72  conf.getParameter<std::vector<edm::ParameterSet> >("vertexCollections");
73 
74  for (std::vector<edm::ParameterSet>::const_iterator algoconf = vertexCollections.begin();
75  algoconf != vertexCollections.end();
76  algoconf++) {
78  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
79  if (fitterAlgorithm == "KalmanVertexFitter") {
80  algorithm.fitter = new KalmanVertexFitter();
81  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
82  algorithm.fitter = new AdaptiveVertexFitter(GeometricAnnealing(algoconf->getParameter<double>("chi2cutoff")));
83  } else {
84  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
85  }
86  algorithm.label = algoconf->getParameter<std::string>("label");
87  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
88  algorithm.useBeamConstraint = algoconf->getParameter<bool>("useBeamConstraint");
89  algorithm.vertexSelector =
90  new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
91  algorithms.push_back(algorithm);
92 
93  produces<reco::VertexCollection>(algorithm.label);
94  }
95  } else {
96  edm::LogWarning("MisConfiguration")
97  << "this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
98 
100  std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
101  if (fitterAlgorithm == "KalmanVertexFitter") {
102  algorithm.fitter = new KalmanVertexFitter();
103  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
104  algorithm.fitter = new AdaptiveVertexFitter();
105  } else {
106  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
107  }
108  algorithm.label = "";
109  algorithm.minNdof = conf.getParameter<double>("minNdof");
110  algorithm.useBeamConstraint = conf.getParameter<bool>("useBeamConstraint");
111 
112  algorithm.vertexSelector = new VertexCompatibleWithBeam(
114  conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
115 
116  algorithms.push_back(algorithm);
117  produces<reco::VertexCollection>(algorithm.label);
118  }
119 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::ParameterSet theConfig
edm::EDGetTokenT< reco::BeamSpot > bsToken
std::vector< algo > algorithms
Common base class.
edm::EDGetTokenT< reco::TrackCollection > trkToken
TrackClusterizerInZ * theTrackClusterizer
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
TrackFilterForPVFindingBase * theTrackFilter
PrimaryVertexProducer::~PrimaryVertexProducer ( )
override

Definition at line 121 of file PrimaryVertexProducer.cc.

References qcdUeDQM_cfi::algorithm, algorithms, theTrackClusterizer, and theTrackFilter.

121  {
122  if (theTrackFilter)
123  delete theTrackFilter;
125  delete theTrackClusterizer;
126  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
127  if (algorithm->fitter)
128  delete algorithm->fitter;
129  if (algorithm->vertexSelector)
130  delete algorithm->vertexSelector;
131  }
132 }
std::vector< algo > algorithms
TrackClusterizerInZ * theTrackClusterizer
TrackFilterForPVFindingBase * theTrackFilter

Member Function Documentation

edm::ParameterSet PrimaryVertexProducer::config ( void  ) const
inline

Definition at line 62 of file PrimaryVertexProducer.h.

References theConfig.

62 { return theConfig; }
edm::ParameterSet theConfig
void PrimaryVertexProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 134 of file PrimaryVertexProducer.cc.

References qcdUeDQM_cfi::algorithm, algorithms, pwdgSkimBPark_cfi::beamSpot, bsToken, TrackClusterizerInZ::clusterize(), bsc_activity_cfg::clusters, gather_cfg::cout, GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::cyy(), GlobalErrorBase< T, ErrorWeightType >::czz(), DEFINE_FWK_MODULE, TransientVertex::degreesOfFreedom(), runTheMatrix::err, VertexState::error(), f4D, fVerbose, edm::EventSetup::get(), edm::Event::getByToken(), edm::HandleBase::isValid(), TransientVertex::isValid(), GlobalErrorBase< T, ErrorWeightType >::matrix(), GlobalErrorBase< T, ErrorWeightType >::matrix4D(), eostools::move(), TransientVertex::originalTracks(), reco::BeamSpot::position(), TransientVertex::position(), TransientVertex::positionError(), edm::Handle< T >::product(), edm::Event::put(), FSQDQM_cfi::pvs, mps_fire::result, reco::BeamSpot::rotatedCovariance3D(), TrackFilterForPVFindingBase::select(), theTrackClusterizer, theTrackFilter, TransientVertex::time(), ntuplemaker::time, TransientVertex::totalChiSquared(), trkTimeResosToken, trkTimesToken, trkToken, findQualityFiles::v, w, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

134  {
135  // get the BeamSpot, it will alwys be needed, even when not used as a constraint
137  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
138  iEvent.getByToken(bsToken, recoBeamSpotHandle);
139  if (recoBeamSpotHandle.isValid()) {
140  beamSpot = *recoBeamSpotHandle;
141  } else {
142  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
143  }
144 
145  bool validBS = true;
146  VertexState beamVertexState(beamSpot);
147  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
148  (beamVertexState.error().czz() <= 0.)) {
149  validBS = false;
150  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
151  }
152 
153  // get RECO tracks from the event
154  // `tks` can be used as a ptr to a reco::TrackCollection
156  iEvent.getByToken(trkToken, tks);
157 
158  // interface RECO tracks to vertex reconstruction
160  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
161  std::vector<reco::TransientTrack> t_tks;
162 
163  if (f4D) {
164  edm::Handle<edm::ValueMap<float> > trackTimesH;
165  edm::Handle<edm::ValueMap<float> > trackTimeResosH;
166  iEvent.getByToken(trkTimesToken, trackTimesH);
167  iEvent.getByToken(trkTimeResosToken, trackTimeResosH);
168  t_tks = (*theB).build(tks, beamSpot, *(trackTimesH.product()), *(trackTimeResosH.product()));
169  } else {
170  t_tks = (*theB).build(tks, beamSpot);
171  }
172  if (fVerbose) {
173  std::cout << "RecoVertex/PrimaryVertexProducer"
174  << "Found: " << t_tks.size() << " reconstructed tracks"
175  << "\n";
176  }
177 
178  // select tracks
179  std::vector<reco::TransientTrack>&& seltks = theTrackFilter->select(t_tks);
180 
181  // clusterize tracks in Z
182  std::vector<std::vector<reco::TransientTrack> >&& clusters = theTrackClusterizer->clusterize(seltks);
183 
184  if (fVerbose) {
185  std::cout << " clustering returned " << clusters.size() << " clusters from " << seltks.size()
186  << " selected tracks" << std::endl;
187  }
188 
189  // vertex fits
190  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
191  auto result = std::make_unique<reco::VertexCollection>();
192  reco::VertexCollection& vColl = (*result);
193 
194  std::vector<TransientVertex> pvs;
195  for (std::vector<std::vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin();
196  iclus != clusters.end();
197  iclus++) {
198  double sumwt = 0.;
199  double sumwt2 = 0.;
200  double sumw = 0.;
201  double meantime = 0.;
202  double vartime = 0.;
203  if (f4D) {
204  for (const auto& tk : *iclus) {
205  const double time = tk.timeExt();
206  const double err = tk.dtErrorExt();
207  const double inverr = err > 0. ? 1.0 / err : 0.;
208  const double w = inverr * inverr;
209  sumwt += w * time;
210  sumwt2 += w * time * time;
211  sumw += w;
212  }
213  meantime = sumwt / sumw;
214  double sumsq = sumwt2 - sumwt * sumwt / sumw;
215  double chisq = iclus->size() > 1 ? sumsq / double(iclus->size() - 1) : sumsq / double(iclus->size());
216  vartime = chisq / sumw;
217  }
218 
220  if (algorithm->useBeamConstraint && validBS && ((*iclus).size() > 1)) {
221  v = algorithm->fitter->vertex(*iclus, beamSpot);
222 
223  if (f4D) {
224  if (v.isValid()) {
225  auto err = v.positionError().matrix4D();
226  err(3, 3) = vartime;
227  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared());
228  }
229  }
230 
231  } else if (!(algorithm->useBeamConstraint) && ((*iclus).size() > 1)) {
232  v = algorithm->fitter->vertex(*iclus);
233 
234  if (f4D) {
235  if (v.isValid()) {
236  auto err = v.positionError().matrix4D();
237  err(3, 3) = vartime;
238  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared());
239  }
240  }
241 
242  } // else: no fit ==> v.isValid()=False
243 
244  if (fVerbose) {
245  if (v.isValid()) {
246  std::cout << "x,y,z";
247  if (f4D)
248  std::cout << ",t";
249  std::cout << "=" << v.position().x() << " " << v.position().y() << " " << v.position().z();
250  if (f4D)
251  std::cout << " " << v.time();
252  std::cout << " cluster size = " << (*iclus).size() << std::endl;
253  } else {
254  std::cout << "Invalid fitted vertex, cluster size=" << (*iclus).size() << std::endl;
255  }
256  }
257 
258  if (v.isValid() && (v.degreesOfFreedom() >= algorithm->minNdof) &&
259  (!validBS || (*(algorithm->vertexSelector))(v, beamVertexState)))
260  pvs.push_back(v);
261  } // end of cluster loop
262 
263  if (fVerbose) {
264  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
265  }
266 
267  if (clusters.size() > 2 && clusters.size() > 2 * pvs.size())
268  edm::LogWarning("PrimaryVertexProducer")
269  << "more than half of candidate vertices lost " << pvs.size() << ' ' << clusters.size();
270 
271  if (pvs.empty() && seltks.size() > 5)
272  edm::LogWarning("PrimaryVertexProducer")
273  << "no vertex found with " << seltks.size() << " tracks and " << clusters.size() << " vertex-candidates";
274 
275  // sort vertices by pt**2 vertex (aka signal vertex tagging)
276  if (pvs.size() > 1) {
277  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
278  }
279 
280  // convert transient vertices returned by the theAlgo to (reco) vertices
281  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin(); iv != pvs.end(); iv++) {
282  reco::Vertex v = *iv;
283  vColl.push_back(v);
284  }
285 
286  if (vColl.empty()) {
287  GlobalError bse(beamSpot.rotatedCovariance3D());
288  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
290  we(0, 0) = 10000;
291  we(1, 1) = 10000;
292  we(2, 2) = 10000;
293  vColl.push_back(reco::Vertex(beamSpot.position(), we, 0., 0., 0));
294  if (fVerbose) {
295  std::cout << "RecoVertex/PrimaryVertexProducer: "
296  << "Beamspot with invalid errors " << bse.matrix() << std::endl;
297  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
298  }
299  } else {
300  vColl.push_back(reco::Vertex(beamSpot.position(), beamSpot.rotatedCovariance3D(), 0., 0., 0));
301  if (fVerbose) {
302  std::cout << "RecoVertex/PrimaryVertexProducer: "
303  << " will put Vertex derived from BeamSpot into Event.\n";
304  }
305  }
306  }
307 
308  if (fVerbose) {
309  int ivtx = 0;
310  for (reco::VertexCollection::const_iterator v = vColl.begin(); v != vColl.end(); ++v) {
311  std::cout << "recvtx " << ivtx++ << "#trk " << std::setw(3) << v->tracksSize() << " chi2 " << std::setw(4)
312  << v->chi2() << " ndof " << std::setw(3) << v->ndof() << " x " << std::setw(6) << v->position().x()
313  << " dx " << std::setw(6) << v->xError() << " y " << std::setw(6) << v->position().y() << " dy "
314  << std::setw(6) << v->yError() << " z " << std::setw(6) << v->position().z() << " dz " << std::setw(6)
315  << v->zError();
316  if (f4D) {
317  std::cout << " t " << std::setw(6) << v->t() << " dt " << std::setw(6) << v->tError();
318  }
319  std::cout << std::endl;
320  }
321  }
322 
323  iEvent.put(std::move(result), algorithm->label);
324  }
325 }
GlobalError positionError() const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
edm::EDGetTokenT< reco::BeamSpot > bsToken
const double w
Definition: UKUtility.cc:23
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
std::vector< algo > algorithms
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
double zError() const
error on z
Definition: Vertex.h:127
float totalChiSquared() const
edm::EDGetTokenT< reco::TrackCollection > trkToken
TrackClusterizerInZ * theTrackClusterizer
T y() const
Definition: PV3DBase.h:60
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const Point & position() const
position
Definition: Vertex.h:113
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
const AlgebraicSymMatrix44 & matrix4D() const
std::vector< reco::TransientTrack > const & originalTracks() const
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
float degreesOfFreedom() const
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:61
double chi2() const
chi-squares
Definition: Vertex.h:102
TrackFilterForPVFindingBase * theTrackFilter
bool isValid() const
Definition: HandleBase.h:70
double ndof() const
Definition: Vertex.h:109
double xError() const
error on x
Definition: Vertex.h:123
T const * product() const
Definition: Handle.h:69
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
double time() const
T get() const
Definition: EventSetup.h:73
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
const Point & position() const
position
Definition: BeamSpot.h:59
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:73
T x() const
Definition: PV3DBase.h:59
bool isValid() const
double yError() const
error on y
Definition: Vertex.h:125
def move(src, dest)
Definition: eostools.py:511
double tError() const
error on t
Definition: Vertex.h:129
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:69
double t() const
t coordinate
Definition: Vertex.h:121

Member Data Documentation

std::vector<algo> PrimaryVertexProducer::algorithms
private

Definition at line 78 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), produce(), and ~PrimaryVertexProducer().

edm::EDGetTokenT<reco::BeamSpot> PrimaryVertexProducer::bsToken
private

Definition at line 83 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

bool PrimaryVertexProducer::f4D
private

Definition at line 88 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

bool PrimaryVertexProducer::fVerbose
private

Definition at line 81 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

edm::ParameterSet PrimaryVertexProducer::theConfig
private

Definition at line 80 of file PrimaryVertexProducer.h.

Referenced by config().

TrackClusterizerInZ* PrimaryVertexProducer::theTrackClusterizer
private

Definition at line 67 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), produce(), and ~PrimaryVertexProducer().

TrackFilterForPVFindingBase* PrimaryVertexProducer::theTrackFilter
private

Definition at line 66 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), produce(), and ~PrimaryVertexProducer().

edm::EDGetTokenT<edm::ValueMap<float> > PrimaryVertexProducer::trkTimeResosToken
private

Definition at line 86 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

edm::EDGetTokenT<edm::ValueMap<float> > PrimaryVertexProducer::trkTimesToken
private

Definition at line 85 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

edm::EDGetTokenT<reco::TrackCollection> PrimaryVertexProducer::trkToken
private

Definition at line 84 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().