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 fRecoveryIteration
 
bool fVerbose
 
double minTrackTimeQuality_
 
edm::EDGetTokenT< reco::VertexCollectionrecoveryVtxToken
 
edm::ParameterSet theConfig
 
TrackClusterizerInZtheTrackClusterizer
 
TrackFilterForPVFindingBasetheTrackFilter
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordtheTTBToken
 
edm::ValueMap< float > trackMTDTimeQualities_
 
edm::EDGetTokenT< edm::ValueMap< float > > trackMTDTimeQualityToken
 
edm::ValueMap< float > trackTimes_
 
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
 
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
 
edm::EDGetTokenT< reco::TrackCollectiontrkToken
 
bool useMVASelection_
 
bool useTransientTrackTime_
 

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 63 of file PrimaryVertexProducer.h.

Constructor & Destructor Documentation

◆ PrimaryVertexProducer()

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

Definition at line 19 of file PrimaryVertexProducer.cc.

References qcdUeDQM_cfi::algorithm, algorithms, bsToken, fRecoveryIteration, fVerbose, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), minTrackTimeQuality_, recoveryVtxToken, AlCaHLTBitMon_QueryRunRegistry::string, theTrackClusterizer, theTrackFilter, trackMTDTimeQualityToken, trkTimeResosToken, trkTimesToken, trkToken, useMVASelection_, useTransientTrackTime_, and HLT_2024v10_cff::vertexCollections.

20  : theTTBToken(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))), theConfig(conf) {
21  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
22  useMVASelection_ = conf.getParameter<bool>("useMVACut");
23 
24  trkToken = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackLabel"));
25  bsToken = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotLabel"));
26  useTransientTrackTime_ = 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("PrimaryVertexProducer: unknown track selection algorithm: " + trackSelectionAlgorithm);
37  }
38 
39  // select and configure the track clusterizer
40  std::string clusteringAlgorithm =
41  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
42  if (clusteringAlgorithm == "gap") {
44  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
45  } else if (clusteringAlgorithm == "DA_vect") {
47  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
48  } else if (clusteringAlgorithm == "DA2D_vect") {
50  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
52  } else {
53  throw VertexException("PrimaryVertexProducer: unknown clustering algorithm: " + clusteringAlgorithm);
54  }
55 
57  trkTimesToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimesLabel"));
58  trkTimeResosToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimeResosLabel"));
60  consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("trackMTDTimeQualityVMapTag"));
61  minTrackTimeQuality_ = conf.getParameter<double>("minTrackTimeQuality");
62  }
63 
64  // select and configure the vertex fitters
65  std::vector<edm::ParameterSet> vertexCollections =
66  conf.getParameter<std::vector<edm::ParameterSet> >("vertexCollections");
67 
68  for (std::vector<edm::ParameterSet>::const_iterator algoconf = vertexCollections.begin();
69  algoconf != vertexCollections.end();
70  algoconf++) {
72 
73  algorithm.label = algoconf->getParameter<std::string>("label");
74 
75  // configure the fitter and selector
76  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
77  if (fitterAlgorithm == "KalmanVertexFitter") {
79  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
80  auto fitter = new AdaptiveVertexFitter(GeometricAnnealing(algoconf->getParameter<double>("chi2cutoff")));
81  algorithm.pv_fitter = new SequentialPrimaryVertexFitterAdapter(fitter);
82  } else if (fitterAlgorithm.empty()) {
83  algorithm.pv_fitter = nullptr;
84  } else if (fitterAlgorithm == "AdaptiveChisquareVertexFitter") {
85  algorithm.pv_fitter = new AdaptiveChisquarePrimaryVertexFitter(algoconf->getParameter<double>("chi2cutoff"),
86  algoconf->getParameter<double>("zcutoff"),
87  algoconf->getParameter<double>("mintrkweight"),
88  false);
89  } else if (fitterAlgorithm == "MultiPrimaryVertexFitter") {
90  algorithm.pv_fitter = new AdaptiveChisquarePrimaryVertexFitter(algoconf->getParameter<double>("chi2cutoff"),
91  algoconf->getParameter<double>("zcutoff"),
92  algoconf->getParameter<double>("mintrkweight"),
93  true);
94  } else if (fitterAlgorithm == "WeightedMeanFitter") {
96  } else {
97  throw VertexException("PrimaryVertexProducer: unknown algorithm: " + fitterAlgorithm);
98  }
99  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
100  algorithm.useBeamConstraint = algoconf->getParameter<bool>("useBeamConstraint");
101  algorithm.vertexSelector =
102  new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
103 
104  // configure separate vertex time reconstruction if applicable
105  // note that the vertex time could, in principle, also come from the clusterizer or the vertex fit
106 
107  const auto& pv_time_conf = algoconf->getParameter<edm::ParameterSet>("vertexTimeParameters");
108  const std::string vertexTimeAlgorithm = pv_time_conf.getParameter<std::string>("algorithm");
109  edm::ConsumesCollector&& collector = consumesCollector();
110 
111  if (vertexTimeAlgorithm.empty()) {
112  algorithm.pv_time_estimator = nullptr;
113  } else if (vertexTimeAlgorithm == "legacy4D") {
114  useTransientTrackTime_ = true;
115  algorithm.pv_time_estimator =
116  new VertexTimeAlgorithmLegacy4D(pv_time_conf.getParameter<edm::ParameterSet>("legacy4D"), collector);
117  } else if (vertexTimeAlgorithm == "fromTracksPID") {
118  algorithm.pv_time_estimator = new VertexTimeAlgorithmFromTracksPID(
119  pv_time_conf.getParameter<edm::ParameterSet>("fromTracksPID"), collector);
120  } else {
121  edm::LogWarning("MisConfiguration") << "unknown vertexTimeParameters.algorithm" << vertexTimeAlgorithm;
122  }
123  algorithms.push_back(algorithm);
124 
125  produces<reco::VertexCollection>(algorithm.label);
126  }
127 
128  //check if this is a recovery iteration
129  fRecoveryIteration = conf.getParameter<bool>("isRecoveryIteration");
130  if (fRecoveryIteration) {
131  if (algorithms.empty()) {
132  throw VertexException("PrimaryVertexProducer: No algorithm specified. ");
133  } else if (algorithms.size() > 1) {
134  throw VertexException(
135  "PrimaryVertexProducer: Running in Recovery mode and more than one algorithm specified. Please "
136  "only one algorithm.");
137  }
138  recoveryVtxToken = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("recoveryVtxCollection"));
139  }
140 }
edm::EDGetTokenT< edm::ValueMap< float > > trackMTDTimeQualityToken
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
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
Log< level::Warning, false > LogWarning

