CMS 3D CMS Logo

IPProducer.h
Go to the documentation of this file.
1 #ifndef RecoBTag_IPProducer
2 #define RecoBTag_IPProducer
3 
4 // system include files
5 #include <algorithm>
6 #include <cmath>
7 #include <functional>
8 #include <iostream>
9 #include <memory>
10 
11 // user include files
22 
28 
32 
41 
47 
48 // system include files
49 #include <memory>
50 
51 // user include files
56 
58 
59 namespace IPProducerHelpers {
60  class FromJTA {
61  public:
64  iC.consumes<reco::JetTracksAssociationCollection>(iConfig.getParameter<edm::InputTag>("jetTracks"))) {}
65  reco::TrackRefVector tracks(const reco::JTATagInfo& it) { return it.tracks(); }
66  std::vector<reco::JTATagInfo> makeBaseVector(const edm::Event& iEvent) {
68  iEvent.getByToken(token_associator, jetTracksAssociation);
69  std::vector<reco::JTATagInfo> bases;
70  size_t i = 0;
71  for (reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin();
72  it != jetTracksAssociation->end();
73  it++, i++) {
74  edm::Ref<reco::JetTracksAssociationCollection> jtaRef(jetTracksAssociation, i);
75  bases.push_back(reco::JTATagInfo(jtaRef));
76  }
77  return bases;
78  }
79 
81  };
83  public:
85  : token_jets(iC.consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>(jets))),
86  token_cands(iC.consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("candidates"))),
87  maxDeltaR(iConfig.getParameter<double>("maxDeltaR")),
88  explicitJTA(iConfig.existsAs<bool>("explicitJTA") ? iConfig.getParameter<bool>("explicitJTA") : false) {}
89 
90  const std::vector<reco::CandidatePtr>& tracks(const reco::JetTagInfo& it) { return m_map[it.jet().key()]; }
91  std::vector<reco::JetTagInfo> makeBaseVector(const edm::Event& iEvent) {
93  iEvent.getByToken(token_jets, jets);
94  std::vector<reco::JetTagInfo> bases;
95 
97  iEvent.getByToken(token_cands, cands);
98  m_map.clear();
99  m_map.resize(jets->size());
100  double maxDeltaR2 = maxDeltaR * maxDeltaR;
101  size_t i = 0;
102  for (edm::View<reco::Jet>::const_iterator it = jets->begin(); it != jets->end(); it++, i++) {
104  bases.push_back(reco::JetTagInfo(jRef));
105  if (explicitJTA) {
106  for (size_t j = 0; j < it->numberOfDaughters(); ++j) {
107  if (it->daughterPtr(j)->bestTrack() != nullptr && it->daughterPtr(j)->charge() != 0) {
108  m_map[i].push_back(it->daughterPtr(j));
109  }
110  }
111  } else {
112  for (size_t j = 0; j < cands->size(); ++j) {
113  if ((*cands)[j].bestTrack() != nullptr && (*cands)[j].charge() != 0 && (*cands)[j].pt() > 0 &&
114  Geom::deltaR2((*cands)[j], (*jets)[i]) < maxDeltaR2) {
115  m_map[i].push_back(cands->ptrAt(j));
116  }
117  }
118  }
119  }
120  return bases;
121  }
122  std::vector<std::vector<reco::CandidatePtr> > m_map;
125  double maxDeltaR;
127  };
128 } // namespace IPProducerHelpers
129 template <class Container, class Base, class Helper>
131 public:
132  typedef std::vector<reco::IPTagInfo<Container, Base> > Product;
133 
134  explicit IPProducer(const edm::ParameterSet&);
135  ~IPProducer() override;
136  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
137 
138  void produce(edm::Event&, const edm::EventSetup&) override;
139 
140 private:
141  void checkEventSetup(const edm::EventSetup& iSetup);
142 
148 
152  std::unique_ptr<HistogramProbabilityEstimator> m_probabilityEstimator;
153  unsigned long long m_calibrationCacheId2D;
154  unsigned long long m_calibrationCacheId3D;
155  bool m_useDB;
156 
159  double m_cutMaxTIP;
160  double m_cutMinPt;
162  double m_cutMaxLIP;
166  Helper m_helper;
167 };
168 
169 // -*- C++ -*-
170 //
171 // Package: IPProducer
172 // Class: IPProducer
173 //
181 //
182 // Original Author: Andrea Rizzi
183 // Created: Thu Apr 6 09:56:23 CEST 2006
184 //
185 //
186 
187 //
188 // constructors and destructor
189 //
190 template <class Container, class Base, class Helper>
192  : m_config(iConfig), m_helper(iConfig, consumesCollector()) {
195 
196  token_primaryVertex = consumes<reco::VertexCollection>(m_config.getParameter<edm::InputTag>("primaryVertex"));
198  esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"));
199  token_calib2D = esConsumes<TrackProbabilityCalibration, BTagTrackProbability2DRcd>();
200  token_calib3D = esConsumes<TrackProbabilityCalibration, BTagTrackProbability3DRcd>();
201 
202  m_computeProbabilities = m_config.getParameter<bool>("computeProbabilities");
203  m_computeGhostTrack = m_config.getParameter<bool>("computeGhostTrack");
204  m_ghostTrackPriorDeltaR = m_config.getParameter<double>("ghostTrackPriorDeltaR");
205  m_cutPixelHits = m_config.getParameter<int>("minimumNumberOfPixelHits");
206  m_cutTotalHits = m_config.getParameter<int>("minimumNumberOfHits");
207  m_cutMaxTIP = m_config.getParameter<double>("maximumTransverseImpactParameter");
208  m_cutMinPt = m_config.getParameter<double>("minimumTransverseMomentum");
209  m_cutMaxChiSquared = m_config.getParameter<double>("maximumChiSquared");
210  m_cutMaxLIP = m_config.getParameter<double>("maximumLongitudinalImpactParameter");
211  m_directionWithTracks = m_config.getParameter<bool>("jetDirectionUsingTracks");
212  m_directionWithGhostTrack = m_config.getParameter<bool>("jetDirectionUsingGhostTrack");
213  m_useTrackQuality = m_config.getParameter<bool>("useTrackQuality");
214 
216  produces<reco::TrackCollection>("ghostTracks");
217  produces<Product>();
218 }
219 
220 template <class Container, class Base, class Helper>
222 
223 //
224 // member functions
225 //
226 // ------------ method called to produce the data ------------
227 template <class Container, class Base, class Helper>
229  // Update probability estimator if event setup is changed
230  if (m_computeProbabilities)
231  checkEventSetup(iSetup);
232 
234  iEvent.getByToken(token_primaryVertex, primaryVertex);
235 
236  edm::ESHandle<TransientTrackBuilder> builder = iSetup.getHandle(token_trackBuilder);
237  // m_algo.setTransientTrackBuilder(builder.product());
238 
239  // output collections
240  auto result = std::make_unique<Product>();
241 
242  std::unique_ptr<reco::TrackCollection> ghostTracks;
243  reco::TrackRefProd ghostTrackRefProd;
244  if (m_computeGhostTrack) {
245  ghostTracks = std::make_unique<reco::TrackCollection>();
246  ghostTrackRefProd = iEvent.getRefBeforePut<reco::TrackCollection>("ghostTracks");
247  }
248 
249  // use first pv of the collection
251  const reco::Vertex* pv = &dummy;
253  if (!primaryVertex->empty()) {
254  pv = &*primaryVertex->begin();
255  // we always use the first vertex (at the moment)
257  } else { // create a dummy PV
259  e(0, 0) = 0.0015 * 0.0015;
260  e(1, 1) = 0.0015 * 0.0015;
261  e(2, 2) = 15. * 15.;
262  reco::Vertex::Point p(0, 0, 0);
263  dummy = reco::Vertex(p, e, 0, 0, 0);
264  }
265 
266  std::vector<Base> baseTagInfos = m_helper.makeBaseVector(iEvent);
267  for (typename std::vector<Base>::const_iterator it = baseTagInfos.begin(); it != baseTagInfos.end(); it++) {
268  Container tracks = m_helper.tracks(*it);
269  math::XYZVector jetMomentum = it->jet()->momentum();
270 
271  if (m_directionWithTracks) {
272  jetMomentum *= 0.5;
273  for (typename Container::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack)
274  if (reco::btag::toTrack(*itTrack)->numberOfValidHits() >= m_cutTotalHits) //minimal quality cuts
275  jetMomentum += (*itTrack)->momentum();
276  }
277 
279  std::vector<reco::TransientTrack> transientTracks;
280 
281  for (typename Container::const_iterator itTrack = tracks.begin(); itTrack != tracks.end(); ++itTrack) {
282  reco::TransientTrack transientTrack = builder->build(*itTrack);
283  const reco::Track& track = transientTrack.track(); //**itTrack;
284  /* cout << " pt " << track.pt() <<
285  " d0 " << fabs(track.d0()) <<
286  " #hit " << track.hitPattern().numberOfValidHits()<<
287  " ipZ " << fabs(track.dz()-pv->z())<<
288  " chi2 " << track.normalizedChi2()<<
289  " #pixel " << track.hitPattern().numberOfValidPixelHits()<< endl;
290 */
291  if (track.pt() > m_cutMinPt &&
292  track.hitPattern().numberOfValidHits() >= m_cutTotalHits && // min num tracker hits
293  track.hitPattern().numberOfValidPixelHits() >= m_cutPixelHits &&
294  track.normalizedChi2() < m_cutMaxChiSquared && std::abs(track.dxy(pv->position())) < m_cutMaxTIP &&
295  std::abs(track.dz(pv->position())) < m_cutMaxLIP) {
296  // std::cout << "selected" << std::endl;
297  selectedTracks.push_back(*itTrack);
298  transientTracks.push_back(transientTrack);
299  }
300  }
301  // std::cout <<"SIZE: " << transientTracks.size() << std::endl;
302  GlobalVector direction(jetMomentum.x(), jetMomentum.y(), jetMomentum.z());
303 
304  std::unique_ptr<reco::GhostTrack> ghostTrack;
305  reco::TrackRef ghostTrackRef;
306  if (m_computeGhostTrack) {
307  reco::GhostTrackFitter fitter;
308  GlobalPoint origin = RecoVertex::convertPos(pv->position());
310  ghostTrack = std::make_unique<reco::GhostTrack>(
311  fitter.fit(origin, error, direction, m_ghostTrackPriorDeltaR, transientTracks));
312 
313  /*
314  if (std::sqrt(jetMomentum.Perp2()) > 30) {
315  double offset = ghostTrack->prediction().lambda(origin);
316  std::cout << "------------------ jet pt " << std::sqrt(jetMomentum.Perp2()) << std::endl;
317  const std::vector<GhostTrackState> *states = &ghostTrack->states();
318  for(std::vector<GhostTrackState>::const_iterator state = states->begin();
319  state != states->end(); ++state) {
320  double dist = state->lambda() - offset;
321  double err = state->lambdaError(ghostTrack->prediction(), error);
322  double ipSig = IPTools::signedImpactParameter3D(state->track(), direction, *pv).second.significance();
323  double axisDist = state->axisDistance(ghostTrack->prediction());
324  std::cout << state->track().impactPointState().freeState()->momentum().perp()
325  << ": " << dist << "/" << err << " [" << (dist / err) << "], ipsig = " << ipSig << ", dist = " << axisDist << ", w = " << state->weight() << std::endl;
326  }
327  }
328 */
329  ghostTrackRef = reco::TrackRef(ghostTrackRefProd, ghostTracks->size());
330  ghostTracks->push_back(*ghostTrack);
331 
332  if (m_directionWithGhostTrack) {
333  const reco::GhostTrackPrediction& pred = ghostTrack->prediction();
334  double lambda = pred.lambda(origin);
337  0,
338  0,
339  0);
340  pv = &dummy;
341  direction = pred.direction();
342  }
343  }
344 
345  std::vector<float> prob2D, prob3D;
346  std::vector<reco::btag::TrackIPData> ipData;
347 
348  for (unsigned int ind = 0; ind < transientTracks.size(); ind++) {
349  const reco::TransientTrack& transientTrack = transientTracks[ind];
350  const reco::Track& track = transientTrack.track();
351 
352  reco::btag::TrackIPData trackIP;
353  trackIP.ip3d = IPTools::signedImpactParameter3D(transientTrack, direction, *pv).second;
354  trackIP.ip2d = IPTools::signedTransverseImpactParameter(transientTrack, direction, *pv).second;
355 
357  IPTools::closestApproachToJet(transientTrack.impactPointState(), *pv, direction, transientTrack.field());
358  if (closest.isValid())
359  trackIP.closestToJetAxis = closest.globalPosition();
360 
361  // TODO: cross check if it is the same using other methods
362  trackIP.distanceToJetAxis = IPTools::jetTrackDistance(transientTrack, direction, *pv).second;
363 
364  if (ghostTrack.get()) {
365  const std::vector<reco::GhostTrackState>& states = ghostTrack->states();
366  std::vector<reco::GhostTrackState>::const_iterator pos =
367  std::find_if(states.begin(), states.end(), [&](auto& arg) { return arg.track() == transientTrack; });
368 
369  if (pos != states.end() && pos->isValid()) {
370  VertexDistance3D dist;
371  const reco::GhostTrackPrediction& pred = ghostTrack->prediction();
372  GlobalPoint p1 = pos->tsos().globalPosition();
373  GlobalError e1 = pos->tsos().cartesianError().position();
374  GlobalPoint p2 = pred.position(pos->lambda());
375  GlobalError e2 = pred.positionError(pos->lambda());
376  trackIP.closestToGhostTrack = p1;
377  trackIP.distanceToGhostTrack = dist.distance(VertexState(p1, e1), VertexState(p2, e2));
378  trackIP.ghostTrackWeight = pos->weight();
379  } else {
380  trackIP.distanceToGhostTrack = Measurement1D(0., -1.);
381  trackIP.ghostTrackWeight = 0.;
382  }
383  } else {
384  trackIP.distanceToGhostTrack = Measurement1D(0., -1.);
385  trackIP.ghostTrackWeight = 1.;
386  }
387 
388  ipData.push_back(trackIP);
389 
390  if (m_computeProbabilities) {
391  //probability with 3D ip
392  std::pair<bool, double> probability = m_probabilityEstimator->probability(
393  m_useTrackQuality, 0, ipData.back().ip3d.significance(), track, *(it->jet()), *pv);
394  prob3D.push_back(probability.first ? probability.second : -1.);
395 
396  //probability with 2D ip
397  probability = m_probabilityEstimator->probability(
398  m_useTrackQuality, 1, ipData.back().ip2d.significance(), track, *(it->jet()), *pv);
399  prob2D.push_back(probability.first ? probability.second : -1.);
400  }
401  }
402 
403  result->push_back(
404  typename Product::value_type(ipData, prob2D, prob3D, selectedTracks, *it, pvRef, direction, ghostTrackRef));
405  }
406 
407  if (m_computeGhostTrack)
408  iEvent.put(std::move(ghostTracks), "ghostTracks");
409  iEvent.put(std::move(result));
410 }
411 
418 
419 template <class Container, class Base, class Helper>
423  unsigned long long cacheId2D = re2D.cacheIdentifier();
424  unsigned long long cacheId3D = re3D.cacheIdentifier();
425 
426  if (cacheId2D != m_calibrationCacheId2D || cacheId3D != m_calibrationCacheId3D) //Calibration changed
427  {
428  //iSetup.get<BTagTrackProbabilityRcd>().get(calib);
429  edm::ESHandle<TrackProbabilityCalibration> calib2DHandle = iSetup.getHandle(token_calib2D);
430  edm::ESHandle<TrackProbabilityCalibration> calib3DHandle = iSetup.getHandle(token_calib3D);
431 
432  const TrackProbabilityCalibration* ca2D = calib2DHandle.product();
433  const TrackProbabilityCalibration* ca3D = calib3DHandle.product();
434 
435  m_probabilityEstimator = std::make_unique<HistogramProbabilityEstimator>(ca3D, ca2D);
436  }
437  m_calibrationCacheId3D = cacheId3D;
438  m_calibrationCacheId2D = cacheId2D;
439 }
440 
441 // Specialized templates used to fill 'descriptions'
442 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
443 template <>
445  edm::ConfigurationDescriptions& descriptions) {
447  desc.add<double>("maximumTransverseImpactParameter", 0.2);
448  desc.add<int>("minimumNumberOfHits", 8);
449  desc.add<double>("minimumTransverseMomentum", 1.0);
450  desc.add<edm::InputTag>("primaryVertex", edm::InputTag("offlinePrimaryVertices"));
451  desc.add<double>("maximumLongitudinalImpactParameter", 17.0);
452  desc.add<bool>("computeGhostTrack", true);
453  desc.add<double>("ghostTrackPriorDeltaR", 0.03);
454  desc.add<edm::InputTag>("jetTracks", edm::InputTag("ak4JetTracksAssociatorAtVertexPF"));
455  desc.add<bool>("jetDirectionUsingGhostTrack", false);
456  desc.add<int>("minimumNumberOfPixelHits", 2);
457  desc.add<bool>("jetDirectionUsingTracks", false);
458  desc.add<bool>("computeProbabilities", true);
459  desc.add<bool>("useTrackQuality", false);
460  desc.add<double>("maximumChiSquared", 5.0);
461  descriptions.addDefault(desc);
462 }
463 
464 template <>
465 inline void
467  edm::ConfigurationDescriptions& descriptions) {
469  desc.add<double>("maximumTransverseImpactParameter", 0.2);
470  desc.add<int>("minimumNumberOfHits", 8);
471  desc.add<double>("minimumTransverseMomentum", 1.0);
472  desc.add<edm::InputTag>("primaryVertex", edm::InputTag("offlinePrimaryVertices"));
473  desc.add<double>("maximumLongitudinalImpactParameter", 17.0);
474  desc.add<bool>("computeGhostTrack", true);
475  desc.add<double>("maxDeltaR", 0.4);
476  desc.add<edm::InputTag>("candidates", edm::InputTag("particleFlow"));
477  desc.add<bool>("jetDirectionUsingGhostTrack", false);
478  desc.add<int>("minimumNumberOfPixelHits", 2);
479  desc.add<bool>("jetDirectionUsingTracks", false);
480  desc.add<bool>("computeProbabilities", true);
481  desc.add<bool>("useTrackQuality", false);
482  desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
483  desc.add<double>("ghostTrackPriorDeltaR", 0.03);
484  desc.add<double>("maximumChiSquared", 5.0);
485  desc.addOptional<bool>("explicitJTA", false);
486  descriptions.addDefault(desc);
487 }
488 
489 #endif
~IPProducer() override
Definition: IPProducer.h:221
const Vector & prediction() const
reco::Vertex::Point convertPos(const GlobalPoint &p)
double m_cutMaxChiSquared
Definition: IPProducer.h:161
const Track & track() const
std::vector< reco::JTATagInfo > makeBaseVector(const edm::Event &iEvent)
Definition: IPProducer.h:66
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool m_directionWithGhostTrack
Definition: IPProducer.h:164
transient_vector_type::const_iterator const_iterator
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
unsigned long long m_calibrationCacheId3D
Definition: IPProducer.h:154
GlobalPoint position(double lambda=0.) const
edm::ESGetToken< TrackProbabilityCalibration, BTagTrackProbability2DRcd > token_calib2D
Definition: IPProducer.h:146
JetTracksAssociation::Container JetTracksAssociationCollection
typedefs for backward compatibility
Measurement1D ip2d
Definition: IPTagInfo.h:31
int closest(std::vector< int > const &vec, int value)
FromJetAndCands(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC, const std::string &jets="jets")
Definition: IPProducer.h:84
double m_cutMaxTIP
Definition: IPProducer.h:159
unsigned long long m_calibrationCacheId2D
Definition: IPProducer.h:153
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:57
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
TrackRefVector tracks(void) const override
returns a list of tracks associated to the jet
Definition: JTATagInfo.h:20
double m_cutMinPt
Definition: IPProducer.h:160
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:81
bool m_useTrackQuality
Definition: IPProducer.h:165
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:182
const_iterator begin() const
edm::ESGetToken< TrackProbabilityCalibration, BTagTrackProbability3DRcd > token_calib3D
Definition: IPProducer.h:147
const reco::Track * toTrack(const reco::TrackBaseRef &t)
Definition: IPTagInfo.h:24
A arg
Definition: Factorize.h:31
int m_cutTotalHits
Definition: IPProducer.h:158
unsigned long long cacheIdentifier() const
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:206
double m_ghostTrackPriorDeltaR
Definition: IPProducer.h:151
std::vector< std::vector< reco::CandidatePtr > > m_map
Definition: IPProducer.h:122
reco::TransientTrack build(const reco::Track *p) const
void checkEventSetup(const edm::EventSetup &iSetup)
Definition: IPProducer.h:420
std::unique_ptr< HistogramProbabilityEstimator > m_probabilityEstimator
Definition: IPProducer.h:152
const GlobalVector direction() const
const std::vector< reco::CandidatePtr > & tracks(const reco::JetTagInfo &it)
Definition: IPProducer.h:90
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
std::vector< reco::IPTagInfo< Container, Base > > Product
Definition: IPProducer.h:132
void addDefault(ParameterSetDescription const &psetDescription)
edm::RefToBase< Jet > jet(void) const override
returns a polymorphic reference to the tagged jet
Definition: JetTagInfo.h:22
GlobalError positionError(double lambda=0.) const
Measurement1D distanceToGhostTrack
Definition: IPTagInfo.h:34
Definition: Jet.py:1
FromJTA(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
Definition: IPProducer.h:62
int m_cutPixelHits
Definition: IPProducer.h:157
const MagneticField * field() const
Helper m_helper
Definition: IPProducer.h:166
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T get() const
Definition: EventSetup.h:79
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
Definition: IPProducer.h:144
reco::TrackRefVector tracks(const reco::JTATagInfo &it)
Definition: IPProducer.h:65
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
Measurement1D ip3d
Definition: IPTagInfo.h:32
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
edm::EDGetTokenT< reco::JetTracksAssociationCollection > token_associator
Definition: IPProducer.h:80
bool m_directionWithTracks
Definition: IPProducer.h:163
const_iterator end() const
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
const edm::ParameterSet & m_config
Definition: IPProducer.h:143
size_t key() const
Definition: RefToBase.h:219
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
GlobalPoint closestToJetAxis
Definition: IPTagInfo.h:29
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
edm::EDGetTokenT< edm::View< reco::Jet > > token_jets
Definition: IPProducer.h:123
double lambda(const GlobalPoint &point) const
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > token_trackBuilder
Definition: IPProducer.h:145
Measurement1D distanceToJetAxis
Definition: IPTagInfo.h:33
void produce(edm::Event &, const edm::EventSetup &) override
Definition: IPProducer.h:228
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
IPProducer(const edm::ParameterSet &)
Definition: IPProducer.h:191
GlobalPoint closestToGhostTrack
Definition: IPTagInfo.h:30
edm::EDGetTokenT< edm::View< reco::Candidate > > token_cands
Definition: IPProducer.h:124
GhostTrack fit(const GlobalPoint &priorPosition, const GlobalError &priorError, const GlobalVector &direction, double coneRadius, const std::vector< TransientTrack > &tracks) const
std::vector< reco::JetTagInfo > makeBaseVector(const edm::Event &iEvent)
Definition: IPProducer.h:91
bool m_useDB
Definition: IPProducer.h:155
bool m_computeProbabilities
Definition: IPProducer.h:149
primaryVertex
hltOfflineBeamSpot for HLTMON
bool m_computeGhostTrack
Definition: IPProducer.h:150
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
double m_cutMaxLIP
Definition: IPProducer.h:162
def move(src, dest)
Definition: eostools.py:511
TrajectoryStateOnSurface impactPointState() const