CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PrimaryVertexProducer.cc
Go to the documentation of this file.
8 
12 
15 
17 
19  : theTTBToken(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))), theConfig(conf) {
20  fVerbose = conf.getUntrackedParameter<bool>("verbose", false);
21 
22  trkToken = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackLabel"));
23  bsToken = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpotLabel"));
24  f4D = false;
25 
26  // select and configure the track selection
27  std::string trackSelectionAlgorithm =
28  conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
29  if (trackSelectionAlgorithm == "filter") {
30  theTrackFilter = new TrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
31  } else if (trackSelectionAlgorithm == "filterWithThreshold") {
33  } else {
34  throw VertexException("PrimaryVertexProducer: unknown track selection algorithm: " + trackSelectionAlgorithm);
35  }
36 
37  // select and configure the track clusterizer
38  std::string clusteringAlgorithm =
39  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
40  if (clusteringAlgorithm == "gap") {
42  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
43  } else if (clusteringAlgorithm == "DA") {
45  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
46  }
47  // provide the vectorized version of the clusterizer, if supported by the build
48  else if (clusteringAlgorithm == "DA_vect") {
50  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
51  } else if (clusteringAlgorithm == "DA2D_vect") {
53  conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
54  f4D = true;
55  }
56 
57  else {
58  throw VertexException("PrimaryVertexProducer: unknown clustering algorithm: " + clusteringAlgorithm);
59  }
60 
61  if (f4D) {
62  trkTimesToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimesLabel"));
63  trkTimeResosToken = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("TrackTimeResosLabel"));
64  }
65 
66  // select and configure the vertex fitters
67  if (conf.exists("vertexCollections")) {
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 {
81  throw VertexException("PrimaryVertexProducer: unknown algorithm: " + fitterAlgorithm);
82  }
83  algorithm.label = algoconf->getParameter<std::string>("label");
84  algorithm.minNdof = algoconf->getParameter<double>("minNdof");
85  algorithm.useBeamConstraint = algoconf->getParameter<bool>("useBeamConstraint");
86  algorithm.vertexSelector =
87  new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
88  algorithms.push_back(algorithm);
89 
90  produces<reco::VertexCollection>(algorithm.label);
91  }
92  } else {
93  edm::LogWarning("MisConfiguration")
94  << "this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
95 
97  std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
98  if (fitterAlgorithm == "KalmanVertexFitter") {
99  algorithm.fitter = new KalmanVertexFitter();
100  } else if (fitterAlgorithm == "AdaptiveVertexFitter") {
101  algorithm.fitter = new AdaptiveVertexFitter();
102  } else {
103  throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);
104  }
105  algorithm.label = "";
106  algorithm.minNdof = conf.getParameter<double>("minNdof");
107  algorithm.useBeamConstraint = conf.getParameter<bool>("useBeamConstraint");
108 
111  conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
112 
113  algorithms.push_back(algorithm);
114  produces<reco::VertexCollection>(algorithm.label);
115  }
116 
117  //check if this is a recovery iteration
118  fRecoveryIteration = conf.getParameter<bool>("isRecoveryIteration");
119  if (fRecoveryIteration) {
120  if (algorithms.empty()) {
121  throw VertexException("PrimaryVertexProducer: No algorithm specified. ");
122  } else if (algorithms.size() > 1) {
123  throw VertexException(
124  "PrimaryVertexProducer: Running in Recovery mode and more than one algorithm specified. Please "
125  "only one algorithm.");
126  }
127  recoveryVtxToken = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("recoveryVtxCollection"));
128  }
129 }
130 
132  if (theTrackFilter)
133  delete theTrackFilter;
135  delete theTrackClusterizer;
136  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
137  if (algorithm->fitter)
138  delete algorithm->fitter;
139  if (algorithm->vertexSelector)
140  delete algorithm->vertexSelector;
141  }
142 }
143 
145  // get the BeamSpot, it will always be needed, even when not used as a constraint
147  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
148  iEvent.getByToken(bsToken, recoBeamSpotHandle);
149  if (recoBeamSpotHandle.isValid()) {
150  beamSpot = *recoBeamSpotHandle;
151  } else {
152  edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
153  }
154 
155  bool validBS = true;
156  VertexState beamVertexState(beamSpot);
157  if ((beamVertexState.error().cxx() <= 0.) || (beamVertexState.error().cyy() <= 0.) ||
158  (beamVertexState.error().czz() <= 0.)) {
159  validBS = false;
160  edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors " << beamVertexState.error().matrix();
161  }
162 
163  //if this is a recovery iteration, check if we already have a valid PV
164  if (fRecoveryIteration) {
165  auto const& oldVertices = iEvent.get(recoveryVtxToken);
166  //look for the first valid (not-BeamSpot) vertex
167  for (auto const& old : oldVertices) {
168  if (!(old.isFake())) {
169  //found a valid vertex, write the first one to the collection and return
170  //otherwise continue with regular vertexing procedure
171  auto result = std::make_unique<reco::VertexCollection>();
172  result->push_back(old);
173  iEvent.put(std::move(result), algorithms.begin()->label);
174  return;
175  }
176  }
177  }
178 
179  // get RECO tracks from the event
180  // `tks` can be used as a ptr to a reco::TrackCollection
182  iEvent.getByToken(trkToken, tks);
183 
184  // interface RECO tracks to vertex reconstruction
185  const auto& theB = &iSetup.getData(theTTBToken);
186  std::vector<reco::TransientTrack> t_tks;
187 
188  if (f4D) {
189  edm::Handle<edm::ValueMap<float> > trackTimesH;
190  edm::Handle<edm::ValueMap<float> > trackTimeResosH;
191  iEvent.getByToken(trkTimesToken, trackTimesH);
192  iEvent.getByToken(trkTimeResosToken, trackTimeResosH);
193  t_tks = (*theB).build(tks, beamSpot, *(trackTimesH.product()), *(trackTimeResosH.product()));
194  } else {
195  t_tks = (*theB).build(tks, beamSpot);
196  }
197  if (fVerbose) {
198  std::cout << "RecoVertex/PrimaryVertexProducer"
199  << "Found: " << t_tks.size() << " reconstructed tracks"
200  << "\n";
201  }
202 
203  // select tracks
204  std::vector<reco::TransientTrack>&& seltks = theTrackFilter->select(t_tks);
205 
206  // clusterize tracks in Z
207  std::vector<std::vector<reco::TransientTrack> >&& clusters = theTrackClusterizer->clusterize(seltks);
208 
209  if (fVerbose) {
210  std::cout << " clustering returned " << clusters.size() << " clusters from " << seltks.size()
211  << " selected tracks" << std::endl;
212  }
213 
214  // vertex fits
215  for (std::vector<algo>::const_iterator algorithm = algorithms.begin(); algorithm != algorithms.end(); algorithm++) {
216  auto result = std::make_unique<reco::VertexCollection>();
217  reco::VertexCollection& vColl = (*result);
218 
219  std::vector<TransientVertex> pvs;
220  for (std::vector<std::vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin();
221  iclus != clusters.end();
222  iclus++) {
223  double sumwt = 0.;
224  double sumwt2 = 0.;
225  double sumw = 0.;
226  double meantime = 0.;
227  double vartime = 0.;
228  if (f4D) {
229  for (const auto& tk : *iclus) {
230  const double time = tk.timeExt();
231  const double err = tk.dtErrorExt();
232  const double inverr = err > 0. ? 1.0 / err : 0.;
233  const double w = inverr * inverr;
234  sumwt += w * time;
235  sumwt2 += w * time * time;
236  sumw += w;
237  }
238  meantime = sumwt / sumw;
239  double sumsq = sumwt2 - sumwt * sumwt / sumw;
240  double chisq = iclus->size() > 1 ? sumsq / double(iclus->size() - 1) : sumsq / double(iclus->size());
241  vartime = chisq / sumw;
242  }
243 
245  if (algorithm->useBeamConstraint && validBS && (iclus->size() > 1)) {
246  v = algorithm->fitter->vertex(*iclus, beamSpot);
247  } else if (!(algorithm->useBeamConstraint) && (iclus->size() > 1)) {
248  v = algorithm->fitter->vertex(*iclus);
249  } // else: no fit ==> v.isValid()=False
250 
251  // 4D vertices: add timing information
252  if (f4D and v.isValid()) {
253  auto err = v.positionError().matrix4D();
254  err(3, 3) = vartime;
255  auto trkWeightMap3d = v.weightMap(); // copy the 3d-fit weights
256  v = TransientVertex(v.position(), meantime, err, v.originalTracks(), v.totalChiSquared(), v.degreesOfFreedom());
257  v.weightMap(trkWeightMap3d);
258  }
259 
260  if (fVerbose) {
261  if (v.isValid()) {
262  std::cout << "x,y,z";
263  if (f4D)
264  std::cout << ",t";
265  std::cout << "=" << v.position().x() << " " << v.position().y() << " " << v.position().z();
266  if (f4D)
267  std::cout << " " << v.time();
268  std::cout << " cluster size = " << (*iclus).size() << std::endl;
269  } else {
270  std::cout << "Invalid fitted vertex, cluster size=" << (*iclus).size() << std::endl;
271  }
272  }
273 
274  if (v.isValid() && (v.degreesOfFreedom() >= algorithm->minNdof) &&
275  (!validBS || (*(algorithm->vertexSelector))(v, beamVertexState)))
276  pvs.push_back(v);
277  } // end of cluster loop
278 
279  if (fVerbose) {
280  std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << pvs.size() << std::endl;
281  }
282 
283  if (clusters.size() > 2 && clusters.size() > 2 * pvs.size())
284  edm::LogWarning("PrimaryVertexProducer")
285  << "more than half of candidate vertices lost " << pvs.size() << ' ' << clusters.size();
286 
287  if (pvs.empty() && seltks.size() > 5)
288  edm::LogWarning("PrimaryVertexProducer")
289  << "no vertex found with " << seltks.size() << " tracks and " << clusters.size() << " vertex-candidates";
290 
291  // sort vertices by pt**2 vertex (aka signal vertex tagging)
292  if (pvs.size() > 1) {
293  sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
294  }
295 
296  // convert transient vertices returned by the theAlgo to (reco) vertices
297  for (std::vector<TransientVertex>::const_iterator iv = pvs.begin(); iv != pvs.end(); iv++) {
298  reco::Vertex v = *iv;
299  vColl.push_back(v);
300  }
301 
302  if (vColl.empty()) {
303  GlobalError bse(beamSpot.rotatedCovariance3D());
304  if ((bse.cxx() <= 0.) || (bse.cyy() <= 0.) || (bse.czz() <= 0.)) {
306  we(0, 0) = 10000;
307  we(1, 1) = 10000;
308  we(2, 2) = 10000;
309  vColl.push_back(reco::Vertex(beamSpot.position(), we, 0., 0., 0));
310  if (fVerbose) {
311  std::cout << "RecoVertex/PrimaryVertexProducer: "
312  << "Beamspot with invalid errors " << bse.matrix() << std::endl;
313  std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
314  }
315  } else {
316  vColl.push_back(reco::Vertex(beamSpot.position(), beamSpot.rotatedCovariance3D(), 0., 0., 0));
317  if (fVerbose) {
318  std::cout << "RecoVertex/PrimaryVertexProducer: "
319  << " will put Vertex derived from BeamSpot into Event.\n";
320  }
321  }
322  }
323 
324  if (fVerbose) {
325  int ivtx = 0;
326  for (reco::VertexCollection::const_iterator v = vColl.begin(); v != vColl.end(); ++v) {
327  std::cout << "recvtx " << ivtx++ << "#trk " << std::setw(3) << v->tracksSize() << " chi2 " << std::setw(4)
328  << v->chi2() << " ndof " << std::setw(3) << v->ndof() << " x " << std::setw(6) << v->position().x()
329  << " dx " << std::setw(6) << v->xError() << " y " << std::setw(6) << v->position().y() << " dy "
330  << std::setw(6) << v->yError() << " z " << std::setw(6) << v->position().z() << " dz " << std::setw(6)
331  << v->zError();
332  if (f4D) {
333  std::cout << " t " << std::setw(6) << v->t() << " dt " << std::setw(6) << v->tError();
334  }
335  std::cout << std::endl;
336  }
337  }
338 
339  iEvent.put(std::move(result), algorithm->label);
340  }
341 }
342 
344  // offlinePrimaryVertices
346  {
348  vpsd1.add<double>("maxDistanceToBeam", 1.0);
349  vpsd1.add<std::string>("algorithm", "AdaptiveVertexFitter");
350  vpsd1.add<bool>("useBeamConstraint", false);
351  vpsd1.add<std::string>("label", "");
352  vpsd1.add<double>("chi2cutoff", 2.5);
353  vpsd1.add<double>("minNdof", 0.0);
354  std::vector<edm::ParameterSet> temp1;
355  temp1.reserve(2);
356  {
357  edm::ParameterSet temp2;
358  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
359  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
360  temp2.addParameter<bool>("useBeamConstraint", false);
361  temp2.addParameter<std::string>("label", "");
362  temp2.addParameter<double>("chi2cutoff", 2.5);
363  temp2.addParameter<double>("minNdof", 0.0);
364  temp1.push_back(temp2);
365  }
366  {
367  edm::ParameterSet temp2;
368  temp2.addParameter<double>("maxDistanceToBeam", 1.0);
369  temp2.addParameter<std::string>("algorithm", "AdaptiveVertexFitter");
370  temp2.addParameter<bool>("useBeamConstraint", true);
371  temp2.addParameter<std::string>("label", "WithBS");
372  temp2.addParameter<double>("chi2cutoff", 2.5);
373  temp2.addParameter<double>("minNdof", 2.0);
374  temp1.push_back(temp2);
375  }
376  desc.addVPSet("vertexCollections", vpsd1, temp1);
377  }
378  desc.addUntracked<bool>("verbose", false);
379  {
382  psd0.add<int>("numTracksThreshold", 0); // HI only
383  desc.add<edm::ParameterSetDescription>("TkFilterParameters", psd0);
384  }
385  desc.add<edm::InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
386  desc.add<edm::InputTag>("TrackLabel", edm::InputTag("generalTracks"));
387  desc.add<edm::InputTag>("TrackTimeResosLabel", edm::InputTag("dummy_default")); // 4D only
388  desc.add<edm::InputTag>("TrackTimesLabel", edm::InputTag("dummy_default")); // 4D only
389 
390  {
392  {
395  psd0.add<edm::ParameterSetDescription>("TkDAClusParameters", psd1);
396 
399  psd0.add<edm::ParameterSetDescription>("TkGapClusParameters", psd2);
400  }
401  psd0.add<std::string>("algorithm", "DA_vect");
402  desc.add<edm::ParameterSetDescription>("TkClusParameters", psd0);
403  }
404 
405  desc.add<bool>("isRecoveryIteration", false);
406  desc.add<edm::InputTag>("recoveryVtxCollection", {""});
407 
408  descriptions.add("primaryVertexProducer", desc);
409 }
410 
411 //define this as a plug-in
GlobalError positionError() const
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
tuple vertexCollections
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
edm::EDGetTokenT< reco::BeamSpot > bsToken
int32_t *__restrict__ iv
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > theTTBToken
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const double w
Definition: UKUtility.cc:23
std::vector< algo > algorithms
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
Common base class.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::VertexCollection > recoveryVtxToken
const AlgebraicSymMatrix33 matrix() const
float totalChiSquared() const
edm::EDGetTokenT< reco::TrackCollection > trkToken
TrackClusterizerInZ * theTrackClusterizer
T y() const
Definition: PV3DBase.h:60
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Log< level::Error, false > LogError
edm::EDGetTokenT< edm::ValueMap< float > > trkTimesToken
const AlgebraicSymMatrix44 & matrix4D() const
tuple result
Definition: mps_fire.py:311
bool getData(T &iHolder) const
Definition: EventSetup.h:128
int iEvent
Definition: GenABIO.cc:224
std::vector< reco::TransientTrack > const & originalTracks() const
edm::EDGetTokenT< edm::ValueMap< float > > trkTimeResosToken
float degreesOfFreedom() const
static void fillPSetDescription(edm::ParameterSetDescription &desc)
VertexCompatibleWithBeam * vertexSelector
virtual std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const =0
GlobalPoint position() const
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
T z() const
Definition: PV3DBase.h:61
def move
Definition: eostools.py:511
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
TrackFilterForPVFindingBase * theTrackFilter
virtual std::vector< reco::TransientTrack > select(const std::vector< reco::TransientTrack > &tracks) const =0
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
T const * product() const
Definition: Handle.h:70
static void fillPSetDescription(edm::ParameterSetDescription &desc)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
double time() const
void produce(edm::Event &, const edm::EventSetup &) override
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
PrimaryVertexProducer(const edm::ParameterSet &)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
GlobalError error() const
Definition: VertexState.h:64
tuple cout
Definition: gather_cfg.py:144
const Point & position() const
position
Definition: BeamSpot.h:59
Covariance3DMatrix rotatedCovariance3D() const
Definition: BeamSpot.cc:73
Log< level::Warning, false > LogWarning
T x() const
Definition: PV3DBase.h:59
bool isValid() const
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
TransientTrackToFloatMap weightMap() const