◆ ~PrimaryVertexProducer()

PrimaryVertexProducer::~PrimaryVertexProducer ( )
override

Definition at line 142 of file PrimaryVertexProducer.cc.

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

142  {
143  if (theTrackFilter)
144  delete theTrackFilter;
146  delete theTrackClusterizer;
147  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
148  if (algorithm->pv_fitter)
149  delete algorithm->pv_fitter;
150  if (algorithm->pv_time_estimator)
151  delete algorithm->pv_time_estimator;
152  if (algorithm->vertexSelector)
153  delete algorithm->vertexSelector;
154  }
155 }
std::vector< algo > algorithms
TrackClusterizerInZ * theTrackClusterizer
TrackFilterForPVFindingBase * theTrackFilter

Member Function Documentation

◆ config()

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

Definition at line 73 of file PrimaryVertexProducer.h.

References theConfig.

73 { return theConfig; }
edm::ParameterSet theConfig

◆ fillDescriptions()

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

Definition at line 359 of file PrimaryVertexProducer.cc.

References edm::ParameterSetDescription::add(), edm::ParameterSet::addParameter(), edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, VertexTimeAlgorithmFromTracksPID::fillPSetDescription(), VertexTimeAlgorithmLegacy4D::fillPSetDescription(), GapClusterizerInZ::fillPSetDescription(), DAClusterizerInZ_vect::fillPSetDescription(), HITrackFilterForPVFinding::fillPSetDescription(), DAClusterizerInZT_vect::fillPSetDescription(), edm::ParameterSetDescription::ifValue(), ProducerED_cfi::InputTag, or, edm::ParameterSetDescription::reserve(), and AlCaHLTBitMon_QueryRunRegistry::string.

