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
 
bool weightFit
 

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 57 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, HLT_2022v15_cff::vertexCollections, and weightFit.

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  weightFit = false;
26 
27  // select and configure the track selection
28  std::string trackSelectionAlgorithm =
29  conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
30  if (trackSelectionAlgorithm == "filter") {
31  theTrackFilter = new TrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
32  } else if (trackSelectionAlgorithm == "filterWithThreshold") {
34  } else {
35  throw VertexException("PrimaryVertexProducer: unknown track selection algorithm: " + trackSelectionAlgorithm);
36  }
37 
38  // select and configure the track clusterizer
39  std::string clusteringAlgorithm =
40  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
41  if (clusteringAlgorithm == "gap") {
43  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
44  } else if (clusteringAlgorithm == "DA") {
46  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
47  }
48  // provide the vectorized version of the clusterizer, if supported by the build
49  else if (clusteringAlgorithm == "DA_vect") {
51  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
52  } else if (clusteringAlgorithm == "DA2D_vect") {
54  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
55  f4D = true;
56  }
57 
58  else {
59  throw VertexException("PrimaryVertexProducer: unknown clustering algorithm: " + clusteringAlgorithm);
60  }
61 
62  if (f4D) {
63  trkTimesToken = consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("TrackTimesLabel"));
64  trkTimeResosToken = consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("TrackTimeResosLabel"));
65  }
66 
67  // select and configure the vertex fitters
68  std::vector<edm::ParameterSet> vertexCollections =
69  conf.getParameter<std::vector<edm::ParameterSet>>("vertexCollections");
70 
71  for (std::vector<edm::ParameterSet>::const_iterator algoconf = vertexCollections.begin();
72  algoconf != vertexCollections.end();
73  algoconf++) {
75  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
76  if (fitterAlgorithm == "KalmanVertexFitter") {
77  algorithm.fitter = new KalmanVertexFitter();
78  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
79  algorithm.fitter = new AdaptiveVertexFitter(GeometricAnnealing(algoconf->getParameter<double>("chi2cutoff")));
80  } else if (fitterAlgorithm == "WeightedMeanFitter") {
81  algorithm.fitter = nullptr;
82  weightFit = true;
83  } else {
84  throw VertexException("PrimaryVertexProducer: 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 
96  //check if this is a recovery iteration
97  fRecoveryIteration = conf.getParameter<bool>("isRecoveryIteration");
98  if (fRecoveryIteration) {
99  if (algorithms.empty()) {
100  throw VertexException("PrimaryVertexProducer: No algorithm specified. ");
101  } else if (algorithms.size() > 1) {
102  throw VertexException(
103  "PrimaryVertexProducer: Running in Recovery mode and more than one algorithm specified. Please "
104  "only one algorithm.");
105  }
106  recoveryVtxToken = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("recoveryVtxCollection"));
107  }
108 }
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 110 of file PrimaryVertexProducer.cc.

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

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

Member Function Documentation

◆ config()

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

Definition at line 67 of file PrimaryVertexProducer.h.

References theConfig.

67 { return theConfig; }
edm::ParameterSet theConfig

◆ fillDescriptions()

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

Definition at line 355 of file PrimaryVertexProducer.cc.

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

355  {
356  // offlinePrimaryVertices
358  {
360  vpsd1.add<double>("maxDistanceToBeam", 1.0);
361  vpsd1.add<std::string>("algorithm", "AdaptiveVertexFitter");
362  vpsd1.add<bool>("useBeamConstraint", false);
363  vpsd1.add<std::string>("label", "");
364  vpsd1.add<double>("chi2cutoff", 2.5);
365  vpsd1.add<double>("minNdof", 0.0);
366  std::vector<edm::ParameterSet> temp1;
367  temp1.reserve(2);
368  {
369  edm::ParameterSet temp2;
370  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
371  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
372  temp2.addParameter<bool>("useBeamConstraint", false);
373  temp2.addParameter<std::string>("label", "");
374  temp2.addParameter<double>("chi2cutoff", 2.5);
375  temp2.addParameter<double>("minNdof", 0.0);
376  temp1.push_back(temp2);
377  }
378  {
379  edm::ParameterSet temp2;
380  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
381  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
382  temp2.addParameter<bool>("useBeamConstraint", true);
383  temp2.addParameter<std::string>("label", "WithBS");
384  temp2.addParameter<double>("chi2cutoff", 2.5);
385  temp2.addParameter<double>("minNdof", 2.0);
386  temp1.push_back(temp2);
387  }
388  desc.addVPSet("vertexCollections", vpsd1, temp1);
389  }
390  desc.addUntracked<bool>("verbose", false);
391  {
394  psd0.add<int>("numTracksThreshold", 0); // HI only
395  psd0.add<int>("maxNumTracksThreshold", 10000000); // HI only
396  psd0.add<double>("minPtTight", 0.0); // HI only
397  desc.add<edm::ParameterSetDescription>("TkFilterParameters", psd0);
398  }
399  desc.add<edm::InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
400  desc.add<edm::InputTag>("TrackLabel", edm::InputTag("generalTracks"));
401  desc.add<edm::InputTag>("TrackTimeResosLabel", edm::InputTag("dummy_default")); // 4D only
402  desc.add<edm::InputTag>("TrackTimesLabel", edm::InputTag("dummy_default")); // 4D only
403 
404  {
406  {
409  psd0.add<edm::ParameterSetDescription>("TkDAClusParameters", psd1);
410 
413  psd0.add<edm::ParameterSetDescription>("TkGapClusParameters", psd2);
414  }
415  psd0.add<std::string>("algorithm", "DA_vect");
416  desc.add<edm::ParameterSetDescription>("TkClusParameters", psd0);
417  }
418 
419  desc.add<bool>("isRecoveryIteration", false);
420  desc.add<edm::InputTag>("recoveryVtxCollection", {""});
421 
422  descriptions.add("primaryVertexProducer", desc);
423 }
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 123 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(), AlCaHLTBitMon_ParallelJobs::p, SiStripOfflineCRack_cfg::p2, hiPixelPairStep_cff::points, edm::Handle< T >::product(), FSQDQM_cfi::pvs, recoveryVtxToken, mps_fire::result, TrackFilterForPVFindingBase::select(), jetUpdater_cfi::sort, WeightedMeanFitter::startError, theTrackClusterizer, theTrackFilter, theTTBToken, protons_cff::time, trkTimeResosToken, trkTimesToken, trkToken, findQualityFiles::v, trackerHitRTTI::vector, w(), WeightedMeanFitter::weightedMeanOutlierRejection(), WeightedMeanFitter::weightedMeanOutlierRejectionBeamSpot(), and weightFit.

123  {
124  // get the BeamSpot, it will always be needed, even when not used as a constraint
126  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
127  iEvent.getByToken(bsToken, recoBeamSpotHandle);
128  if (recoBeamSpotHandle.isValid()) {
129  beamSpot = *recoBeamSpotHandle;
130  } else {
131  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
132  }
133 
134  bool validBS = true;
135  VertexState beamVertexState(beamSpot);
136  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
137  (beamVertexState.error().czz() <= 0.)) {
138  validBS = false;
139  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
140  }
141 
142  //if this is a recovery iteration, check if we already have a valid PV
143  if (fRecoveryIteration) {
144  auto const& oldVertices = iEvent.get(recoveryVtxToken);
145  //look for the first valid (not-BeamSpot) vertex
146  for (auto const& old : oldVertices) {
147  if (!(old.isFake())) {
148  //found a valid vertex, write the first one to the collection and return
149  //otherwise continue with regular vertexing procedure
150  auto result = std::make_unique<reco::VertexCollection>();
151  result->push_back(old);
152  iEvent.put(std::move(result), algorithms.begin()->label);
153  return;
154  }
155  }
156  }
157 
158  // get RECO tracks from the event
159  // `tks` can be used as a ptr to a reco::TrackCollection
161  iEvent.getByToken(trkToken, tks);
162 
163  // interface RECO tracks to vertex reconstruction
164  const auto& theB = &iSetup.getData(theTTBToken);
165  std::vector<reco::TransientTrack> t_tks;
166 
167  if (f4D) {
169  edm::Handle<edm::ValueMap<float>> trackTimeResosH;
170  iEvent.getByToken(trkTimesToken, trackTimesH);
171  iEvent.getByToken(trkTimeResosToken, trackTimeResosH);
172  t_tks = (*theB).build(tks, beamSpot, *(trackTimesH.product()), *(trackTimeResosH.product()));
173  } else {
174  t_tks = (*theB).build(tks, beamSpot);
175  }
176  if (fVerbose) {
177  std::cout << "RecoVertex/PrimaryVertexProducer"
178  << "Found: " << t_tks.size() << " reconstructed tracks"
179  << "\n";
180  }
181 
182  // select tracks
183  std::vector<reco::TransientTrack>&& seltks = theTrackFilter->select(t_tks);
184 
185  // clusterize tracks in Z
186  std::vector<std::vector<reco::TransientTrack>>&& clusters = theTrackClusterizer->clusterize(seltks);
187 
188  if (fVerbose) {
189  std::cout << " clustering returned " << clusters.size() << " clusters from " << seltks.size()
190  << " selected tracks" << std::endl;
191  }
192 
193  // vertex fits
194  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
195  auto result = std::make_unique<reco::VertexCollection>();
196  reco::VertexCollection& vColl = (*result);
197 
198  std::vector<TransientVertex> pvs;
199  for (std::vector<std::vector<reco::TransientTrack>>::const_iterator iclus = clusters.begin();
200  iclus != clusters.end();
201  iclus++) {
202  double sumwt = 0.;
203  double sumwt2 = 0.;
204  double sumw = 0.;
205  double meantime = 0.;
206  double vartime = 0.;
207  if (f4D) {
208  for (const auto& tk : *iclus) {
209  const double time = tk.timeExt();
210  const double err = tk.dtErrorExt();
211  const double inverr = err > 0. ? 1.0 / err : 0.;
212  const double w = inverr * inverr;
213  sumwt += w * time;
214  sumwt2 += w * time * time;
215  sumw += w;
216  }
217  meantime = sumwt / sumw;
218  double sumsq = sumwt2 - sumwt * sumwt / sumw;
219  double chisq = iclus->size() > 1 ? sumsq / double(iclus->size() - 1) : sumsq / double(iclus->size());
220  vartime = chisq / sumw;
221  }
222 
224  if (algorithm->fitter) {
225  if (algorithm->useBeamConstraint && validBS && (iclus->size() > 1)) {
226  v = algorithm->fitter->vertex(*iclus, beamSpot);
227  } else if (!(algorithm->useBeamConstraint) && (iclus->size() > 1)) {
228  v = algorithm->fitter->vertex(*iclus);
229  } // else: no fit ==> v.isValid()=False
230  } else if (weightFit) {
231  std::vector<std::pair<GlobalPoint, GlobalPoint>> points;
232  if (algorithm->useBeamConstraint && validBS && (iclus->size() > 1)) {
233  for (const auto& itrack : *iclus) {
234  GlobalPoint p = itrack.stateAtBeamLine().trackStateAtPCA().position();
235  GlobalPoint err(itrack.stateAtBeamLine().transverseImpactParameter().error(),
236  itrack.stateAtBeamLine().transverseImpactParameter().error(),
237  itrack.track().dzError());
238  std::pair<GlobalPoint, GlobalPoint> p2(p, err);
239  points.push_back(p2);
240  }
241 
243  if ((v.positionError().matrix())(2, 2) != (WeightedMeanFitter::startError * WeightedMeanFitter::startError))
244  pvs.push_back(v);
245  } else if (!(algorithm->useBeamConstraint) && (iclus->size() > 1)) {
246  for (const auto& itrack : *iclus) {
247  GlobalPoint p = itrack.impactPointState().globalPosition();
248  GlobalPoint err(itrack.track().dxyError(), itrack.track().dxyError(), itrack.track().dzError());
249  std::pair<GlobalPoint, GlobalPoint> p2(p, err);
250  points.push_back(p2);
251  }
252 
254  if ((v.positionError().matrix())(2, 2) != (WeightedMeanFitter::startError * WeightedMeanFitter::startError))
255  pvs.push_back(v); //FIX with constants
256  }
257  } else
258  throw VertexException(
259  "PrimaryVertexProducer: Something went wrong. You are not using the weighted mean fit and no algorithm was "
260  "selected.");
261 
262  // 4D vertices: add timing information
263  if (f4D and v.isValid()) {
264  auto err = v.positionError().matrix4D();
265  err(3, 3) = vartime;
266  auto trkWeightMap3d = v.weightMap(); // copy the 3d-fit weights
267  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared(), v.degreesOfFreedom());
268  v.weightMap(trkWeightMap3d);
269  }
270 
271  if (fVerbose) {
272  if (v.isValid()) {
273  std::cout << "x,y,z";
274  if (f4D)
275  std::cout << ",t";
276  std::cout << "=" << v.position().x() << " " << v.position().y() << " " << v.position().z();
277  if (f4D)
278  std::cout << " " << v.time();
279  std::cout << " cluster size = " << (*iclus).size() << std::endl;
280  } else {
281  std::cout << "Invalid fitted vertex, cluster size=" << (*iclus).size() << std::endl;
282  }
283  }
284 
285  //for weightFit we have already pushed it above (no timing infomration anyway)
286  if (v.isValid() && not weightFit && (v.degreesOfFreedom() >= algorithm->minNdof) &&
287  (!validBS || (*(algorithm->vertexSelector))(v, beamVertexState)))
288  pvs.push_back(v);
289  } // end of cluster loop
290 
291  if (fVerbose) {
292  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
293  }
294 
295  if (clusters.size() > 2 && clusters.size() > 2 * pvs.size())
296  edm::LogWarning("PrimaryVertexProducer")
297  << "more than half of candidate vertices lost " << pvs.size() << ' ' << clusters.size();
298 
299  if (pvs.empty() && seltks.size() > 5)
300  edm::LogWarning("PrimaryVertexProducer")
301  << "no vertex found with " << seltks.size() << " tracks and " << clusters.size() << " vertex-candidates";
302 
303  // sort vertices by pt**2 vertex (aka signal vertex tagging)
304  if (pvs.size() > 1) {
305  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
306  }
307 
308  // convert transient vertices returned by the theAlgo to (reco) vertices
309  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin(); iv != pvs.end(); iv++) {
310  reco::Vertex v = *iv;
311  vColl.push_back(v);
312  }
313 
314  if (vColl.empty()) {
315  GlobalError bse(beamSpot.rotatedCovariance3D());
316  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
318  we(0, 0) = 10000;
319  we(1, 1) = 10000;
320  we(2, 2) = 10000;
321  vColl.push_back(reco::Vertex(beamSpot.position(), we, 0., 0., 0));
322  if (fVerbose) {
323  std::cout << "RecoVertex/PrimaryVertexProducer: "
324  << "Beamspot with invalid errors " << bse.matrix() << std::endl;
325  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
326  }
327  } else {
328  vColl.push_back(reco::Vertex(beamSpot.position(), beamSpot.rotatedCovariance3D(), 0., 0., 0));
329  if (fVerbose) {
330  std::cout << "RecoVertex/PrimaryVertexProducer: "
331  << " will put Vertex derived from BeamSpot into Event.\n";
332  }
333  }
334  }
335 
336  if (fVerbose) {
337  int ivtx = 0;
338  for (reco::VertexCollection::const_iterator v = vColl.begin(); v != vColl.end(); ++v) {
339  std::cout << "recvtx " << ivtx++ << "#trk " << std::setw(3) << v->tracksSize() << " chi2 " << std::setw(4)
340  << v->chi2() << " ndof " << std::setw(3) << v->ndof() << " x " << std::setw(6) << v->position().x()
341  << " dx " << std::setw(6) << v->xError() << " y " << std::setw(6) << v->position().y() << " dy "
342  << std::setw(6) << v->yError() << " z " << std::setw(6) << v->position().z() << " dz " << std::setw(6)
343  << v->zError();
344  if (f4D) {
345  std::cout << " t " << std::setw(6) << v->t() << " dt " << std::setw(6) << v->tError();
346  }
347  std::cout << std::endl;
348  }
349  }
350 
351  iEvent.put(std::move(result), algorithm->label);
352  }
353 }
edm::EDGetTokenT< reco::BeamSpot > bsToken
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
int32_t *__restrict__ iv
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > theTTBToken
std::vector< algo > algorithms
T w() const
Common base class.
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
TransientVertex weightedMeanOutlierRejectionBeamSpot(const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< reco::TransientTrack > iclus, const reco::BeamSpot &beamSpot)
TrackFilterForPVFindingBase * theTrackFilter
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
TransientVertex weightedMeanOutlierRejection(const std::vector< std::pair< GlobalPoint, GlobalPoint >> &points, std::vector< reco::TransientTrack > iclus)
constexpr float startError

