CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static 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
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

std::vector< algoalgorithms
 
edm::EDGetTokenT< reco::BeamSpotbsToken
 
bool f4D
 
bool fRecoveryIteration
 
bool fVerbose
 
edm::EDGetTokenT< reco::VertexCollectionrecoveryVtxToken
 
edm::ParameterSet theConfig
 
TrackClusterizerInZtheTrackClusterizer
 
TrackFilterForPVFindingBasetheTrackFilter
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordtheTTBToken
 
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<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::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::PrimaryVertexProducer ( const edm::ParameterSet conf)

Definition at line 18 of file PrimaryVertexProducer.cc.

References qcdUeDQM_cfi::algorithm, algorithms, bsToken, f4D, fRecoveryIteration, fVerbose, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), recoveryVtxToken, AlCaHLTBitMon_QueryRunRegistry::string, theTrackClusterizer, theTrackFilter, trkTimeResosToken, trkTimesToken, trkToken, and HLT_2022v12_cff::vertexCollections.

19  : theTTBToken(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))), theConfig(conf) {
20  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
21 
22  trkToken = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackLabel"));
23  bsToken = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotLabel"));
24  f4D = false;
25 
26  // select and configure the track selection
27  std::string trackSelectionAlgorithm =
28  conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
29  if (trackSelectionAlgorithm == "filter") {
30  theTrackFilter = new TrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
31  } else if (trackSelectionAlgorithm == "filterWithThreshold") {
33  } else {
34  throw VertexException("PrimaryVertexProducer: unknown track selection algorithm: " + trackSelectionAlgorithm);
35  }
36 
37  // select and configure the track clusterizer
38  std::string clusteringAlgorithm =
39  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
40  if (clusteringAlgorithm == "gap") {
42  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
43  } else if (clusteringAlgorithm == "DA") {
45  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
46  }
47  // provide the vectorized version of the clusterizer, if supported by the build
48  else if (clusteringAlgorithm == "DA_vect") {
50  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
51  } else if (clusteringAlgorithm == "DA2D_vect") {
53  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
54  f4D = true;
55  }
56 
57  else {
58  throw VertexException("PrimaryVertexProducer: 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  // select and configure the vertex fitters
67  std::vector<edm::ParameterSet> vertexCollections =
68  conf.getParameter<std::vector<edm::ParameterSet> >("vertexCollections");
69 
70  for (std::vector<edm::ParameterSet>::const_iterator algoconf = vertexCollections.begin();
71  algoconf != vertexCollections.end();
72  algoconf++) {
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("PrimaryVertexProducer: 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 =
86  new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
87  algorithms.push_back(algorithm);
88 
89  produces<reco::VertexCollection>(algorithm.label);
90  }
91 
92  //check if this is a recovery iteration
93  fRecoveryIteration = conf.getParameter<bool>("isRecoveryIteration");
94  if (fRecoveryIteration) {
95  if (algorithms.empty()) {
96  throw VertexException("PrimaryVertexProducer: No algorithm specified. ");
97  } else if (algorithms.size() > 1) {
98  throw VertexException(
99  "PrimaryVertexProducer: Running in Recovery mode and more than one algorithm specified. Please "
100  "only one algorithm.");
101  }
102  recoveryVtxToken = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("recoveryVtxCollection"));
103  }
104 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ParameterSet theConfig
edm::EDGetTokenT< reco::BeamSpot > bsToken
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > theTTBToken
std::vector< algo > algorithms
Common base class.
edm::EDGetTokenT< reco::VertexCollection > recoveryVtxToken
edm::EDGetTokenT< reco::TrackCollection > trkToken
TrackClusterizerInZ * theTrackClusterizer
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
TrackFilterForPVFindingBase * theTrackFilter

◆ ~PrimaryVertexProducer()

PrimaryVertexProducer::~PrimaryVertexProducer ( )
override

Definition at line 106 of file PrimaryVertexProducer.cc.

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

106  {
107  if (theTrackFilter)
108  delete theTrackFilter;
110  delete theTrackClusterizer;
111  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
112  if (algorithm->fitter)
113  delete algorithm->fitter;
114  if (algorithm->vertexSelector)
115  delete algorithm->vertexSelector;
116  }
117 }
std::vector< algo > algorithms
TrackClusterizerInZ * theTrackClusterizer
TrackFilterForPVFindingBase * theTrackFilter

Member Function Documentation

◆ config()

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

Definition at line 66 of file PrimaryVertexProducer.h.

References theConfig.

66 { return theConfig; }
edm::ParameterSet theConfig

◆ fillDescriptions()

void PrimaryVertexProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 318 of file PrimaryVertexProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ParameterSet::addParameter(), submitPVResolutionJobs::desc, GapClusterizerInZ::fillPSetDescription(), TrackFilterForPVFinding::fillPSetDescription(), DAClusterizerInZT_vect::fillPSetDescription(), HLT_2022v12_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

318  {
319  // offlinePrimaryVertices
321  {
323  vpsd1.add<double>("maxDistanceToBeam", 1.0);
324  vpsd1.add<std::string>("algorithm", "AdaptiveVertexFitter");
325  vpsd1.add<bool>("useBeamConstraint", false);
326  vpsd1.add<std::string>("label", "");
327  vpsd1.add<double>("chi2cutoff", 2.5);
328  vpsd1.add<double>("minNdof", 0.0);
329  std::vector<edm::ParameterSet> temp1;
330  temp1.reserve(2);
331  {
332  edm::ParameterSet temp2;
333  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
334  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
335  temp2.addParameter<bool>("useBeamConstraint", false);
336  temp2.addParameter<std::string>("label", "");
337  temp2.addParameter<double>("chi2cutoff", 2.5);
338  temp2.addParameter<double>("minNdof", 0.0);
339  temp1.push_back(temp2);
340  }
341  {
342  edm::ParameterSet temp2;
343  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
344  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
345  temp2.addParameter<bool>("useBeamConstraint", true);
346  temp2.addParameter<std::string>("label", "WithBS");
347  temp2.addParameter<double>("chi2cutoff", 2.5);
348  temp2.addParameter<double>("minNdof", 2.0);
349  temp1.push_back(temp2);
350  }
351  desc.addVPSet("vertexCollections", vpsd1, temp1);
352  }
353  desc.addUntracked<bool>("verbose", false);
354  {
357  psd0.add<int>("numTracksThreshold", 0); // HI only
358  desc.add<edm::ParameterSetDescription>("TkFilterParameters", psd0);
359  }
360  desc.add<edm::InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
361  desc.add<edm::InputTag>("TrackLabel", edm::InputTag("generalTracks"));
362  desc.add<edm::InputTag>("TrackTimeResosLabel", edm::InputTag("dummy_default")); // 4D only
363  desc.add<edm::InputTag>("TrackTimesLabel", edm::InputTag("dummy_default")); // 4D only
364 
365  {
367  {
370  psd0.add<edm::ParameterSetDescription>("TkDAClusParameters", psd1);
371 
374  psd0.add<edm::ParameterSetDescription>("TkGapClusParameters", psd2);
375  }
376  psd0.add<std::string>("algorithm", "DA_vect");
377  desc.add<edm::ParameterSetDescription>("TkClusParameters", psd0);
378  }
379 
380  desc.add<bool>("isRecoveryIteration", false);
381  desc.add<edm::InputTag>("recoveryVtxCollection", {""});
382 
383  descriptions.add("primaryVertexProducer", desc);
384 }
static void fillPSetDescription(edm::ParameterSetDescription &desc)
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillPSetDescription(edm::ParameterSetDescription &desc)

◆ produce()

void PrimaryVertexProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 119 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(), submitPVResolutionJobs::err, VertexState::error(), f4D, fRecoveryIteration, fVerbose, edm::EventSetup::getData(), iEvent, edm::HandleBase::isValid(), gpuVertexFinder::iv, GlobalErrorBase< T, ErrorWeightType >::matrix(), eostools::move(), edm::Handle< T >::product(), FSQDQM_cfi::pvs, recoveryVtxToken, mps_fire::result, TrackFilterForPVFindingBase::select(), jetUpdater_cfi::sort, theTrackClusterizer, theTrackFilter, theTTBToken, protons_cff::time, trkTimeResosToken, trkTimesToken, trkToken, findQualityFiles::v, trackerHitRTTI::vector, and w().

119  {
120  // get the BeamSpot, it will always be needed, even when not used as a constraint
122  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
123  iEvent.getByToken(bsToken, recoBeamSpotHandle);
124  if (recoBeamSpotHandle.isValid()) {
125  beamSpot = *recoBeamSpotHandle;
126  } else {
127  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
128  }
129 
130  bool validBS = true;
131  VertexState beamVertexState(beamSpot);
132  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
133  (beamVertexState.error().czz() <= 0.)) {
134  validBS = false;
135  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
136  }
137 
138  //if this is a recovery iteration, check if we already have a valid PV
139  if (fRecoveryIteration) {
140  auto const& oldVertices = iEvent.get(recoveryVtxToken);
141  //look for the first valid (not-BeamSpot) vertex
142  for (auto const& old : oldVertices) {
143  if (!(old.isFake())) {
144  //found a valid vertex, write the first one to the collection and return
145  //otherwise continue with regular vertexing procedure
146  auto result = std::make_unique<reco::VertexCollection>();
147  result->push_back(old);
148  iEvent.put(std::move(result), algorithms.begin()->label);
149  return;
150  }
151  }
152  }
153 
154  // get RECO tracks from the event
155  // `tks` can be used as a ptr to a reco::TrackCollection
157  iEvent.getByToken(trkToken, tks);
158 
159  // interface RECO tracks to vertex reconstruction
160  const auto& theB = &iSetup.getData(theTTBToken);
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  } else if (!(algorithm->useBeamConstraint) && (iclus->size() > 1)) {
223  v = algorithm->fitter->vertex(*iclus);
224  } // else: no fit ==> v.isValid()=False
225 
226  // 4D vertices: add timing information
227  if (f4D and v.isValid()) {
228  auto err = v.positionError().matrix4D();
229  err(3, 3) = vartime;
230  auto trkWeightMap3d = v.weightMap(); // copy the 3d-fit weights
231  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared(), v.degreesOfFreedom());
232  v.weightMap(trkWeightMap3d);
233  }
234 
235  if (fVerbose) {
236  if (v.isValid()) {
237  std::cout << "x,y,z";
238  if (f4D)
239  std::cout << ",t";
240  std::cout << "=" << v.position().x() << " " << v.position().y() << " " << v.position().z();
241  if (f4D)
242  std::cout << " " << v.time();
243  std::cout << " cluster size = " << (*iclus).size() << std::endl;
244  } else {
245  std::cout << "Invalid fitted vertex, cluster size=" << (*iclus).size() << std::endl;
246  }
247  }
248 
249  if (v.isValid() && (v.degreesOfFreedom() >= algorithm->minNdof) &&
250  (!validBS || (*(algorithm->vertexSelector))(v, beamVertexState)))
251  pvs.push_back(v);
252  } // end of cluster loop
253 
254  if (fVerbose) {
255  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
256  }
257 
258  if (clusters.size() > 2 && clusters.size() > 2 * pvs.size())
259  edm::LogWarning("PrimaryVertexProducer")
260  << "more than half of candidate vertices lost " << pvs.size() << ' ' << clusters.size();
261 
262  if (pvs.empty() && seltks.size() > 5)
263  edm::LogWarning("PrimaryVertexProducer")
264  << "no vertex found with " << seltks.size() << " tracks and " << clusters.size() << " vertex-candidates";
265 
266  // sort vertices by pt**2 vertex (aka signal vertex tagging)
267  if (pvs.size() > 1) {
268  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
269  }
270 
271  // convert transient vertices returned by the theAlgo to (reco) vertices
272  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin(); iv != pvs.end(); iv++) {
273  reco::Vertex v = *iv;
274  vColl.push_back(v);
275  }
276 
277  if (vColl.empty()) {
278  GlobalError bse(beamSpot.rotatedCovariance3D());
279  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
281  we(0, 0) = 10000;
282  we(1, 1) = 10000;
283  we(2, 2) = 10000;
284  vColl.push_back(reco::Vertex(beamSpot.position(), we, 0., 0., 0));
285  if (fVerbose) {
286  std::cout << "RecoVertex/PrimaryVertexProducer: "
287  << "Beamspot with invalid errors " << bse.matrix() << std::endl;
288  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
289  }
290  } else {
291  vColl.push_back(reco::Vertex(beamSpot.position(), beamSpot.rotatedCovariance3D(), 0., 0., 0));
292  if (fVerbose) {
293  std::cout << "RecoVertex/PrimaryVertexProducer: "
294  << " will put Vertex derived from BeamSpot into Event.\n";
295  }
296  }
297  }
298 
299  if (fVerbose) {
300  int ivtx = 0;
301  for (reco::VertexCollection::const_iterator v = vColl.begin(); v != vColl.end(); ++v) {
302  std::cout << "recvtx " << ivtx++ << "#trk " << std::setw(3) << v->tracksSize() << " chi2 " << std::setw(4)
303  << v->chi2() << " ndof " << std::setw(3) << v->ndof() << " x " << std::setw(6) << v->position().x()
304  << " dx " << std::setw(6) << v->xError() << " y " << std::setw(6) << v->position().y() << " dy "
305  << std::setw(6) << v->yError() << " z " << std::setw(6) << v->position().z() << " dz " << std::setw(6)
306  << v->zError();
307  if (f4D) {
308  std::cout << " t " << std::setw(6) << v->t() << " dt " << std::setw(6) << v->tError();
309  }
310  std::cout << std::endl;
311  }
312  }
313 
314  iEvent.put(std::move(result), algorithm->label);
315  }
316 }
edm::EDGetTokenT< reco::BeamSpot > bsToken
int32_t *__restrict__ iv
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > theTTBToken
std::vector< algo > algorithms
T w() const
edm::EDGetTokenT< reco::VertexCollection > recoveryVtxToken
T const * product() const
Definition: Handle.h:70
edm::EDGetTokenT< reco::TrackCollection > trkToken
TrackClusterizerInZ * theTrackClusterizer
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Log< level::Error, false > LogError
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
TrackFilterForPVFindingBase * theTrackFilter
bool getData(T &iHolder) const
Definition: EventSetup.h:122
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
bool isValid() const
Definition: HandleBase.h:70
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ algorithms

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

