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 hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () 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 56 of file PrimaryVertexProducer.h.

Constructor & Destructor Documentation

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

Definition at line 21 of file PrimaryVertexProducer.cc.

References electronCleaner_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, HIPixelAdaptiveVertex_cfi::vertexCollections, and PrimaryVertexProducer::algo::vertexSelector.

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

Definition at line 116 of file PrimaryVertexProducer.cc.

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

116  {
117  if (theTrackFilter) delete theTrackFilter;
119  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
120  if (algorithm->fitter) delete algorithm->fitter;
121  if (algorithm->vertexSelector) delete algorithm->vertexSelector;
122  }
123 }
TrackClusterizerInZ * theTrackClusterizer
TrackFilterForPVFindingBase * theTrackFilter
std::vector< algo > algorithms

Member Function Documentation

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

Definition at line 64 of file PrimaryVertexProducer.h.

References theConfig.

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

Definition at line 127 of file PrimaryVertexProducer.cc.

References electronCleaner_cfi::algorithm, algorithms, ecalDrivenElectronSeedsParameters_cff::beamSpot, bsToken, TrackClusterizerInZ::clusterize(), fastPrimaryVertexProducer_cfi::clusters, gather_cfg::cout, GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::cyy(), GlobalErrorBase< T, ErrorWeightType >::czz(), DEFINE_FWK_MODULE, TransientVertex::degreesOfFreedom(), 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(), 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().

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

Member Data Documentation

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

Definition at line 80 of file PrimaryVertexProducer.h.

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

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

Definition at line 85 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

bool PrimaryVertexProducer::f4D
private

Definition at line 90 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

bool PrimaryVertexProducer::fVerbose
private

Definition at line 83 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

edm::ParameterSet PrimaryVertexProducer::theConfig
private

Definition at line 82 of file PrimaryVertexProducer.h.

Referenced by config().

TrackClusterizerInZ* PrimaryVertexProducer::theTrackClusterizer
private

Definition at line 69 of file PrimaryVertexProducer.h.

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

TrackFilterForPVFindingBase* PrimaryVertexProducer::theTrackFilter
private

Definition at line 68 of file PrimaryVertexProducer.h.

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

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

Definition at line 88 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

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

Definition at line 87 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

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

Definition at line 86 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().