CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
IPProducer.h
Go to the documentation of this file.
1 #ifndef RecoBTag_IPProducer
2 #define RecoBTag_IPProducer
3 
4 // system include files
5 #include <cmath>
6 #include <memory>
7 #include <iostream>
8 #include <algorithm>
9 
10 #include "boost/bind.hpp"
11 
12 // user include files
20 
26 
30 
39 
43 
44 // system include files
45 #include <memory>
46 
47 // user include files
52 
54 using namespace std;
55 using namespace reco;
56 using namespace edm;
57 using boost::bind;
58 
59 namespace IPProducerHelpers {
60  class FromJTA{
61  public:
62  FromJTA(const edm::ParameterSet& iConfig, edm::ConsumesCollector && iC) : token_associator(iC.consumes<reco::JetTracksAssociationCollection>(iConfig.getParameter<edm::InputTag>("jetTracks")))
63  {}
65  {
66  return it.tracks();
67  }
68  std::vector<reco::JTATagInfo> makeBaseVector(edm::Event& iEvent){
70  iEvent.getByToken(token_associator, jetTracksAssociation);
71  std::vector<reco::JTATagInfo> bases;
72  size_t i = 0;
74  jetTracksAssociation->begin();
75  it != jetTracksAssociation->end(); it++, i++) {
76  Ref<JetTracksAssociationCollection> jtaRef(jetTracksAssociation, i);
77  bases.push_back(reco::JTATagInfo(jtaRef));
78  }
79  return bases;
80  }
81 
83  };
85  public:
86  FromJetAndCands(const edm::ParameterSet& iConfig, edm::ConsumesCollector && iC): token_jets(iC.consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"))),
87  token_cands(iC.consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("candidates"))), maxDeltaR(iConfig.getParameter<double>("maxDeltaR")),
88  explicitJTA(iConfig.existsAs<bool>("explicitJTA") ? iConfig.getParameter<bool>("explicitJTA") : false) {}
89 
90  std::vector<reco::CandidatePtr> tracks(edm::Event&,const reco::JetTagInfo & it)
91  {
92  return m_map[it.jet().key()];
93  }
94  std::vector<reco::JetTagInfo> makeBaseVector(edm::Event& iEvent){
96  iEvent.getByToken(token_jets, jets);
97  std::vector<reco::JetTagInfo> bases;
98 
100  iEvent.getByToken(token_cands, cands);
101  m_map.clear();
102  m_map.resize(jets->size());
103  double maxDeltaR2 = maxDeltaR*maxDeltaR;
104  size_t i = 0;
105  for(edm::View<reco::Jet>::const_iterator it = jets->begin();
106  it != jets->end(); it++, i++) {
107  edm::RefToBase<reco::Jet> jRef(jets, i);
108  bases.push_back(reco::JetTagInfo(jRef));
109  if( explicitJTA )
110  {
111  for(size_t j=0;j<it->numberOfDaughters();++j) {
112  if( it->daughterPtr(j)->bestTrack()!=0 && it->daughterPtr(j)->charge() !=0 ){
113  m_map[i].push_back(it->daughterPtr(j));
114  }
115  }
116  }
117  else
118  {
119  for(size_t j=0;j<cands->size();++j) {
120  if( (*cands)[j].bestTrack()!=0 && Geom::deltaR2((*cands)[j].p4(),(*jets)[i].p4()) < maxDeltaR2 && (*cands)[j].charge() !=0 ){
121  m_map[i].push_back(cands->ptrAt(j));
122  }
123  }
124  }
125  }
126  return bases;
127  }
128  std::vector<std::vector<reco::CandidatePtr> > m_map;
131  double maxDeltaR;
133  };
134 }
135 template <class Container, class Base, class Helper>
137  public:
138  typedef std::vector<reco::IPTagInfo<Container,Base> > Product;
139 
140 
141  explicit IPProducer(const edm::ParameterSet&);
142  ~IPProducer();
143 
144  virtual void produce(edm::Event&, const edm::EventSetup&);
145  private:
146  void checkEventSetup(const edm::EventSetup & iSetup);
147 
150 
154  std::auto_ptr<HistogramProbabilityEstimator> m_probabilityEstimator;
155  unsigned long long m_calibrationCacheId2D;
156  unsigned long long m_calibrationCacheId3D;
157  bool m_useDB;
158 
161  double m_cutMaxTIP;
162  double m_cutMinPt;
164  double m_cutMaxLIP;
168  Helper m_helper;
169 };
170 
171 // -*- C++ -*-
172 //
173 // Package: IPProducer
174 // Class: IPProducer
175 //
183 //
184 // Original Author: Andrea Rizzi
185 // Created: Thu Apr 6 09:56:23 CEST 2006
186 //
187 //
188 
189 
190 
191 
192 
193 //
194 // constructors and destructor
195 //
196 template <class Container, class Base, class Helper> IPProducer<Container,Base,Helper>::IPProducer(const edm::ParameterSet& iConfig) :
197  m_config(iConfig),m_helper(iConfig,consumesCollector())
198 {
201 
202  token_primaryVertex = consumes<reco::VertexCollection>(m_config.getParameter<edm::InputTag>("primaryVertex"));
203 
204  m_computeProbabilities = m_config.getParameter<bool>("computeProbabilities");
205  m_computeGhostTrack = m_config.getParameter<bool>("computeGhostTrack");
206  m_ghostTrackPriorDeltaR = m_config.getParameter<double>("ghostTrackPriorDeltaR");
207  m_cutPixelHits = m_config.getParameter<int>("minimumNumberOfPixelHits");
208  m_cutTotalHits = m_config.getParameter<int>("minimumNumberOfHits");
209  m_cutMaxTIP = m_config.getParameter<double>("maximumTransverseImpactParameter");
210  m_cutMinPt = m_config.getParameter<double>("minimumTransverseMomentum");
211  m_cutMaxChiSquared = m_config.getParameter<double>("maximumChiSquared");
212  m_cutMaxLIP = m_config.getParameter<double>("maximumLongitudinalImpactParameter");
213  m_directionWithTracks = m_config.getParameter<bool>("jetDirectionUsingTracks");
214  m_directionWithGhostTrack = m_config.getParameter<bool>("jetDirectionUsingGhostTrack");
215  m_useTrackQuality = m_config.getParameter<bool>("useTrackQuality");
216 
217  if (m_computeGhostTrack)
218  produces<reco::TrackCollection>("ghostTracks");
219  produces<Product>();
220 }
221 
222 template <class Container, class Base, class Helper> IPProducer<Container,Base,Helper>::~IPProducer()
223 {
224 }
225 
226 //
227 // member functions
228 //
229 // ------------ method called to produce the data ------------
230 template <class Container, class Base, class Helper> void
232 {
233  // Update probability estimator if event setup is changed
234  if (m_computeProbabilities)
235  checkEventSetup(iSetup);
236 
237 
239  iEvent.getByToken(token_primaryVertex, primaryVertex);
240 
242  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
243  // m_algo.setTransientTrackBuilder(builder.product());
244 
245  // output collections
246  auto_ptr<Product> result(new Product);
247 
248  auto_ptr<reco::TrackCollection> ghostTracks;
249  TrackRefProd ghostTrackRefProd;
250  if (m_computeGhostTrack) {
251  ghostTracks.reset(new reco::TrackCollection);
252  ghostTrackRefProd = iEvent.getRefBeforePut<TrackCollection>("ghostTracks");
253  }
254 
255  // use first pv of the collection
256  Vertex dummy;
257  const Vertex *pv = &dummy;
259  if (primaryVertex->size() != 0) {
260  pv = &*primaryVertex->begin();
261  // we always use the first vertex (at the moment)
263  } else { // create a dummy PV
265  e(0, 0) = 0.0015 * 0.0015;
266  e(1, 1) = 0.0015 * 0.0015;
267  e(2, 2) = 15. * 15.;
268  Vertex::Point p(0, 0, 0);
269  dummy = Vertex(p, e, 0, 0, 0);
270  }
271 
272  std::vector<Base> baseTagInfos = m_helper.makeBaseVector(iEvent);
273  for(typename std::vector<Base>::const_iterator it = baseTagInfos.begin(); it != baseTagInfos.end(); it++) {
274  Container tracks = m_helper.tracks(iEvent,*it);
275  math::XYZVector jetMomentum = it->jet()->momentum();
276 
277  if (m_directionWithTracks) {
278  jetMomentum *= 0.5;
279  for(typename Container::const_iterator itTrack = tracks.begin();
280  itTrack != tracks.end(); ++itTrack)
281  if (reco::btag::toTrack(*itTrack)->numberOfValidHits() >= m_cutTotalHits) //minimal quality cuts
282  jetMomentum += (*itTrack)->momentum();
283  }
284 
286  vector<TransientTrack> transientTracks;
287 
288  for(typename Container::const_iterator itTrack = tracks.begin();
289  itTrack != tracks.end(); ++itTrack) {
290  TransientTrack transientTrack = builder->build(*itTrack);
291  const Track & track = transientTrack.track(); //**itTrack;
292  /* cout << " pt " << track.pt() <<
293  " d0 " << fabs(track.d0()) <<
294  " #hit " << track.hitPattern().numberOfValidHits()<<
295  " ipZ " << fabs(track.dz()-pv->z())<<
296  " chi2 " << track.normalizedChi2()<<
297  " #pixel " << track.hitPattern().numberOfValidPixelHits()<< endl;
298 */
299  if (track.pt() > m_cutMinPt &&
300  track.hitPattern().numberOfValidHits() >= m_cutTotalHits && // min num tracker hits
301  track.hitPattern().numberOfValidPixelHits() >= m_cutPixelHits &&
302  track.normalizedChi2() < m_cutMaxChiSquared &&
303  std::abs(track.dxy(pv->position())) < m_cutMaxTIP &&
304  std::abs(track.dz(pv->position())) < m_cutMaxLIP) {
305 // std::cout << "selected" << std::endl;
306  selectedTracks.push_back(*itTrack);
307  transientTracks.push_back(transientTrack);
308  }
309  }
310 // std::cout <<"SIZE: " << transientTracks.size() << std::endl;
311  GlobalVector direction(jetMomentum.x(), jetMomentum.y(), jetMomentum.z());
312 
313  auto_ptr<GhostTrack> ghostTrack;
314  TrackRef ghostTrackRef;
315  if (m_computeGhostTrack) {
316  GhostTrackFitter fitter;
319  ghostTrack.reset(new GhostTrack(fitter.fit(origin, error, direction,
320  m_ghostTrackPriorDeltaR,
321  transientTracks)));
322 
323 /*
324  if (std::sqrt(jetMomentum.Perp2()) > 30) {
325  double offset = ghostTrack->prediction().lambda(origin);
326  std::cout << "------------------ jet pt " << std::sqrt(jetMomentum.Perp2()) << std::endl;
327  const std::vector<GhostTrackState> *states = &ghostTrack->states();
328  for(std::vector<GhostTrackState>::const_iterator state = states->begin();
329  state != states->end(); ++state) {
330  double dist = state->lambda() - offset;
331  double err = state->lambdaError(ghostTrack->prediction(), error);
332  double ipSig = IPTools::signedImpactParameter3D(state->track(), direction, *pv).second.significance();
333  double axisDist = state->axisDistance(ghostTrack->prediction());
334  std::cout << state->track().impactPointState().freeState()->momentum().perp()
335  << ": " << dist << "/" << err << " [" << (dist / err) << "], ipsig = " << ipSig << ", dist = " << axisDist << ", w = " << state->weight() << std::endl;
336  }
337  }
338 */
339  ghostTrackRef = TrackRef(ghostTrackRefProd, ghostTracks->size());
340  ghostTracks->push_back(*ghostTrack);
341 
342  if (m_directionWithGhostTrack) {
343  const GhostTrackPrediction &pred = ghostTrack->prediction();
344  double lambda = pred.lambda(origin);
345  dummy = Vertex(RecoVertex::convertPos(pred.position(lambda)),
347  0, 0, 0);
348  pv = &dummy;
349  direction = pred.direction();
350  }
351  }
352 
353  vector<float> prob2D, prob3D;
354  vector<reco::btag::TrackIPData> ipData;
355 
356  for(unsigned int ind = 0; ind < transientTracks.size(); ind++) {
357  const TransientTrack &transientTrack = transientTracks[ind];
358  const Track & track = transientTrack.track();
359 
360  reco::btag::TrackIPData trackIP;
361  trackIP.ip3d = IPTools::signedImpactParameter3D(transientTrack, direction, *pv).second;
362  trackIP.ip2d = IPTools::signedTransverseImpactParameter(transientTrack, direction, *pv).second;
363 
364  TrajectoryStateOnSurface closest =
366  *pv, direction,
367  transientTrack.field());
368  if (closest.isValid())
369  trackIP.closestToJetAxis = closest.globalPosition();
370 
371  // TODO: cross check if it is the same using other methods
372  trackIP.distanceToJetAxis = IPTools::jetTrackDistance(transientTrack, direction, *pv).second;
373 
374  if (ghostTrack.get()) {
375  const std::vector<GhostTrackState> &states = ghostTrack->states();
376  std::vector<GhostTrackState>::const_iterator pos =
377  std::find_if(states.begin(), states.end(),
378  bind(equal_to<TransientTrack>(),
379  bind(&GhostTrackState::track, _1),
380  transientTrack));
381 
382  if (pos != states.end() && pos->isValid()) {
383  VertexDistance3D dist;
384  const GhostTrackPrediction &pred = ghostTrack->prediction();
385  GlobalPoint p1 = pos->tsos().globalPosition();
386  GlobalError e1 = pos->tsos().cartesianError().position();
387  GlobalPoint p2 = pred.position(pos->lambda());
388  GlobalError e2 = pred.positionError(pos->lambda());
389  trackIP.closestToGhostTrack = p1;
390  trackIP.distanceToGhostTrack = dist.distance(VertexState(p1, e1),
391  VertexState(p2, e2));
392  trackIP.ghostTrackWeight = pos->weight();
393  } else {
394  trackIP.distanceToGhostTrack = Measurement1D(-1. -1.);
395  trackIP.ghostTrackWeight = 0.;
396  }
397  } else {
398  trackIP.distanceToGhostTrack = Measurement1D(-1. -1.);
399  trackIP.ghostTrackWeight = 1.;
400  }
401 
402  ipData.push_back(trackIP);
403 
404  if (m_computeProbabilities) {
405  //probability with 3D ip
406  pair<bool,double> probability = m_probabilityEstimator->probability(m_useTrackQuality, 0,ipData.back().ip3d.significance(),track,*(it->jet()),*pv);
407  prob3D.push_back(probability.first ? probability.second : -1.);
408 
409  //probability with 2D ip
410  probability = m_probabilityEstimator->probability(m_useTrackQuality,1,ipData.back().ip2d.significance(),track,*(it->jet()),*pv);
411  prob2D.push_back(probability.first ? probability.second : -1.);
412  }
413  }
414 
415  result->push_back(typename Product::value_type(ipData, prob2D, prob3D, selectedTracks,
416  *it, pvRef, direction, ghostTrackRef));
417  }
418 
419  if (m_computeGhostTrack)
420  iEvent.put(ghostTracks, "ghostTracks");
421  iEvent.put(result);
422 }
423 
424 
431 
432 template <class Container, class Base, class Helper> void IPProducer<Container,Base,Helper>::checkEventSetup(const EventSetup & iSetup)
433  {
434  using namespace edm;
435  using namespace edm::eventsetup;
436 
437  const EventSetupRecord & re2D= iSetup.get<BTagTrackProbability2DRcd>();
438  const EventSetupRecord & re3D= iSetup.get<BTagTrackProbability3DRcd>();
439  unsigned long long cacheId2D= re2D.cacheIdentifier();
440  unsigned long long cacheId3D= re3D.cacheIdentifier();
441 
442  if(cacheId2D!=m_calibrationCacheId2D || cacheId3D!=m_calibrationCacheId3D ) //Calibration changed
443  {
444  //iSetup.get<BTagTrackProbabilityRcd>().get(calib);
446  iSetup.get<BTagTrackProbability2DRcd>().get(calib2DHandle);
448  iSetup.get<BTagTrackProbability3DRcd>().get(calib3DHandle);
449 
450  const TrackProbabilityCalibration * ca2D= calib2DHandle.product();
451  const TrackProbabilityCalibration * ca3D= calib3DHandle.product();
452 
453  m_probabilityEstimator.reset(new HistogramProbabilityEstimator(ca3D,ca2D));
454 
455  }
456  m_calibrationCacheId3D=cacheId3D;
457  m_calibrationCacheId2D=cacheId2D;
458 }
459 
460 
461 #endif
virtual Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const
reco::Vertex::Point convertPos(const GlobalPoint &p)
T getParameter(std::string const &) const
unsigned long long cacheIdentifier() const
double m_cutMaxChiSquared
Definition: IPProducer.h:163
int i
Definition: DBlmapReader.cc:9
bool m_directionWithGhostTrack
Definition: IPProducer.h:166
virtual edm::RefToBase< Jet > jet(void) const
returns a polymorphic reference to the tagged jet
Definition: JetTagInfo.h:22
EcalChannelStatus Container
double lambda(const GlobalPoint &point) const
transient_vector_type::const_iterator const_iterator
unsigned long long m_calibrationCacheId3D
Definition: IPProducer.h:156
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:609
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
int numberOfValidHits() const
Definition: HitPattern.h:734
Base class for all types of Jets.
Definition: Jet.h:20
Measurement1D ip2d
Definition: IPTagInfo.h:31
double m_cutMaxTIP
Definition: IPProducer.h:161
unsigned long long m_calibrationCacheId2D
Definition: IPProducer.h:155
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:50
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
double m_cutMinPt
Definition: IPProducer.h:162
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:71
GlobalPoint globalPosition() const
bool m_useTrackQuality
Definition: IPProducer.h:167
TrajectoryStateOnSurface closestApproachToJet(const TrajectoryStateOnSurface &state, const reco::Vertex &vertex, const GlobalVector &aJetDirection, const MagneticField *field)
Definition: IPTools.cc:177
const MagneticField * field() const
std::auto_ptr< HistogramProbabilityEstimator > m_probabilityEstimator
Definition: IPProducer.h:154
const reco::Track * toTrack(const reco::TrackBaseRef &t)
Definition: IPTagInfo.h:24
std::vector< reco::CandidatePtr > tracks(edm::Event &, const reco::JetTagInfo &it)
Definition: IPProducer.h:90
const Point & position() const
position
Definition: Vertex.h:106
int m_cutTotalHits
Definition: IPProducer.h:160
std::vector< reco::IPTagInfo< Container, Base > > Product
Definition: IPProducer.h:138
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:200
double m_ghostTrackPriorDeltaR
Definition: IPProducer.h:153
GhostTrack fit(const GlobalPoint &priorPosition, const GlobalError &priorError, const GlobalVector &direction, double coneRadius, const std::vector< TransientTrack > &tracks) const
std::vector< reco::JetTagInfo > makeBaseVector(edm::Event &iEvent)
Definition: IPProducer.h:94
std::vector< std::vector< reco::CandidatePtr > > m_map
Definition: IPProducer.h:128
virtual TrackRefVector tracks(void) const
returns a list of tracks associated to the jet
Definition: JTATagInfo.h:21
void checkEventSetup(const edm::EventSetup &iSetup)
Definition: IPProducer.h:432
int iEvent
Definition: GenABIO.cc:230
Measurement1D distanceToGhostTrack
Definition: IPTagInfo.h:34
FromJTA(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
Definition: IPProducer.h:62
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
size_t key() const
Definition: RefToBase.h:229
double p4[4]
Definition: TauolaWrapper.h:92
int m_cutPixelHits
Definition: IPProducer.h:159
Helper m_helper
Definition: IPProducer.h:168
reco::TrackRefVector tracks(edm::Event &, const reco::JTATagInfo &it)
Definition: IPProducer.h:64
vector< PseudoJet > jets
double pt() const
track transverse momentum
Definition: TrackBase.h:669
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
Definition: IPProducer.h:149
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
edm::EDGetTokenT< JetTracksAssociationCollection > token_associator
Definition: IPProducer.h:82
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:868
Measurement1D ip3d
Definition: IPTagInfo.h:32
Container::value_type value_type
double p2[4]
Definition: TauolaWrapper.h:90
bool m_directionWithTracks
Definition: IPProducer.h:165
RefProd< PROD > getRefBeforePut()
Definition: Event.h:133
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:657
const edm::ParameterSet & m_config
Definition: IPProducer.h:148
GlobalPoint closestToJetAxis
Definition: IPTagInfo.h:29
const Track & track() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:384
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const T & get() const
Definition: EventSetup.h:55
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:19
edm::EDGetTokenT< edm::View< reco::Jet > > token_jets
Definition: IPProducer.h:129
T const * product() const
Definition: ESHandle.h:86
Error error() const
return SMatrix
Definition: Vertex.h:129
Measurement1D distanceToJetAxis
Definition: IPTagInfo.h:33
const Vector & prediction() const
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: IPProducer.h:231
double p1[4]
Definition: TauolaWrapper.h:89
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
FromJetAndCands(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
Definition: IPProducer.h:86
IPProducer(const edm::ParameterSet &)
Definition: IPProducer.h:196
GlobalPoint position(double lambda=0.) const
GlobalPoint closestToGhostTrack
Definition: IPTagInfo.h:30
edm::EDGetTokenT< edm::View< reco::Candidate > > token_cands
Definition: IPProducer.h:130
int numberOfValidPixelHits() const
Definition: HitPattern.h:749
bool m_useDB
Definition: IPProducer.h:157
TrajectoryStateOnSurface impactPointState() const
bool m_computeProbabilities
Definition: IPProducer.h:151
volatile std::atomic< bool > shutdown_flag false
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:639
bool m_computeGhostTrack
Definition: IPProducer.h:152
double m_cutMaxLIP
Definition: IPProducer.h:164
std::vector< reco::JTATagInfo > makeBaseVector(edm::Event &iEvent)
Definition: IPProducer.h:68
GlobalError positionError(double lambda=0.) const
const GlobalVector direction() const