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 
89  std::vector<reco::CandidatePtr> tracks(edm::Event&,const reco::JetTagInfo & it)
90  {
91  return m_map[it.jet().key()];
92  }
93  std::vector<reco::JetTagInfo> makeBaseVector(edm::Event& iEvent){
95  iEvent.getByToken(token_jets, jets);
96  std::vector<reco::JetTagInfo> bases;
97 
99  iEvent.getByToken(token_cands, cands);
100  m_map.clear();
101  m_map.resize(jets->size());
102  size_t i = 0;
103  for(edm::View<reco::Jet>::const_iterator it = jets->begin();
104  it != jets->end(); it++, i++) {
105  edm::RefToBase<reco::Jet> jRef(jets, i);
106  bases.push_back(jRef);
107  //FIXME: add deltaR or any other requirement here
108  for(size_t j=0;j<cands->size();j++) {
109  if((*cands)[j].bestTrack()!=0 && ROOT::Math::VectorUtil::DeltaR((*cands)[j].p4(),(*jets)[i].p4()) < maxDeltaR && (*cands)[j].charge() !=0 ){
110  m_map[i].push_back(cands->ptrAt(j));
111  }
112  }
113  }
114  return bases;
115  }
116  std::vector<std::vector<reco::CandidatePtr> > m_map;
119  double maxDeltaR;
120  };
121 }
122 template <class Container, class Base, class Helper>
124  public:
125  typedef std::vector<reco::IPTagInfo<Container,Base> > Product;
126 
127 
128  explicit IPProducer(const edm::ParameterSet&);
129  ~IPProducer();
130 
131  virtual void produce(edm::Event&, const edm::EventSetup&);
132  private:
133  void checkEventSetup(const edm::EventSetup & iSetup);
134 
137 
141  std::auto_ptr<HistogramProbabilityEstimator> m_probabilityEstimator;
142  unsigned long long m_calibrationCacheId2D;
143  unsigned long long m_calibrationCacheId3D;
144  bool m_useDB;
145 
148  double m_cutMaxTIP;
149  double m_cutMinPt;
151  double m_cutMaxLIP;
155  Helper m_helper;
156 };
157 
158 // -*- C++ -*-
159 //
160 // Package: IPProducer
161 // Class: IPProducer
162 //
170 //
171 // Original Author: Andrea Rizzi
172 // Created: Thu Apr 6 09:56:23 CEST 2006
173 //
174 //
175 
176 
177 
178 
179 
180 //
181 // constructors and destructor
182 //
183 template <class Container, class Base, class Helper> IPProducer<Container,Base,Helper>::IPProducer(const edm::ParameterSet& iConfig) :
184  m_config(iConfig),m_helper(iConfig,consumesCollector())
185 {
188 
189  token_primaryVertex = consumes<reco::VertexCollection>(m_config.getParameter<edm::InputTag>("primaryVertex"));
190 
191  m_computeProbabilities = m_config.getParameter<bool>("computeProbabilities");
192  m_computeGhostTrack = m_config.getParameter<bool>("computeGhostTrack");
193  m_ghostTrackPriorDeltaR = m_config.getParameter<double>("ghostTrackPriorDeltaR");
194  m_cutPixelHits = m_config.getParameter<int>("minimumNumberOfPixelHits");
195  m_cutTotalHits = m_config.getParameter<int>("minimumNumberOfHits");
196  m_cutMaxTIP = m_config.getParameter<double>("maximumTransverseImpactParameter");
197  m_cutMinPt = m_config.getParameter<double>("minimumTransverseMomentum");
198  m_cutMaxChiSquared = m_config.getParameter<double>("maximumChiSquared");
199  m_cutMaxLIP = m_config.getParameter<double>("maximumLongitudinalImpactParameter");
200  m_directionWithTracks = m_config.getParameter<bool>("jetDirectionUsingTracks");
201  m_directionWithGhostTrack = m_config.getParameter<bool>("jetDirectionUsingGhostTrack");
202  m_useTrackQuality = m_config.getParameter<bool>("useTrackQuality");
203 
204  if (m_computeGhostTrack)
205  produces<reco::TrackCollection>("ghostTracks");
206  produces<Product>();
207 }
208 
209 template <class Container, class Base, class Helper> IPProducer<Container,Base,Helper>::~IPProducer()
210 {
211 }
212 
213 //
214 // member functions
215 //
216 // ------------ method called to produce the data ------------
217 template <class Container, class Base, class Helper> void
219 {
220  // Update probability estimator if event setup is changed
221  if (m_computeProbabilities)
222  checkEventSetup(iSetup);
223 
224 
226  iEvent.getByToken(token_primaryVertex, primaryVertex);
227 
229  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
230  // m_algo.setTransientTrackBuilder(builder.product());
231 
232  // output collections
233  auto_ptr<Product> result(new Product);
234 
235  auto_ptr<reco::TrackCollection> ghostTracks;
236  TrackRefProd ghostTrackRefProd;
237  if (m_computeGhostTrack) {
238  ghostTracks.reset(new reco::TrackCollection);
239  ghostTrackRefProd = iEvent.getRefBeforePut<TrackCollection>("ghostTracks");
240  }
241 
242  // use first pv of the collection
243  Vertex dummy;
244  const Vertex *pv = &dummy;
246  if (primaryVertex->size() != 0) {
247  pv = &*primaryVertex->begin();
248  // we always use the first vertex (at the moment)
250  } else { // create a dummy PV
252  e(0, 0) = 0.0015 * 0.0015;
253  e(1, 1) = 0.0015 * 0.0015;
254  e(2, 2) = 15. * 15.;
255  Vertex::Point p(0, 0, 0);
256  dummy = Vertex(p, e, 0, 0, 0);
257  }
258 
259  std::vector<Base> baseTagInfos = m_helper.makeBaseVector(iEvent);
260  for(typename std::vector<Base>::const_iterator it = baseTagInfos.begin(); it != baseTagInfos.end(); it++) {
261  Container tracks = m_helper.tracks(iEvent,*it);
262  math::XYZVector jetMomentum = it->jet()->momentum();
263 
264  if (m_directionWithTracks) {
265  jetMomentum *= 0.5;
266  for(typename Container::const_iterator itTrack = tracks.begin();
267  itTrack != tracks.end(); ++itTrack)
268  if (reco::btag::toTrack(*itTrack)->numberOfValidHits() >= m_cutTotalHits) //minimal quality cuts
269  jetMomentum += (*itTrack)->momentum();
270  }
271 
273  vector<TransientTrack> transientTracks;
274 
275  for(typename Container::const_iterator itTrack = tracks.begin();
276  itTrack != tracks.end(); ++itTrack) {
277  TransientTrack transientTrack = builder->build(*itTrack);
278  const Track & track = transientTrack.track(); //**itTrack;
279  /* cout << " pt " << track.pt() <<
280  " d0 " << fabs(track.d0()) <<
281  " #hit " << track.hitPattern().numberOfValidHits()<<
282  " ipZ " << fabs(track.dz()-pv->z())<<
283  " chi2 " << track.normalizedChi2()<<
284  " #pixel " << track.hitPattern().numberOfValidPixelHits()<< endl;
285 */
286  if (track.pt() > m_cutMinPt &&
287  track.hitPattern().numberOfValidHits() >= m_cutTotalHits && // min num tracker hits
288  track.hitPattern().numberOfValidPixelHits() >= m_cutPixelHits &&
289  track.normalizedChi2() < m_cutMaxChiSquared &&
290  std::abs(track.dxy(pv->position())) < m_cutMaxTIP &&
291  std::abs(track.dz(pv->position())) < m_cutMaxLIP) {
292 // std::cout << "selected" << std::endl;
293  selectedTracks.push_back(*itTrack);
294  transientTracks.push_back(transientTrack);
295  }
296  }
297 // std::cout <<"SIZE: " << transientTracks.size() << std::endl;
298  GlobalVector direction(jetMomentum.x(), jetMomentum.y(), jetMomentum.z());
299 
300  auto_ptr<GhostTrack> ghostTrack;
301  TrackRef ghostTrackRef;
302  if (m_computeGhostTrack) {
303  GhostTrackFitter fitter;
306  ghostTrack.reset(new GhostTrack(fitter.fit(origin, error, direction,
307  m_ghostTrackPriorDeltaR,
308  transientTracks)));
309 
310 /*
311  if (std::sqrt(jetMomentum.Perp2()) > 30) {
312  double offset = ghostTrack->prediction().lambda(origin);
313  std::cout << "------------------ jet pt " << std::sqrt(jetMomentum.Perp2()) << std::endl;
314  const std::vector<GhostTrackState> *states = &ghostTrack->states();
315  for(std::vector<GhostTrackState>::const_iterator state = states->begin();
316  state != states->end(); ++state) {
317  double dist = state->lambda() - offset;
318  double err = state->lambdaError(ghostTrack->prediction(), error);
319  double ipSig = IPTools::signedImpactParameter3D(state->track(), direction, *pv).second.significance();
320  double axisDist = state->axisDistance(ghostTrack->prediction());
321  std::cout << state->track().impactPointState().freeState()->momentum().perp()
322  << ": " << dist << "/" << err << " [" << (dist / err) << "], ipsig = " << ipSig << ", dist = " << axisDist << ", w = " << state->weight() << std::endl;
323  }
324  }
325 */
326  ghostTrackRef = TrackRef(ghostTrackRefProd, ghostTracks->size());
327  ghostTracks->push_back(*ghostTrack);
328 
329  if (m_directionWithGhostTrack) {
330  const GhostTrackPrediction &pred = ghostTrack->prediction();
331  double lambda = pred.lambda(origin);
332  dummy = Vertex(RecoVertex::convertPos(pred.position(lambda)),
334  0, 0, 0);
335  pv = &dummy;
336  direction = pred.direction();
337  }
338  }
339 
340  vector<float> prob2D, prob3D;
341  vector<reco::btag::TrackIPData> ipData;
342 
343  for(unsigned int ind = 0; ind < transientTracks.size(); ind++) {
344  const TransientTrack &transientTrack = transientTracks[ind];
345  const Track & track = transientTrack.track();
346 
347  reco::btag::TrackIPData trackIP;
348  trackIP.ip3d = IPTools::signedImpactParameter3D(transientTrack, direction, *pv).second;
349  trackIP.ip2d = IPTools::signedTransverseImpactParameter(transientTrack, direction, *pv).second;
350 
351  TrajectoryStateOnSurface closest =
353  *pv, direction,
354  transientTrack.field());
355  if (closest.isValid())
356  trackIP.closestToJetAxis = closest.globalPosition();
357 
358  // TODO: cross check if it is the same using other methods
359  trackIP.distanceToJetAxis = IPTools::jetTrackDistance(transientTrack, direction, *pv).second;
360 
361  if (ghostTrack.get()) {
362  const std::vector<GhostTrackState> &states = ghostTrack->states();
363  std::vector<GhostTrackState>::const_iterator pos =
364  std::find_if(states.begin(), states.end(),
365  bind(equal_to<TransientTrack>(),
366  bind(&GhostTrackState::track, _1),
367  transientTrack));
368 
369  if (pos != states.end() && pos->isValid()) {
370  VertexDistance3D dist;
371  const 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),
378  VertexState(p2, e2));
379  trackIP.ghostTrackWeight = pos->weight();
380  } else {
381  trackIP.distanceToGhostTrack = Measurement1D(-1. -1.);
382  trackIP.ghostTrackWeight = 0.;
383  }
384  } else {
385  trackIP.distanceToGhostTrack = Measurement1D(-1. -1.);
386  trackIP.ghostTrackWeight = 1.;
387  }
388 
389  ipData.push_back(trackIP);
390 
391  if (m_computeProbabilities) {
392  //probability with 3D ip
393  pair<bool,double> probability = m_probabilityEstimator->probability(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(m_useTrackQuality,1,ipData.back().ip2d.significance(),track,*(it->jet()),*pv);
398  prob2D.push_back(probability.first ? probability.second : -1.);
399  }
400  }
401 
402  result->push_back(typename Product::value_type(ipData, prob2D, prob3D, selectedTracks,
403  *it, pvRef, direction, ghostTrackRef));
404  }
405 
406  if (m_computeGhostTrack)
407  iEvent.put(ghostTracks, "ghostTracks");
408  iEvent.put(result);
409 }
410 
411 
418 
419 template <class Container, class Base, class Helper> void IPProducer<Container,Base,Helper>::checkEventSetup(const EventSetup & iSetup)
420  {
421  using namespace edm;
422  using namespace edm::eventsetup;
423 
424  const EventSetupRecord & re2D= iSetup.get<BTagTrackProbability2DRcd>();
425  const EventSetupRecord & re3D= iSetup.get<BTagTrackProbability3DRcd>();
426  unsigned long long cacheId2D= re2D.cacheIdentifier();
427  unsigned long long cacheId3D= re3D.cacheIdentifier();
428 
429  if(cacheId2D!=m_calibrationCacheId2D || cacheId3D!=m_calibrationCacheId3D ) //Calibration changed
430  {
431  //iSetup.get<BTagTrackProbabilityRcd>().get(calib);
433  iSetup.get<BTagTrackProbability2DRcd>().get(calib2DHandle);
435  iSetup.get<BTagTrackProbability3DRcd>().get(calib3DHandle);
436 
437  const TrackProbabilityCalibration * ca2D= calib2DHandle.product();
438  const TrackProbabilityCalibration * ca3D= calib3DHandle.product();
439 
440  m_probabilityEstimator.reset(new HistogramProbabilityEstimator(ca3D,ca2D));
441 
442  }
443  m_calibrationCacheId3D=cacheId3D;
444  m_calibrationCacheId2D=cacheId2D;
445 }
446 
447 
448 #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:150
int i
Definition: DBlmapReader.cc:9
bool m_directionWithGhostTrack
Definition: IPProducer.h:153
const reco::Track * toTrack(const reco::TrackRef &t)
Definition: IPTagInfo.h:24
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:143
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:537
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
int numberOfValidHits() const
Definition: HitPattern.h:737
Base class for all types of Jets.
Definition: Jet.h:20
Measurement1D ip2d
Definition: IPTagInfo.h:30
double m_cutMaxTIP
Definition: IPProducer.h:148
unsigned long long m_calibrationCacheId2D
Definition: IPProducer.h:142
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:149
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:154
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:141
std::vector< reco::CandidatePtr > tracks(edm::Event &, const reco::JetTagInfo &it)
Definition: IPProducer.h:89
const Point & position() const
position
Definition: Vertex.h:106
int m_cutTotalHits
Definition: IPProducer.h:147
std::vector< reco::IPTagInfo< Container, Base > > Product
Definition: IPProducer.h:125
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:140
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:93
std::vector< std::vector< reco::CandidatePtr > > m_map
Definition: IPProducer.h:116
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:419
int iEvent
Definition: GenABIO.cc:230
Measurement1D distanceToGhostTrack
Definition: IPTagInfo.h:33
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:146
Helper m_helper
Definition: IPProducer.h:155
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:597
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:136
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:796
Measurement1D ip3d
Definition: IPTagInfo.h:31
Container::value_type value_type
double p2[4]
Definition: TauolaWrapper.h:90
bool m_directionWithTracks
Definition: IPProducer.h:152
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:585
const edm::ParameterSet & m_config
Definition: IPProducer.h:135
GlobalPoint closestToJetAxis
Definition: IPTagInfo.h:28
const Track & track() const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:364
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:117
T const * product() const
Definition: ESHandle.h:86
Error error() const
return SMatrix
Definition: Vertex.h:129
Measurement1D distanceToJetAxis
Definition: IPTagInfo.h:32
const Vector & prediction() const
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: IPProducer.h:218
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:183
GlobalPoint position(double lambda=0.) const
GlobalPoint closestToGhostTrack
Definition: IPTagInfo.h:29
edm::EDGetTokenT< edm::View< reco::Candidate > > token_cands
Definition: IPProducer.h:118
int numberOfValidPixelHits() const
Definition: HitPattern.h:752
bool m_useDB
Definition: IPProducer.h:144
TrajectoryStateOnSurface impactPointState() const
bool m_computeProbabilities
Definition: IPProducer.h:138
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:567
bool m_computeGhostTrack
Definition: IPProducer.h:139
double m_cutMaxLIP
Definition: IPProducer.h:151
std::vector< reco::JTATagInfo > makeBaseVector(edm::Event &iEvent)
Definition: IPProducer.h:68
GlobalError positionError(double lambda=0.) const
const GlobalVector direction() const