Member Data Documentation

◆ algorithms

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

Definition at line 85 of file PrimaryVertexProducer.h.

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

◆ bsToken

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

Definition at line 93 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ f4D

bool PrimaryVertexProducer::f4D
private

Definition at line 98 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ fRecoveryIteration

bool PrimaryVertexProducer::fRecoveryIteration
private

Definition at line 90 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ fVerbose

bool PrimaryVertexProducer::fVerbose
private

Definition at line 88 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ recoveryVtxToken

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

Definition at line 91 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ theConfig

edm::ParameterSet PrimaryVertexProducer::theConfig
private

Definition at line 87 of file PrimaryVertexProducer.h.

Referenced by config().

◆ theTrackClusterizer

TrackClusterizerInZ* PrimaryVertexProducer::theTrackClusterizer
private

Definition at line 74 of file PrimaryVertexProducer.h.

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

◆ theTrackFilter

TrackFilterForPVFindingBase* PrimaryVertexProducer::theTrackFilter
private

Definition at line 73 of file PrimaryVertexProducer.h.

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

◆ theTTBToken

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

Definition at line 71 of file PrimaryVertexProducer.h.

Referenced by produce().

◆ trkTimeResosToken

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

Definition at line 96 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ trkTimesToken

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

Definition at line 95 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ trkToken

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

Definition at line 94 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ weightFit

bool PrimaryVertexProducer::weightFit
private

Definition at line 99 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().