Definition at line 84 of file PrimaryVertexProducer.h.

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

◆ bsToken

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

Definition at line 92 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ f4D

bool PrimaryVertexProducer::f4D
private

Definition at line 97 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ fRecoveryIteration

bool PrimaryVertexProducer::fRecoveryIteration
private

Definition at line 89 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ fVerbose

bool PrimaryVertexProducer::fVerbose
private

Definition at line 87 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ recoveryVtxToken

edm::EDGetTokenT<reco::VertexCollection> PrimaryVertexProducer::recoveryVtxToken
private

Definition at line 90 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ theConfig

edm::ParameterSet PrimaryVertexProducer::theConfig
private

Definition at line 86 of file PrimaryVertexProducer.h.

Referenced by config().

◆ theTrackClusterizer

TrackClusterizerInZ* PrimaryVertexProducer::theTrackClusterizer
private

Definition at line 73 of file PrimaryVertexProducer.h.

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

◆ theTrackFilter

TrackFilterForPVFindingBase* PrimaryVertexProducer::theTrackFilter
private

Definition at line 72 of file PrimaryVertexProducer.h.

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

◆ theTTBToken

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> PrimaryVertexProducer::theTTBToken
private

Definition at line 70 of file PrimaryVertexProducer.h.

Referenced by produce().

◆ trkTimeResosToken

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

Definition at line 95 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ trkTimesToken

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

Definition at line 94 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ trkToken

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

Definition at line 93 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().