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 
121  if (theTrackFilter)
122  delete theTrackFilter;
124  delete theTrackClusterizer;
125  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
126  if (algorithm->fitter)
127  delete algorithm->fitter;
128  if (algorithm->vertexSelector)
129  delete algorithm->vertexSelector;
130  }
131 }
132 
134  // get the BeamSpot, it will alwys be needed, even when not used as a constraint
136  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
137  iEvent.getByToken(bsToken, recoBeamSpotHandle);
138  if (recoBeamSpotHandle.isValid()) {
139  beamSpot = *recoBeamSpotHandle;
140  } else {
141  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
142  }
143 
144  bool validBS = true;
145  VertexState beamVertexState(beamSpot);
146  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
147  (beamVertexState.error().czz() <= 0.)) {
148  validBS = false;
149  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
150  }
151 
152  // get RECO tracks from the event
153  // `tks` can be used as a ptr to a reco::TrackCollection
155  iEvent.getByToken(trkToken, tks);
156 
157  // interface RECO tracks to vertex reconstruction
159  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
160  std::vector<reco::TransientTrack> t_tks;
161 
162  if (f4D) {
163  edm::Handle<edm::ValueMap<float> > trackTimesH;
164  edm::Handle<edm::ValueMap<float> > trackTimeResosH;
165  iEvent.getByToken(trkTimesToken, trackTimesH);
166  iEvent.getByToken(trkTimeResosToken, trackTimeResosH);
167  t_tks = (*theB).build(tks, beamSpot, *(trackTimesH.product()), *(trackTimeResosH.product()));
168  } else {
169  t_tks = (*theB).build(tks, beamSpot);
170  }
171  if (fVerbose) {
172  std::cout << "RecoVertex/PrimaryVertexProducer"
173  << "Found: " << t_tks.size() << " reconstructed tracks"
174  << "\n";
175  }
176 
177  // select tracks
178  std::vector<reco::TransientTrack>&& seltks = theTrackFilter->select(t_tks);
179 
180  // clusterize tracks in Z
181  std::vector<std::vector<reco::TransientTrack> >&& clusters = theTrackClusterizer->clusterize(seltks);
182 
183  if (fVerbose) {
184  std::cout << " clustering returned " << clusters.size() << " clusters from " << seltks.size()
185  << " selected tracks" << std::endl;
186  }
187 
188  // vertex fits
189  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
190  auto result = std::make_unique<reco::VertexCollection>();
191  reco::VertexCollection& vColl = (*result);
192 
193  std::vector<TransientVertex> pvs;
194  for (std::vector<std::vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin();
195  iclus != clusters.end();
196  iclus++) {
197  double sumwt = 0.;
198  double sumwt2 = 0.;
199  double sumw = 0.;
200  double meantime = 0.;
201  double vartime = 0.;
202  if (f4D) {
203  for (const auto& tk : *iclus) {
204  const double time = tk.timeExt();
205  const double err = tk.dtErrorExt();
206  const double inverr = err > 0. ? 1.0 / err : 0.;
207  const double w = inverr * inverr;
208  sumwt += w * time;
209  sumwt2 += w * time * time;
210  sumw += w;
211  }
212  meantime = sumwt / sumw;
213  double sumsq = sumwt2 - sumwt * sumwt / sumw;
214  double chisq = iclus->size() > 1 ? sumsq / double(iclus->size() - 1) : sumsq / double(iclus->size());
215  vartime = chisq / sumw;
216  }
217 
219  if (algorithm->useBeamConstraint && validBS && ((*iclus).size() > 1)) {
220  v = algorithm->fitter->vertex(*iclus, beamSpot);
221 
222  if (f4D) {
223  if (v.isValid()) {
224  auto err = v.positionError().matrix4D();
225  err(3, 3) = vartime;
226  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared());
227  }
228  }
229 
230  } else if (!(algorithm->useBeamConstraint) && ((*iclus).size() > 1)) {
231  v = algorithm->fitter->vertex(*iclus);
232 
233  if (f4D) {
234  if (v.isValid()) {
235  auto err = v.positionError().matrix4D();
236  err(3, 3) = vartime;
237  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared());
238  }
239  }
240 
241  } // else: no fit ==> v.isValid()=False
242 
243  if (fVerbose) {
244  if (v.isValid()) {
245  std::cout << "x,y,z";
246  if (f4D)
247  std::cout << ",t";
248  std::cout << "=" << v.position().x() << " " << v.position().y() << " " << v.position().z();
249  if (f4D)
250  std::cout << " " << v.time();
251  std::cout << " cluster size = " << (*iclus).size() << std::endl;
252  } else {
253  std::cout << "Invalid fitted vertex, cluster size=" << (*iclus).size() << std::endl;
254  }
255  }
256 
257  if (v.isValid() && (v.degreesOfFreedom() >= algorithm->minNdof) &&
258  (!validBS || (*(algorithm->vertexSelector))(v, beamVertexState)))
259  pvs.push_back(v);
260  } // end of cluster loop
261 
262  if (fVerbose) {
263  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
264  }
265 
266  if (clusters.size() > 2 && clusters.size() > 2 * pvs.size())
267  edm::LogWarning("PrimaryVertexProducer")
268  << "more than half of candidate vertices lost " << pvs.size() << ' ' << clusters.size();
269 
270  if (pvs.empty() && seltks.size() > 5)
271  edm::LogWarning("PrimaryVertexProducer")
272  << "no vertex found with " << seltks.size() << " tracks and " << clusters.size() << " vertex-candidates";
273 
274  // sort vertices by pt**2 vertex (aka signal vertex tagging)
275  if (pvs.size() > 1) {
276  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
277  }
278 
279  // convert transient vertices returned by the theAlgo to (reco) vertices
280  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin(); iv != pvs.end(); iv++) {
281  reco::Vertex v = *iv;
282  vColl.push_back(v);
283  }
284 
285  if (vColl.empty()) {
286  GlobalError bse(beamSpot.rotatedCovariance3D());
287  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
289  we(0, 0) = 10000;
290  we(1, 1) = 10000;
291  we(2, 2) = 10000;
292  vColl.push_back(reco::Vertex(beamSpot.position(), we, 0., 0., 0));
293  if (fVerbose) {
294  std::cout << "RecoVertex/PrimaryVertexProducer: "
295  << "Beamspot with invalid errors " << bse.matrix() << std::endl;
296  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
297  }
298  } else {
299  vColl.push_back(reco::Vertex(beamSpot.position(), beamSpot.rotatedCovariance3D(), 0., 0., 0));
300  if (fVerbose) {
301  std::cout << "RecoVertex/PrimaryVertexProducer: "
302  << " will put Vertex derived from BeamSpot into Event.\n";
303  }
304  }
305  }
306 
307  if (fVerbose) {
308  int ivtx = 0;
309  for (reco::VertexCollection::const_iterator v = vColl.begin(); v != vColl.end(); ++v) {
310  std::cout << "recvtx " << ivtx++ << "#trk " << std::setw(3) << v->tracksSize() << " chi2 " << std::setw(4)
311  << v->chi2() << " ndof " << std::setw(3) << v->ndof() << " x " << std::setw(6) << v->position().x()
312  << " dx " << std::setw(6) << v->xError() << " y " << std::setw(6) << v->position().y() << " dy "
313  << std::setw(6) << v->yError() << " z " << std::setw(6) << v->position().z() << " dz " << std::setw(6)
314  << v->zError();
315  if (f4D) {
316  std::cout << " t " << std::setw(6) << v->t() << " dt " << std::setw(6) << v->tError();
317  }
318  std::cout << std::endl;
319  }
320  }
321 
322  iEvent.put(std::move(result), algorithm->label);
323  }
324 }
325 
327  // offlinePrimaryVertices
329  {
331  vpsd1.add<double>("maxDistanceToBeam", 1.0);
332  vpsd1.add<std::string>("algorithm", "AdaptiveVertexFitter");
333  vpsd1.add<bool>("useBeamConstraint", false);
334  vpsd1.add<std::string>("label", "");
335  vpsd1.add<double>("chi2cutoff", 2.5);
336  vpsd1.add<double>("minNdof", 0.0);
337  std::vector<edm::ParameterSet> temp1;
338  temp1.reserve(2);
339  {
340  edm::ParameterSet temp2;
341  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
342  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
343  temp2.addParameter<bool>("useBeamConstraint", false);
344  temp2.addParameter<std::string>("label", "");
345  temp2.addParameter<double>("chi2cutoff", 2.5);
346  temp2.addParameter<double>("minNdof", 0.0);
347  temp1.push_back(temp2);
348  }
349  {
350  edm::ParameterSet temp2;
351  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
352  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
353  temp2.addParameter<bool>("useBeamConstraint", true);
354  temp2.addParameter<std::string>("label", "WithBS");
355  temp2.addParameter<double>("chi2cutoff", 2.5);
356  temp2.addParameter<double>("minNdof", 2.0);
357  temp1.push_back(temp2);
358  }
359  desc.addVPSet("vertexCollections", vpsd1, temp1);
360  }
361  desc.addUntracked<bool>("verbose", false);
362  {
364  psd0.add<double>("maxNormalizedChi2", 10.0);
365  psd0.add<double>("minPt", 0.0);
366  psd0.add<std::string>("algorithm", "filter");
367  psd0.add<double>("maxEta", 2.4);
368  psd0.add<double>("maxD0Significance", 4.0);
369  psd0.add<double>("maxD0Error", 1.0);
370  psd0.add<double>("maxDzError", 1.0);
371  psd0.add<std::string>("trackQuality", "any");
372  psd0.add<int>("minPixelLayersWithHits", 2);
373  psd0.add<int>("minSiliconLayersWithHits", 5);
374  psd0.add<int>("numTracksThreshold", 0); // HI only
375  desc.add<edm::ParameterSetDescription>("TkFilterParameters", psd0);
376  }
377  desc.add<edm::InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
378  desc.add<edm::InputTag>("TrackLabel", edm::InputTag("generalTracks"));
379  desc.add<edm::InputTag>("TrackTimeResosLabel", edm::InputTag("dummy_default")); // 4D only
380  desc.add<edm::InputTag>("TrackTimesLabel", edm::InputTag("dummy_default")); // 4D only
381  {
383  {
385  psd1.addUntracked<bool>("verbose", false);
386  psd1.addUntracked<double>("zdumpcenter", 0.);
387  psd1.addUntracked<double>("zdumpwidth", 20.);
388  psd1.addUntracked<bool>("use_vdt", false); // obsolete, appears in HLT configs
389  psd1.add<double>("d0CutOff", 3.0);
390  psd1.add<double>("Tmin", 2.0);
391  psd1.add<double>("delta_lowT", 0.001);
392  psd1.add<double>("zmerge", 0.01);
393  psd1.add<double>("dzCutOff", 3.0);
394  psd1.add<double>("Tpurge", 2.0);
395  psd1.add<int>("convergence_mode", 0);
396  psd1.add<double>("delta_highT", 0.01);
397  psd1.add<double>("Tstop", 0.5);
398  psd1.add<double>("coolingFactor", 0.6);
399  psd1.add<double>("vertexSize", 0.006);
400  psd1.add<double>("uniquetrkweight", 0.8);
401  psd1.add<double>("zrange", 4.0);
402 
403  psd1.add<double>("tmerge", 0.01); // 4D only
404  psd1.add<double>("dtCutOff", 4.); // 4D only
405  psd1.add<double>("t0Max", 1.0); // 4D only
406  psd1.add<double>("vertexSizeTime", 0.008); // 4D only
407 
408  psd0.add<edm::ParameterSetDescription>("TkDAClusParameters", psd1);
409 
411  psd2.add<double>("zSeparation", 1.0);
412  psd0.add<edm::ParameterSetDescription>("TkGapClusParameters", psd2);
413  }
414  psd0.add<std::string>("algorithm", "DA_vect");
415  desc.add<edm::ParameterSetDescription>("TkClusParameters", psd0);
416  }
417 
418  descriptions.add("primaryVertexProducer", desc);
419 }
420 
421 //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
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
PrimaryVertexProducer::trkTimeResosToken
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
Definition: PrimaryVertexProducer.h:88
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:87
TransientTrack.h
findQualityFiles.v
v
Definition: findQualityFiles.py:179
edm::Handle< reco::BeamSpot >
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:73
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
reco::BeamSpot
Definition: BeamSpot.h:21
TransientTrackRecord
Definition: TransientTrackRecord.h:11
PrimaryVertexProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PrimaryVertexProducer.cc:326
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
edm::LogWarning
Definition: MessageLogger.h:141
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:674
edm::ParameterSetDescription::addUntracked
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:100
TransientTrackBuilder.h
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
PrimaryVertexProducer::f4D
bool f4D
Definition: PrimaryVertexProducer.h:90
edm::ParameterSet
Definition: ParameterSet.h:36
edm::LogError
Definition: MessageLogger.h:183
edm::ParameterSetDescription::addVPSet
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
Definition: ParameterSetDescription.h:149
runTheMatrix.err
err
Definition: runTheMatrix.py:288
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:124
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
GlobalErrorBase< double, ErrorMatrixTag >
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
VertexDistanceXY
Definition: VertexDistanceXY.h:11
TransientTrackRecord.h
get
#define get
FSQDQM_cfi.pvs
pvs
Definition: FSQDQM_cfi.py:12
HLT_2018_cff.vertexCollections
vertexCollections
Definition: HLT_2018_cff.py:33346
PrimaryVertexProducer::trkToken
edm::EDGetTokenT< reco::TrackCollection > trkToken
Definition: PrimaryVertexProducer.h:86
InputTag.h
VertexFwd.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
TransientVertex.h
DAClusterizerInZ_vect
Definition: DAClusterizerInZ_vect.h:20
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:120
PrimaryVertexProducer::PrimaryVertexProducer
PrimaryVertexProducer(const edm::ParameterSet &)
Definition: PrimaryVertexProducer.cc:20
VertexState
Definition: VertexState.h:13
PrimaryVertexProducer::bsToken
edm::EDGetTokenT< reco::BeamSpot > bsToken
Definition: PrimaryVertexProducer.h:85
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:133
PrimaryVertexProducer::algorithms
std::vector< algo > algorithms
Definition: PrimaryVertexProducer.h:80
mps_fire.result
result
Definition: mps_fire.py:303
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
ntuplemaker.time
time
Definition: ntuplemaker.py:310
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
reco::Vertex
Definition: Vertex.h:35
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22