CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SecondaryVertexProducer.cc
Go to the documentation of this file.
1 #include <functional>
2 #include <algorithm>
3 #include <iterator>
4 #include <cstddef>
5 #include <string>
6 #include <vector>
7 #include <map>
8 #include <set>
9 
10 #include <boost/iterator/transform_iterator.hpp>
11 
20 
24 
28 
36 
40 
46 
47 using namespace reco;
48 
49 namespace {
50  template<typename T>
51  struct RefToBaseLess : public std::binary_function<edm::RefToBase<T>,
52  edm::RefToBase<T>,
53  bool> {
54  inline bool operator()(const edm::RefToBase<T> &r1,
55  const edm::RefToBase<T> &r2) const
56  {
57  return r1.id() < r2.id() ||
58  (r1.id() == r2.id() && r1.key() < r2.key());
59  }
60  };
61 }
62 
64  public:
65  explicit SecondaryVertexProducer(const edm::ParameterSet &params);
67 
68  virtual void produce(edm::Event &event, const edm::EventSetup &es);
69 
70  private:
72  CONSTRAINT_NONE = 0,
77  CONSTRAINT_PV_PRIMARIES_IN_FIT
78  };
79 
80  static ConstraintType getConstraintType(const std::string &name);
81 
94 };
95 
98 {
99  if (name == "None")
100  return CONSTRAINT_NONE;
101  else if (name == "BeamSpot")
102  return CONSTRAINT_BEAMSPOT;
103  else if (name == "BeamSpot+PVPosition")
104  return CONSTRAINT_PV_BEAMSPOT_SIZE;
105  else if (name == "BeamSpotZ+PVErrorScaledXY")
106  return CONSTRAINT_PV_BS_Z_ERRORS_SCALED;
107  else if (name == "PVErrorScaled")
108  return CONSTRAINT_PV_ERROR_SCALED;
109  else if (name == "BeamSpot+PVTracksInFit")
110  return CONSTRAINT_PV_PRIMARIES_IN_FIT;
111  else
112  throw cms::Exception("InvalidArgument")
113  << "SecondaryVertexProducer: ``constraint'' parameter "
114  "value \"" << name << "\" not understood."
115  << std::endl;
116 }
117 
119 getGhostTrackFitType(const std::string &name)
120 {
121  if (name == "AlwaysWithGhostTrack")
123  else if (name == "SingleTracksWithGhostTrack")
125  else if (name == "RefitGhostTrackWithVertices")
127  else
128  throw cms::Exception("InvalidArgument")
129  << "SecondaryVertexProducer: ``fitType'' "
130  "parameter value \"" << name << "\" for "
131  "GhostTrackVertexFinder settings not "
132  "understood." << std::endl;
133 }
134 
136  const edm::ParameterSet &params) :
137  trackIPTagInfoLabel(params.getParameter<edm::InputTag>("trackIPTagInfos")),
138  sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
139  trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
140  constraint(getConstraintType(params.getParameter<std::string>("constraint"))),
141  constraintScaling(1.0),
142  vtxRecoPSet(params.getParameter<edm::ParameterSet>("vertexReco")),
143  useGhostTrack(vtxRecoPSet.getParameter<std::string>("finder") == "gtvr"),
144  withPVError(params.getParameter<bool>("usePVError")),
145  minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
146  vertexFilter(params.getParameter<edm::ParameterSet>("vertexCuts")),
147  vertexSorting(params.getParameter<edm::ParameterSet>("vertexSelection"))
148 {
151  constraintScaling = params.getParameter<double>("pvErrorScaling");
152 
157  beamSpotTag = params.getParameter<edm::InputTag>("beamSpotTag");
158 
159  produces<SecondaryVertexTagInfoCollection>();
160 }
161 
163 {
164 }
165 
166 namespace {
167  struct SVBuilder :
168  public std::unary_function<const Vertex&, SecondaryVertex> {
169 
170  SVBuilder(const reco::Vertex &pv,
171  const GlobalVector &direction,
172  bool withPVError) :
173  pv(pv), direction(direction),
174  withPVError(withPVError) {}
175 
176  SecondaryVertex operator () (const reco::Vertex &sv) const
177  { return SecondaryVertex(pv, sv, direction, withPVError); }
178 
179  const Vertex &pv;
180  const GlobalVector &direction;
181  bool withPVError;
182  };
183 
184  struct SVFilter :
185  public std::unary_function<const SecondaryVertex&, bool> {
186 
187  SVFilter(const VertexFilter &filter, const Vertex &pv,
188  const GlobalVector &direction) :
189  filter(filter), pv(pv), direction(direction) {}
190 
191  inline bool operator () (const SecondaryVertex &sv) const
192  { return !filter(pv, sv, direction); }
193 
194  const VertexFilter &filter;
195  const Vertex &pv;
196  const GlobalVector &direction;
197  };
198 
199 } // anonynmous namespace
200 
202  const edm::EventSetup &es)
203 {
205  RefToBaseLess<Track> > TransientTrackMap;
206 
208  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
209  trackBuilder);
210 
212  event.getByLabel(trackIPTagInfoLabel, trackIPTagInfos);
213 
214  edm::Handle<BeamSpot> beamSpot;
215  unsigned int bsCovSrc[7] = { 0, };
216  double sigmaZ = 0.0, beamWidth = 0.0;
217  switch(constraint) {
219  event.getByLabel(beamSpotTag,beamSpot);
220  bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = bsCovSrc[6] = 1;
221  sigmaZ = beamSpot->sigmaZ();
222  beamWidth = beamSpot->BeamWidthX();
223  break;
224 
226  event.getByLabel(beamSpotTag,beamSpot);
227  bsCovSrc[0] = bsCovSrc[1] = 2;
228  bsCovSrc[3] = bsCovSrc[4] = bsCovSrc[5] = 1;
229  sigmaZ = beamSpot->sigmaZ();
230  break;
231 
233  bsCovSrc[0] = bsCovSrc[1] = bsCovSrc[2] = 2;
234  break;
235 
236  case CONSTRAINT_BEAMSPOT:
238  event.getByLabel(beamSpotTag,beamSpot);
239  break;
240 
241  default:
242  /* nothing */;
243  }
244 
245  std::auto_ptr<ConfigurableVertexReconstructor> vertexReco;
246  std::auto_ptr<GhostTrackVertexFinder> vertexRecoGT;
247  if (useGhostTrack)
248  vertexRecoGT.reset(new GhostTrackVertexFinder(
249  vtxRecoPSet.getParameter<double>("maxFitChi2"),
250  vtxRecoPSet.getParameter<double>("mergeThreshold"),
251  vtxRecoPSet.getParameter<double>("primcut"),
252  vtxRecoPSet.getParameter<double>("seccut"),
253  getGhostTrackFitType(vtxRecoPSet.getParameter<std::string>("fitType"))));
254  else
255  vertexReco.reset(
257 
258  TransientTrackMap primariesMap;
259 
260  // result secondary vertices
261 
262  std::auto_ptr<SecondaryVertexTagInfoCollection>
263  tagInfos(new SecondaryVertexTagInfoCollection);
264 
265  for(TrackIPTagInfoCollection::const_iterator iterJets =
266  trackIPTagInfos->begin(); iterJets != trackIPTagInfos->end();
267  ++iterJets) {
268 
269  std::vector<SecondaryVertexTagInfo::IndexedTrackData> trackData;
270 
271  const Vertex &pv = *iterJets->primaryVertex();
272 
273  std::set<TransientTrack> primaries;
275  for(Vertex::trackRef_iterator iter = pv.tracks_begin();
276  iter != pv.tracks_end(); ++iter) {
277  TransientTrackMap::iterator pos =
278  primariesMap.lower_bound(*iter);
279 
280  if (pos != primariesMap.end() &&
281  pos->first == *iter)
282  primaries.insert(pos->second);
283  else {
284  TransientTrack track =
285  trackBuilder->build(
286  iter->castTo<TrackRef>());
287  primariesMap.insert(pos,
288  std::make_pair(*iter, track));
289  primaries.insert(track);
290  }
291  }
292  }
293 
294  edm::RefToBase<Jet> jetRef = iterJets->jet();
295 
296  GlobalVector jetDir(jetRef->momentum().x(),
297  jetRef->momentum().y(),
298  jetRef->momentum().z());
299 
300  std::vector<std::size_t> indices =
301  iterJets->sortedIndexes(sortCriterium);
302 
303  TrackRefVector trackRefs = iterJets->sortedTracks(indices);
304 
305  const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
306  iterJets->impactParameterData();
307 
308  // build transient tracks used for vertex reconstruction
309 
310  std::vector<TransientTrack> fitTracks;
311  std::vector<GhostTrackState> gtStates;
312  std::auto_ptr<GhostTrackPrediction> gtPred;
313  if (useGhostTrack)
314  gtPred.reset(new GhostTrackPrediction(
315  *iterJets->ghostTrack()));
316 
317  for(unsigned int i = 0; i < indices.size(); i++) {
318  typedef SecondaryVertexTagInfo::IndexedTrackData IndexedTrackData;
319 
320  const TrackRef &trackRef = trackRefs[i];
321 
322  trackData.push_back(IndexedTrackData());
323  trackData.back().first = indices[i];
324 
325  // select tracks for SV finder
326 
327  if (!trackSelector(*trackRef, ipData[indices[i]], *jetRef,
329  pv.position()))) {
330  trackData.back().second.svStatus =
331  SecondaryVertexTagInfo::TrackData::trackSelected;
332  continue;
333  }
334 
335  TransientTrackMap::const_iterator pos =
336  primariesMap.find(
337  TrackBaseRef(trackRef));
338  TransientTrack fitTrack;
339  if (pos != primariesMap.end()) {
340  primaries.erase(pos->second);
341  fitTrack = pos->second;
342  } else
343  fitTrack = trackBuilder->build(trackRef);
344  fitTracks.push_back(fitTrack);
345 
346  trackData.back().second.svStatus =
347  SecondaryVertexTagInfo::TrackData::trackUsedForVertexFit;
348 
349  if (useGhostTrack) {
350  GhostTrackState gtState(fitTrack);
351  GlobalPoint pos =
352  ipData[indices[i]].closestToGhostTrack;
353  gtState.linearize(*gtPred, true,
354  gtPred->lambda(pos));
355  gtState.setWeight(ipData[indices[i]].ghostTrackWeight);
356  gtStates.push_back(gtState);
357  }
358  }
359 
360  std::auto_ptr<GhostTrack> ghostTrack;
361  if (useGhostTrack)
362  ghostTrack.reset(new GhostTrack(
366  GlobalVector(
367  iterJets->ghostTrack()->px(),
368  iterJets->ghostTrack()->py(),
369  iterJets->ghostTrack()->pz()),
370  0.05),
371  *gtPred, gtStates,
372  iterJets->ghostTrack()->chi2(),
373  iterJets->ghostTrack()->ndof()));
374 
375  // perform actual vertex finding
376 
377  std::vector<TransientVertex> fittedSVs;
378  switch(constraint) {
379  case CONSTRAINT_NONE:
380  if (useGhostTrack)
381  fittedSVs = vertexRecoGT->vertices(
382  pv, *ghostTrack);
383  else
384  fittedSVs = vertexReco->vertices(fitTracks);
385  break;
386 
387  case CONSTRAINT_BEAMSPOT:
388  if (useGhostTrack)
389  fittedSVs = vertexRecoGT->vertices(
390  pv, *beamSpot, *ghostTrack);
391  else
392  fittedSVs = vertexReco->vertices(fitTracks,
393  *beamSpot);
394  break;
395 
400  for(unsigned int i = 0; i < 7; i++) {
401  unsigned int covSrc = bsCovSrc[i];
402  for(unsigned int j = 0; j < 7; j++) {
403  double v=0.0;
404  if (!covSrc || bsCovSrc[j] != covSrc)
405  v = 0.0;
406  else if (covSrc == 1)
407  v = beamSpot->covariance(i, j);
408  else if (j<3 && i<3)
409  v = pv.covariance(i, j) *
411  cov(i, j) = v;
412  }
413  }
414 
415  BeamSpot bs(pv.position(), sigmaZ,
416  beamSpot.isValid() ? beamSpot->dxdz() : 0.,
417  beamSpot.isValid() ? beamSpot->dydz() : 0.,
418  beamWidth, cov, BeamSpot::Unknown);
419 
420  if (useGhostTrack)
421  fittedSVs = vertexRecoGT->vertices(
422  pv, bs, *ghostTrack);
423  else
424  fittedSVs = vertexReco->vertices(fitTracks, bs);
425  } break;
426 
428  std::vector<TransientTrack> primaries_(
429  primaries.begin(), primaries.end());
430  if (useGhostTrack)
431  fittedSVs = vertexRecoGT->vertices(
432  pv, *beamSpot, primaries_,
433  *ghostTrack);
434  else
435  fittedSVs = vertexReco->vertices(
436  primaries_, fitTracks,
437  *beamSpot);
438  } break;
439  }
440 
441  // build combined SV information and filter
442 
443  std::vector<SecondaryVertex> SVs;
444  SVBuilder svBuilder(pv, jetDir, withPVError);
445  std::remove_copy_if(boost::make_transform_iterator(
446  fittedSVs.begin(), svBuilder),
447  boost::make_transform_iterator(
448  fittedSVs.end(), svBuilder),
449  std::back_inserter(SVs),
450  SVFilter(vertexFilter, pv, jetDir));
451 
452  // clean up now unneeded collections
453 
454  gtPred.reset();
455  ghostTrack.reset();
456  gtStates.clear();
457  fitTracks.clear();
458  fittedSVs.clear();
459 
460  // sort SVs by importance
461 
462  std::vector<unsigned int> vtxIndices = vertexSorting(SVs);
463 
464  std::vector<SecondaryVertexTagInfo::VertexData> svData;
465 
466  svData.resize(vtxIndices.size());
467  for(unsigned int idx = 0; idx < vtxIndices.size(); idx++) {
468  const SecondaryVertex &sv = SVs[vtxIndices[idx]];
469 
470  svData[idx].vertex = sv;
471  svData[idx].dist2d = sv.dist2d();
472  svData[idx].dist3d = sv.dist3d();
473  svData[idx].direction =
474  GlobalVector(sv.x() - pv.x(),
475  sv.y() - pv.y(),
476  sv.z() - pv.z());
477 
478  // mark tracks successfully used in vertex fit
479 
480  for(Vertex::trackRef_iterator iter = sv.tracks_begin();
481  iter != sv.tracks_end(); iter++) {
482  if (sv.trackWeight(*iter) < minTrackWeight)
483  continue;
484 
486  std::find(trackRefs.begin(), trackRefs.end(),
487  iter->castTo<TrackRef>());
488  if (pos == trackRefs.end())
489  throw cms::Exception("TrackNotFound")
490  << "Could not find track from secondary "
491  "vertex in original tracks."
492  << std::endl;
493 
494  unsigned int index = pos - trackRefs.begin();
495  trackData[index].second.svStatus =
497  ((unsigned int)SecondaryVertexTagInfo::TrackData::trackAssociatedToVertex + idx);
498  }
499  }
500 
501  // fill result into tag infos
502 
503  tagInfos->push_back(
505  trackData, svData, SVs.size(),
506  TrackIPTagInfoRef(trackIPTagInfos,
507  iterJets - trackIPTagInfos->begin())));
508  }
509 
510  event.put(tagInfos);
511 }
512 
513 //define this as a plug-in
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:32
reco::Vertex::Point convertPos(const GlobalPoint &p)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
reco::TrackIPTagInfo::SortCriteria getCriterium(const std::string &name)
Definition: TrackSorting.cc:11
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
static ConstraintType getConstraintType(const std::string &name)
Measurement1D dist3d() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double y() const
y coordinate
Definition: Vertex.h:97
SecondaryVertexProducer(const edm::ParameterSet &params)
const edm::InputTag trackIPTagInfoLabel
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
virtual Vector momentum() const
spatial momentum vector
virtual void produce(edm::Event &event, const edm::EventSetup &es)
void setWeight(double weight)
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:110
ProductID id() const
Definition: RefToBase.h:220
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:243
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
const Point & position() const
position
Definition: Vertex.h:93
TrackIPTagInfo::SortCriteria sortCriterium
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:238
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
int j
Definition: DBlmapReader.cc:9
double z() const
y coordinate
Definition: Vertex.h:99
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:76
size_t key() const
Definition: RefToBase.h:228
Measurement1D dist2d() const
double x() const
x coordinate
Definition: Vertex.h:95
std::pair< unsigned int, TrackData > IndexedTrackData
tuple filter
USE THIS FOR SKIMMED TRACKS process.p = cms.Path(process.hltLevel1GTSeed*process.skimming*process.offlineBeamSpot*process.TrackRefitter2) OTHERWISE USE THIS.
Definition: align_tpl.py:86
bool linearize(const GhostTrackPrediction &pred, bool initial=false, double lambda=0.)
const T & get() const
Definition: EventSetup.h:55
Error error() const
return SMatrix
Definition: Vertex.h:116
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:60
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
mathSSE::Vec4< T > v
static GhostTrackVertexFinder::FitType getGhostTrackFitType(const std::string &name)
Global3DVector GlobalVector
Definition: GlobalVector.h:10
tuple trackSelector
Tracks selection.
Definition: valSkim_cff.py:4