CMS 3D CMS Logo

PrimaryVertexProducer.cc
Go to the documentation of this file.
8 
12 
17 
19 
21  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
22 
23  trkToken = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackLabel"));
24  bsToken = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotLabel"));
25  f4D = 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("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " +
36  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") {
47  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") {
52  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
53  } else if (clusteringAlgorithm == "DA2D_vect") {
55  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
56  f4D = true;
57  }
58 
59  else {
60  throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);
61  }
62 
63  if (f4D) {
64  trkTimesToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimesLabel"));
65  trkTimeResosToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimeResosLabel"));
66  }
67 
68  // select and configure the vertex fitters
69  if (conf.exists("vertexCollections")) {
70  std::vector<edm::ParameterSet> vertexCollections =
71  conf.getParameter<std::vector<edm::ParameterSet> >("vertexCollections");
72 
73  for (std::vector<edm::ParameterSet>::const_iterator algoconf = vertexCollections.begin();
74  algoconf != vertexCollections.end();
75  algoconf++) {
77  std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
78  if (fitterAlgorithm == "KalmanVertexFitter") {
79  algorithm.fitter = new KalmanVertexFitter();
80  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
81  algorithm.fitter = new AdaptiveVertexFitter(GeometricAnnealing(algoconf->getParameter<double>("chi2cutoff")));
82  } else {
83  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
84  }
85  algorithm.label = algoconf->getParameter<std::string>("label");
86  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
87  algorithm.useBeamConstraint = algoconf->getParameter<bool>("useBeamConstraint");
88  algorithm.vertexSelector =
89  new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
90  algorithms.push_back(algorithm);
91 
92  produces<reco::VertexCollection>(algorithm.label);
93  }
94  } else {
95  edm::LogWarning("MisConfiguration")
96  << "this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
97 
99  std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
100  if (fitterAlgorithm == "KalmanVertexFitter") {
101  algorithm.fitter = new KalmanVertexFitter();
102  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
103  algorithm.fitter = new AdaptiveVertexFitter();
104  } else {
105  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
106  }
107  algorithm.label = "";
108  algorithm.minNdof = conf.getParameter<double>("minNdof");
109  algorithm.useBeamConstraint = conf.getParameter<bool>("useBeamConstraint");
110 
111  algorithm.vertexSelector = new VertexCompatibleWithBeam(
113  conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
114 
115  algorithms.push_back(algorithm);
116  produces<reco::VertexCollection>(algorithm.label);
117  }
118 
119  //check if this is a recovery iteration
120  fRecoveryIteration = conf.getParameter<bool>("isRecoveryIteration");
121  if (fRecoveryIteration) {
122  if (algorithms.empty()) {
123  throw VertexException("PrimaryVertexProducerAlgorithm: No algorithm specified. ");
124  } else if (algorithms.size() > 1) {
125  throw VertexException(
126  "PrimaryVertexProducerAlgorithm: Running in Recovery mode and more than one algorithm specified. Please "
127  "only one algorithm.");
128  }
129  recoveryVtxToken = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("recoveryVtxCollection"));
130  }
131 }
132 
134  if (theTrackFilter)
135  delete theTrackFilter;
137  delete theTrackClusterizer;
138  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
139  if (algorithm->fitter)
140  delete algorithm->fitter;
141  if (algorithm->vertexSelector)
142  delete algorithm->vertexSelector;
143  }
144 }
145 
147  // get the BeamSpot, it will alwys be needed, even when not used as a constraint
149  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
150  iEvent.getByToken(bsToken, recoBeamSpotHandle);
151  if (recoBeamSpotHandle.isValid()) {
152  beamSpot = *recoBeamSpotHandle;
153  } else {
154  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
155  }
156 
157  bool validBS = true;
158  VertexState beamVertexState(beamSpot);
159  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
160  (beamVertexState.error().czz() <= 0.)) {
161  validBS = false;
162  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
163  }
164 
165  //if this is a recovery iteration, check if we already have a valid PV
166  if (fRecoveryIteration) {
167  auto const& oldVertices = iEvent.get(recoveryVtxToken);
168  //look for the first valid (not-BeamSpot) vertex
169  for (auto const& old : oldVertices) {
170  if (!(old.isFake())) {
171  //found a valid vertex, write the first one to the collection and return
172  //otherwise continue with regular vertexing procedure
173  auto result = std::make_unique<reco::VertexCollection>();
174  result->push_back(old);
175  iEvent.put(std::move(result), algorithms.begin()->label);
176  return;
177  }
178  }
179  }
180 
181  // get RECO tracks from the event
182  // `tks` can be used as a ptr to a reco::TrackCollection
184  iEvent.getByToken(trkToken, tks);
185 
186  // interface RECO tracks to vertex reconstruction
188  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
189  std::vector<reco::TransientTrack> t_tks;
190 
191  if (f4D) {
192  edm::Handle<edm::ValueMap<float> > trackTimesH;
193  edm::Handle<edm::ValueMap<float> > trackTimeResosH;
194  iEvent.getByToken(trkTimesToken, trackTimesH);
195  iEvent.getByToken(trkTimeResosToken, trackTimeResosH);
196  t_tks = (*theB).build(tks, beamSpot, *(trackTimesH.product()), *(trackTimeResosH.product()));
197  } else {
198  t_tks = (*theB).build(tks, beamSpot);
199  }
200  if (fVerbose) {
201  std::cout << "RecoVertex/PrimaryVertexProducer"
202  << "Found: " << t_tks.size() << " reconstructed tracks"
203  << "\n";
204  }
205 
206  // select tracks
207  std::vector<reco::TransientTrack>&& seltks = theTrackFilter->select(t_tks);
208 
209  // clusterize tracks in Z
210  std::vector<std::vector<reco::TransientTrack> >&& clusters = theTrackClusterizer->clusterize(seltks);
211 
212  if (fVerbose) {
213  std::cout << " clustering returned " << clusters.size() << " clusters from " << seltks.size()
214  << " selected tracks" << std::endl;
215  }
216 
217  // vertex fits
218  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
219  auto result = std::make_unique<reco::VertexCollection>();
220  reco::VertexCollection& vColl = (*result);
221 
222  std::vector<TransientVertex> pvs;
223  for (std::vector<std::vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin();
224  iclus != clusters.end();
225  iclus++) {
226  double sumwt = 0.;
227  double sumwt2 = 0.;
228  double sumw = 0.;
229  double meantime = 0.;
230  double vartime = 0.;
231  if (f4D) {
232  for (const auto& tk : *iclus) {
233  const double time = tk.timeExt();
234  const double err = tk.dtErrorExt();
235  const double inverr = err > 0. ? 1.0 / err : 0.;
236  const double w = inverr * inverr;
237  sumwt += w * time;
238  sumwt2 += w * time * time;
239  sumw += w;
240  }
241  meantime = sumwt / sumw;
242  double sumsq = sumwt2 - sumwt * sumwt / sumw;
243  double chisq = iclus->size() > 1 ? sumsq / double(iclus->size() - 1) : sumsq / double(iclus->size());
244  vartime = chisq / sumw;
245  }
246 
248  if (algorithm->useBeamConstraint && validBS && ((*iclus).size() > 1)) {
249  v = algorithm->fitter->vertex(*iclus, beamSpot);
250 
251  if (f4D) {
252  if (v.isValid()) {
253  auto err = v.positionError().matrix4D();
254  err(3, 3) = vartime;
255  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared());
256  }
257  }
258 
259  } else if (!(algorithm->useBeamConstraint) && ((*iclus).size() > 1)) {
260  v = algorithm->fitter->vertex(*iclus);
261 
262  if (f4D) {
263  if (v.isValid()) {
264  auto err = v.positionError().matrix4D();
265  err(3, 3) = vartime;
266  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared());
267  }
268  }
269 
270  } // else: no fit ==> v.isValid()=False
271 
272  if (fVerbose) {
273  if (v.isValid()) {
274  std::cout << "x,y,z";
275  if (f4D)
276  std::cout << ",t";
277  std::cout << "=" << v.position().x() << " " << v.position().y() << " " << v.position().z();
278  if (f4D)
279  std::cout << " " << v.time();
280  std::cout << " cluster size = " << (*iclus).size() << std::endl;
281  } else {
282  std::cout << "Invalid fitted vertex, cluster size=" << (*iclus).size() << std::endl;
283  }
284  }
285 
286  if (v.isValid() && (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 }
354 
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  {
393  psd0.add<double>("maxNormalizedChi2", 10.0);
394  psd0.add<double>("minPt", 0.0);
395  psd0.add<std::string>("algorithm", "filter");
396  psd0.add<double>("maxEta", 2.4);
397  psd0.add<double>("maxD0Significance", 4.0);
398  psd0.add<double>("maxD0Error", 1.0);
399  psd0.add<double>("maxDzError", 1.0);
400  psd0.add<std::string>("trackQuality", "any");
401  psd0.add<int>("minPixelLayersWithHits", 2);
402  psd0.add<int>("minSiliconLayersWithHits", 5);
403  psd0.add<int>("numTracksThreshold", 0); // HI only
404  desc.add<edm::ParameterSetDescription>("TkFilterParameters", psd0);
405  }
406  desc.add<edm::InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
407  desc.add<edm::InputTag>("TrackLabel", edm::InputTag("generalTracks"));
408  desc.add<edm::InputTag>("TrackTimeResosLabel", edm::InputTag("dummy_default")); // 4D only
409  desc.add<edm::InputTag>("TrackTimesLabel", edm::InputTag("dummy_default")); // 4D only
410  {
412  {
414  psd1.addUntracked<bool>("verbose", false);
415  psd1.addUntracked<double>("zdumpcenter", 0.);
416  psd1.addUntracked<double>("zdumpwidth", 20.);
417  psd1.addUntracked<bool>("use_vdt", false); // obsolete, appears in HLT configs
418  psd1.add<double>("d0CutOff", 3.0);
419  psd1.add<double>("Tmin", 2.0);
420  psd1.add<double>("delta_lowT", 0.001);
421  psd1.add<double>("zmerge", 0.01);
422  psd1.add<double>("dzCutOff", 3.0);
423  psd1.add<double>("Tpurge", 2.0);
424  psd1.add<int>("convergence_mode", 0);
425  psd1.add<double>("delta_highT", 0.01);
426  psd1.add<double>("Tstop", 0.5);
427  psd1.add<double>("coolingFactor", 0.6);
428  psd1.add<double>("vertexSize", 0.006);
429  psd1.add<double>("uniquetrkweight", 0.8);
430  psd1.add<double>("zrange", 4.0);
431 
432  psd1.add<double>("tmerge", 0.01); // 4D only
433  psd1.add<double>("dtCutOff", 4.); // 4D only
434  psd1.add<double>("t0Max", 1.0); // 4D only
435  psd1.add<double>("vertexSizeTime", 0.008); // 4D only
436 
437  psd0.add<edm::ParameterSetDescription>("TkDAClusParameters", psd1);
438 
440  psd2.add<double>("zSeparation", 1.0);
441  psd0.add<edm::ParameterSetDescription>("TkGapClusParameters", psd2);
442  }
443  psd0.add<std::string>("algorithm", "DA_vect");
444  desc.add<edm::ParameterSetDescription>("TkClusParameters", psd0);
445  }
446  desc.add<bool>("isRecoveryIteration", false);
447  desc.add<edm::InputTag>("recoveryVtxCollection", {""});
448 
449  descriptions.add("primaryVertexProducer", desc);
450 }
451 
452 //define this as a plug-in
AdaptiveVertexFitter
Definition: AdaptiveVertexFitter.h:29
TrackFilterForPVFindingBase::select
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
Handle.h
AlgebraicSymMatrix33
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Definition: AlgebraicROOTObjects.h:21
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
pwdgSkimBPark_cfi.beamSpot
beamSpot
Definition: pwdgSkimBPark_cfi.py:5
MessageLogger.h
edm::Handle::product
T const * product() const
Definition: Handle.h:70
ESHandle.h
VertexException
Common base class.
Definition: VertexException.h:12
reco::VertexCollection
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
gather_cfg.cout
cout
Definition: gather_cfg.py:144
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89353
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
PrimaryVertexProducer::trkTimeResosToken
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
Definition: PrimaryVertexProducer.h:91
PrimaryVertexProducer::recoveryVtxToken
edm::EDGetTokenT< reco::VertexCollection > recoveryVtxToken
Definition: PrimaryVertexProducer.h:86
GlobalErrorBase::matrix
const AlgebraicSymMatrix33 matrix() const
Definition: GlobalErrorBase.h:121
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
PrimaryVertexProducer
Definition: PrimaryVertexProducer.h:54
PrimaryVertexProducer::trkTimesToken
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
Definition: PrimaryVertexProducer.h:90
TransientTrack.h
findQualityFiles.v
v
Definition: findQualityFiles.py:179
edm::Handle< reco::BeamSpot >
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
PrimaryVertexProducer::theTrackClusterizer
TrackClusterizerInZ * theTrackClusterizer
Definition: PrimaryVertexProducer.h:69
PrimaryVertexProducer::algo
Definition: PrimaryVertexProducer.h:72
PrimaryVertexProducer.h
MakerMacros.h
GlobalErrorBase::cyy
T cyy() const
Definition: GlobalErrorBase.h:101
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
VertexDistanceXY.h
BeamSpot.h
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
w
const double w
Definition: UKUtility.cc:23
PrimaryVertexProducer::fRecoveryIteration
bool fRecoveryIteration
Definition: PrimaryVertexProducer.h:85
reco::BeamSpot
Definition: BeamSpot.h:21
TransientTrackRecord
Definition: TransientTrackRecord.h:11
PrimaryVertexProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PrimaryVertexProducer.cc:355
edm::ESHandle< TransientTrackBuilder >
VertexState::error
GlobalError error() const
Definition: VertexState.h:64
PrimaryVertexProducer::theTrackFilter
TrackFilterForPVFindingBase * theTrackFilter
Definition: PrimaryVertexProducer.h:68
GlobalErrorBase::cxx
T cxx() const
Definition: GlobalErrorBase.h:97
DAClusterizerInZT_vect
Definition: DAClusterizerInZT_vect.h:22
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
GeometricAnnealing.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
edm::ParameterSet::exists
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: ParameterSet.cc:681
edm::ParameterSetDescription::addUntracked
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:100
TransientTrackBuilder.h
PrimaryVertexProducer::f4D
bool f4D
Definition: PrimaryVertexProducer.h:93
edm::ParameterSet
Definition: ParameterSet.h:47
TrackClusterizerInZ::clusterize
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
edm::ParameterSet::addParameter
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
qcdUeDQM_cfi.algorithm
algorithm
Definition: qcdUeDQM_cfi.py:32
iEvent
int iEvent
Definition: GenABIO.cc:224
VertexCompatibleWithBeam
Definition: VertexCompatibleWithBeam.h:15
GeometricAnnealing
Definition: GeometricAnnealing.h:7
DAClusterizerInZ
Definition: DAClusterizerInZ.h:18
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
GlobalErrorBase< double, ErrorMatrixTag >
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
VertexDistanceXY
Definition: VertexDistanceXY.h:11
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
TransientTrackRecord.h
submitPVResolutionJobs.err
err
Definition: submitPVResolutionJobs.py:85
get
#define get
FSQDQM_cfi.pvs
pvs
Definition: FSQDQM_cfi.py:12
PrimaryVertexProducer::trkToken
edm::EDGetTokenT< reco::TrackCollection > trkToken
Definition: PrimaryVertexProducer.h:89
InputTag.h
VertexFwd.h
TransientVertex.h
DAClusterizerInZ_vect
Definition: DAClusterizerInZ_vect.h:20
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
PrimaryVertexProducer::fVerbose
bool fVerbose
Definition: PrimaryVertexProducer.h:83
HITrackFilterForPVFinding
Definition: HITrackFilterForPVFinding.h:13
PrimaryVertexProducer::~PrimaryVertexProducer
~PrimaryVertexProducer() override
Definition: PrimaryVertexProducer.cc:133
PrimaryVertexProducer::PrimaryVertexProducer
PrimaryVertexProducer(const edm::ParameterSet &)
Definition: PrimaryVertexProducer.cc:20
VertexState
Definition: VertexState.h:13
HLT_FULL_cff.vertexCollections
vertexCollections
Definition: HLT_FULL_cff.py:34864
PrimaryVertexProducer::bsToken
edm::EDGetTokenT< reco::BeamSpot > bsToken
Definition: PrimaryVertexProducer.h:88
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
GapClusterizerInZ
Definition: GapClusterizerInZ.h:14
GlobalErrorBase::czz
T czz() const
Definition: GlobalErrorBase.h:107
TrackFilterForPVFinding
Definition: TrackFilterForPVFinding.h:14
VertexHigherPtSquared
Definition: VertexHigherPtSquared.h:13
PrimaryVertexProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PrimaryVertexProducer.cc:146
PrimaryVertexProducer::algorithms
std::vector< algo > algorithms
Definition: PrimaryVertexProducer.h:80
mps_fire.result
result
Definition: mps_fire.py:311
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
ntuplemaker.time
time
Definition: ntuplemaker.py:310
edm::Event
Definition: Event.h:73
edm::Log
Definition: MessageLogger.h:70
edm::InputTag
Definition: InputTag.h:15
reco::Vertex
Definition: Vertex.h:35
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22