359  {
360  edm::ParameterSetDescription psd_pv_time;
361  {
364  psd_pv_time.add<edm::ParameterSetDescription>("fromTracksPID", psd1);
365 
368  psd_pv_time.add<edm::ParameterSetDescription>("legacy4D", psd2);
369  }
370  psd_pv_time.add<std::string>("algorithm", ""); // default = none
371 
372  // vertex collections
374  {
376  vpsd1.add<double>("maxDistanceToBeam", 1.0);
377  vpsd1.add<std::string>("algorithm", "AdaptiveVertexFitter");
378  vpsd1.add<bool>("useBeamConstraint", false);
379  vpsd1.add<std::string>("label", "");
380  vpsd1.add<double>("chi2cutoff", 2.5);
381  vpsd1.add<double>("zcutoff", 1.0);
382  vpsd1.add<double>("mintrkweight", 0.0);
383  vpsd1.add<double>("minNdof", 0.0);
384  vpsd1.add<edm::ParameterSetDescription>("vertexTimeParameters", psd_pv_time);
385 
386  // two default values : with- and without beam constraint
387  std::vector<edm::ParameterSet> temp1;
388  temp1.reserve(2);
389  {
390  edm::ParameterSet temp2;
391  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
392  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
393  temp2.addParameter<bool>("useBeamConstraint", false);
394  temp2.addParameter<std::string>("label", "");
395  temp2.addParameter<double>("chi2cutoff", 2.5);
396  temp2.addParameter<double>("zcutoff", 1.0);
397  temp2.addParameter<double>("mintrkweight", 0.);
398  temp2.addParameter<double>("minNdof", 0.0);
399  edm::ParameterSet temp_vertexTime;
400  temp_vertexTime.addParameter<std::string>("algorithm", "");
401  temp2.addParameter<edm::ParameterSet>("vertexTimeParameters", temp_vertexTime);
402  temp1.push_back(temp2);
403  }
404  {
405  edm::ParameterSet temp2;
406  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
407  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
408  temp2.addParameter<bool>("useBeamConstraint", true);
409  temp2.addParameter<std::string>("label", "WithBS");
410  temp2.addParameter<double>("chi2cutoff", 2.5);
411  temp2.addParameter<double>("zcutoff", 1.0);
412  temp2.addParameter<double>("mintrkweight", 0.);
413  temp2.addParameter<double>("minNdof", 2.0);
414  edm::ParameterSet temp_vertexTime;
415  temp_vertexTime.addParameter<std::string>("algorithm", "");
416  temp2.addParameter<edm::ParameterSet>("vertexTimeParameters", temp_vertexTime);
417  temp1.push_back(temp2);
418  }
419  desc.addVPSet("vertexCollections", vpsd1, temp1);
420  }
421  desc.addUntracked<bool>("verbose", false);
422  {
424  HITrackFilterForPVFinding::fillPSetDescription(psd0); // extension of TrackFilterForPVFinding
425  desc.add<edm::ParameterSetDescription>("TkFilterParameters", psd0);
426  }
427  desc.add<edm::InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
428  desc.add<edm::InputTag>("TrackLabel", edm::InputTag("generalTracks"));
429  desc.add<edm::InputTag>("TrackTimeResosLabel", edm::InputTag("dummy_default")); // 4D only
430  desc.add<edm::InputTag>("TrackTimesLabel", edm::InputTag("dummy_default")); // 4D only
431  desc.add<edm::InputTag>("trackMTDTimeQualityVMapTag", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA")); // 4D only
432 
433  {
435  {
438 
441 
444 
445  psd0.ifValue(
446  edm::ParameterDescription<std::string>("algorithm", "DA_vect", true),
447  "DA_vect" >> edm::ParameterDescription<edm::ParameterSetDescription>("TkDAClusParameters", psd1, true) or
448  "DA2D_vect" >>
449  edm::ParameterDescription<edm::ParameterSetDescription>("TkDAClusParameters", psd2, true) or
450  "gap" >> edm::ParameterDescription<edm::ParameterSetDescription>("TkGapClusParameters", psd3, true));
451  }
452  desc.add<edm::ParameterSetDescription>("TkClusParameters", psd0);
453  }
454 
455  desc.add<bool>("isRecoveryIteration", false);
456  desc.add<edm::InputTag>("recoveryVtxCollection", {""});
457  desc.add<bool>("useMVACut", false);
458  desc.add<double>("minTrackTimeQuality", 0.8);
459 
460  descriptions.addWithDefaultLabel(desc);
461 }
ParameterDescriptionNode * ifValue(ParameterDescription< T > const &switchParameter, std::unique_ptr< ParameterDescriptionCases< T >> cases)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
void reserve(SetDescriptionEntries::size_type n)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:136
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
static void fillPSetDescription(edm::ParameterSetDescription &desc)

◆ produce()

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

Definition at line 157 of file PrimaryVertexProducer.cc.

References qcdUeDQM_cfi::algorithm, algorithms, pwdgSkimBPark_cfi::beamSpot, bsToken, bsc_activity_cfg::clusters, gather_cfg::cout, GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::cyy(), GlobalErrorBase< T, ErrorWeightType >::czz(), VertexState::error(), alignBH_cfg::fixed, fRecoveryIteration, fVerbose, edm::EventSetup::getData(), mps_fire::i, iEvent, edm::HandleBase::isValid(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::iv, GlobalErrorBase< T, ErrorWeightType >::matrix(), SiStripPI::max, minTrackTimeQuality_, eostools::move(), PrimaryVertexProducer::algo::pv_time_estimator, FSQDQM_cfi::pvs, recoveryVtxToken, mps_fire::result, TrackFilterForPVFindingBase::select(), VertexTimeAlgorithmBase::setEvent(), jetUpdater_cfi::sort, theTrackClusterizer, theTrackFilter, theTTBToken, trackMTDTimeQualities_, trackMTDTimeQualityToken, trackTimes_, trkTimeResosToken, trkTimesToken, trkToken, useMVASelection_, useTransientTrackTime_, findQualityFiles::v, and TrackClusterizerInZ::vertices().

157  {
158  // get the BeamSpot, it will always be needed, even when not used as a constraint
160  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
161  iEvent.getByToken(bsToken, recoBeamSpotHandle);
162  if (recoBeamSpotHandle.isValid()) {
163  beamSpot = *recoBeamSpotHandle;
164  } else {
165  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
166  }
167 
168  bool validBS = true;
169  VertexState beamVertexState(beamSpot);
170  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
171  (beamVertexState.error().czz() <= 0.)) {
172  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
173  validBS = false;
174  }
175 
176  //if this is a recovery iteration, check if we already have a valid PV
177  if (fRecoveryIteration) {
178  auto const& oldVertices = iEvent.get(recoveryVtxToken);
179  //look for the first valid (not-BeamSpot) vertex
180  for (auto const& old : oldVertices) {
181  if (!(old.isFake())) {
182  //found a valid vertex, write the first one to the collection and return
183  //otherwise continue with regular vertexing procedure
184  auto result = std::make_unique<reco::VertexCollection>();
185  result->push_back(old);
186  iEvent.put(std::move(result), algorithms.begin()->label);
187  return;
188  }
189  }
190  }
191 
192  // get RECO tracks from the event
193  // `tks` can be used as a ptr to a reco::TrackCollection
195  iEvent.getByToken(trkToken, tks);
196 
197  // mechanism to put the beamspot if the track collection is empty
198  if (!tks.isValid()) {
199  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
200  auto result = std::make_unique<reco::VertexCollection>();
201  reco::VertexCollection& vColl = (*result);
202 
203  GlobalError bse(beamSpot.rotatedCovariance3D());
204  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
206  we(0, 0) = 10000;
207  we(1, 1) = 10000;
208  we(2, 2) = 10000;
209  vColl.push_back(reco::Vertex(beamSpot.position(), we, 0., 0., 0));
210  if (fVerbose) {
211  std::cout << "RecoVertex/PrimaryVertexProducer: "
212  << "Beamspot with invalid errors " << bse.matrix() << std::endl;
213  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
214  }
215  } else {
216  vColl.push_back(reco::Vertex(beamSpot.position(), beamSpot.rotatedCovariance3D(), 0., 0., 0));
217  if (fVerbose) {
218  std::cout << "RecoVertex/PrimaryVertexProducer: "
219  << " will put Vertex derived from BeamSpot into Event.\n";
220  }
221  }
222  iEvent.put(std::move(result), algorithm->label);
223  }
224 
225  return; // early return
226  }
227 
228  for (auto& algo : algorithms) {
229  if (algo.pv_time_estimator) {
230  algo.pv_time_estimator->setEvent(iEvent, iSetup);
231  }
232  }
233 
234  // interface RECO tracks to vertex reconstruction
235  const auto& theB = &iSetup.getData(theTTBToken);
236  std::vector<reco::TransientTrack> t_tks;
237 
239  auto const& trackTimeResos_ = iEvent.get(trkTimeResosToken);
240  auto trackTimes_ = iEvent.get(trkTimesToken);
241 
242  if (useMVASelection_) {
244 
245  for (unsigned int i = 0; i < (*tks).size(); i++) {
246  const reco::TrackRef ref(tks, i);
247  auto const trkTimeQuality = trackMTDTimeQualities_[ref];
248  if (trkTimeQuality < minTrackTimeQuality_) {
250  }
251  }
252  t_tks = (*theB).build(tks, beamSpot, trackTimes_, trackTimeResos_);
253  } else {
254  t_tks = (*theB).build(tks, beamSpot, trackTimes_, trackTimeResos_);
255  }
256  } else {
257  t_tks = (*theB).build(tks, beamSpot);
258  }
259 
260  // select tracks
261  std::vector<reco::TransientTrack>&& seltks = theTrackFilter->select(t_tks);
262 
263  // clusterize tracks in Z
264  std::vector<TransientVertex>&& clusters = theTrackClusterizer->vertices(seltks);
265 
266  if (fVerbose) {
267  edm::LogPrint("PrimaryVertexProducer")
268  << "Clustering returned " << clusters.size() << " clusters from " << seltks.size() << " selected tracks";
269  }
270 
271  // vertex fits
272  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
273  auto result = std::make_unique<reco::VertexCollection>();
274  reco::VertexCollection& vColl = (*result);
275  std::vector<TransientVertex> pvs;
276 
277  if (algorithm->pv_fitter == nullptr) {
278  pvs = clusters;
279  } else {
280  pvs = algorithm->pv_fitter->fit(seltks, clusters, beamSpot, algorithm->useBeamConstraint);
281  }
282 
283  if (algorithm->pv_time_estimator != nullptr) {
284  algorithm->pv_time_estimator->fill_vertex_times(pvs);
285  }
286 
287  // sort vertices by pt**2 vertex
288  if (pvs.size() > 1) {
289  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
290  }
291 
292  // select and convert transient vertices to (reco) vertices
293  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin(); iv != pvs.end(); iv++) {
294  if (iv->isValid() && (iv->degreesOfFreedom() >= algorithm->minNdof)) {
295  reco::Vertex v = *iv;
296  if (!validBS || ((*(algorithm->vertexSelector))(v, beamVertexState))) {
297  vColl.push_back(v);
298  }
299  }
300  }
301 
302  if (fVerbose) {
303  edm::LogPrint("PrimaryVertexProducer") << "PrimaryVertexProducer \"" << algorithm->label << "\" contains "
304  << pvs.size() << " reco::Vertex candidates";
305  }
306 
307  if (clusters.size() > 2 && clusters.size() > 2 * pvs.size()) {
308  edm::LogWarning("PrimaryVertexProducer")
309  << "More than 50% of candidate vertices lost (" << pvs.size() << " out of " << clusters.size() << ")";
310  }
311 
312  if (pvs.empty() && seltks.size() > 5) {
313  edm::LogWarning("PrimaryVertexProducer")
314  << "No vertex found with " << seltks.size() << " tracks and " << clusters.size() << " vertex candidates";
315  }
316 
317  if (vColl.empty()) {
318  GlobalError bse(beamSpot.rotatedCovariance3D());
319  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
321  we(0, 0) = 10000;
322  we(1, 1) = 10000;
323  we(2, 2) = 10000;
324  vColl.push_back(reco::Vertex(beamSpot.position(), we, 0., 0., 0));
325  edm::LogWarning("PrimaryVertexProducer") << "Zero recostructed vertices, will put reco::Vertex derived from "
326  "dummy/fake BeamSpot into Event, BeamSpot has invalid errors: "
327  << bse.matrix();
328  } else {
329  vColl.push_back(reco::Vertex(beamSpot.position(), beamSpot.rotatedCovariance3D(), 0., 0., 0));
330  if (fVerbose) {
331  edm::LogWarning("PrimaryVertexProducer")
332  << "Zero recostructed vertices, will put reco::Vertex derived from BeamSpot into Event.";
333  }
334  }
335  }
336 
337  if (fVerbose) {
338  int ivtx = 0;
339  for (reco::VertexCollection::const_iterator v = vColl.begin(); v != vColl.end(); ++v) {
340  edm::LogPrint("PrimaryVertexProducer")
341  << "recvtx " << std::setw(3) << std::fixed << ivtx++ << " #trk " << std::setw(3) << v->tracksSize()
342  << " chi2 " << std::setw(5) << std::setprecision(1) << v->chi2() << " ndof " << std::setw(5)
343  << std::setprecision(1) << v->ndof() << " x " << std::setw(7) << std::setprecision(4) << v->position().x()
344  << " dx " << std::setw(6) << std::setprecision(4) << v->xError() << " y " << std::setw(7)
345  << std::setprecision(4) << v->position().y() << " dy " << std::setw(6) << std::setprecision(4)
346  << v->yError() << " z " << std::setw(8) << std::setprecision(4) << v->position().z() << " dz "
347  << std::setw(6) << std::setprecision(4) << v->zError();
348  if (v->tError() > 0) {
349  edm::LogPrint("PrimaryVertexProducer") << " t " << std::setw(6) << std::setprecision(3) << v->t() << " dt "
350  << std::setw(6) << std::setprecision(3) << v->tError();
351  }
352  }
353  }
354 
355  iEvent.put(std::move(result), algorithm->label);
356  }
357 }
edm::EDGetTokenT< edm::ValueMap< float > > trackMTDTimeQualityToken
edm::EDGetTokenT< reco::BeamSpot > bsToken
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > theTTBToken
std::vector< algo > algorithms
edm::EDGetTokenT< reco::VertexCollection > recoveryVtxToken
edm::EDGetTokenT< reco::TrackCollection > trkToken
TrackClusterizerInZ * theTrackClusterizer
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::ValueMap< float > trackTimes_
Log< level::Error, false > LogError
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
int iEvent
Definition: GenABIO.cc:224
edm::ValueMap< float > trackMTDTimeQualities_
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
TrackFilterForPVFindingBase * theTrackFilter
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
Log< level::Warning, true > LogPrint
bool isValid() const
Definition: HandleBase.h:70
virtual std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const =0
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ algorithms

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

Definition at line 92 of file PrimaryVertexProducer.h.

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

◆ bsToken

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

Definition at line 100 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ fRecoveryIteration

bool PrimaryVertexProducer::fRecoveryIteration
private

Definition at line 97 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ fVerbose

bool PrimaryVertexProducer::fVerbose
private

Definition at line 95 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ minTrackTimeQuality_

double PrimaryVertexProducer::minTrackTimeQuality_
private

Definition at line 110 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ recoveryVtxToken

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

Definition at line 98 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ theConfig

edm::ParameterSet PrimaryVertexProducer::theConfig
private

Definition at line 94 of file PrimaryVertexProducer.h.

Referenced by config().

◆ theTrackClusterizer

TrackClusterizerInZ* PrimaryVertexProducer::theTrackClusterizer
private

Definition at line 80 of file PrimaryVertexProducer.h.

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

◆ theTrackFilter

TrackFilterForPVFindingBase* PrimaryVertexProducer::theTrackFilter
private

Definition at line 79 of file PrimaryVertexProducer.h.

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

◆ theTTBToken

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

Definition at line 77 of file PrimaryVertexProducer.h.

Referenced by produce().

◆ trackMTDTimeQualities_

edm::ValueMap<float> PrimaryVertexProducer::trackMTDTimeQualities_
private

Definition at line 108 of file PrimaryVertexProducer.h.

Referenced by produce().

◆ trackMTDTimeQualityToken

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

Definition at line 104 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ trackTimes_

edm::ValueMap<float> PrimaryVertexProducer::trackTimes_
private

Definition at line 109 of file PrimaryVertexProducer.h.

Referenced by produce().

◆ trkTimeResosToken

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

Definition at line 103 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ trkTimesToken

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

Definition at line 102 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ trkToken

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

Definition at line 101 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ useMVASelection_

bool PrimaryVertexProducer::useMVASelection_
private

Definition at line 107 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

◆ useTransientTrackTime_

bool PrimaryVertexProducer::useTransientTrackTime_
private

Definition at line 106 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().