CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GhostTrackComputer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cstddef>
3 #include <string>
4 #include <cmath>
5 #include <vector>
6 
7 #include <Math/VectorUtil.h>
8 
11 
24 
30 
32 
33 using namespace reco;
34 
36 {
37  edm::ParameterSet psetCopy(pset);
38  psetCopy.addParameter<double>("jetDeltaRMax", 99999.0);
39  return psetCopy;
40 }
41 
43  charmCut(params.getParameter<double>("charmCut")),
44  sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
45  trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
46  trackNoDeltaRSelector(dropDeltaR(params.getParameter<edm::ParameterSet>("trackSelection"))),
47  minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
48  trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter"))
49 {
50 }
51 
55  const reco::Jet &jet,
56  const GlobalPoint &pv) const
57 {
59  trackIPTagInfo.selectedTracks();
60  const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
61  trackIPTagInfo.impactParameterData();
62  std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
63 
64  TrackKinematics kin;
65  for(std::vector<std::size_t>::const_iterator iter = indices.begin();
66  iter != indices.end(); ++iter) {
67  const TrackIPTagInfo::TrackIPData &data = ipData[*iter];
68  const Track &track = *tracks[*iter];
69 
70  if (!trackNoDeltaRSelector(track, data, jet, pv))
71  continue;
72 
73  kin.add(track);
74  if (kin.vectorSum().M() > charmCut)
75  return data;
76  }
77 
78  static const TrackIPTagInfo::TrackIPData dummy = {
79  GlobalPoint(),
80  GlobalPoint(),
81  Measurement1D(-1.0, 1.0),
82  Measurement1D(-1.0, 1.0),
83  Measurement1D(-1.0, 1.0),
84  Measurement1D(-1.0, 1.0),
85  0.
86  };
87  return dummy;
88 }
89 
90 static double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
91 {
92  double momPar = dir.Dot(track);
93  double energy = std::sqrt(track.Mag2() +
94  ROOT::Math::Square(ParticleMasses::piPlus));
95 
96  return 0.5 * std::log((energy + momPar) / (energy - momPar));
97 }
98 
99 static void addMeas(std::pair<double, double> &sum, Measurement1D meas)
100 {
101  double weight = 1. / meas.error();
102  weight *= weight;
103  sum.first += weight * meas.value();
104  sum.second += weight;
105 }
106 
109  const SecondaryVertexTagInfo &svInfo) const
110 {
111  using namespace ROOT::Math;
112 
113  edm::RefToBase<Jet> jet = ipInfo.jet();
114  math::XYZVector jetDir = jet->momentum().Unit();
115  bool havePv = ipInfo.primaryVertex().isNonnull();
116  GlobalPoint pv;
117  if (havePv)
118  pv = GlobalPoint(ipInfo.primaryVertex()->x(),
119  ipInfo.primaryVertex()->y(),
120  ipInfo.primaryVertex()->z());
121 
123 
124  TaggingVariableList vars;
125 
126  vars.insert(btau::jetPt, jet->pt(), true);
127  vars.insert(btau::jetEta, jet->eta(), true);
128 
129  TrackKinematics allKinematics;
130  TrackKinematics vertexKinematics;
131  TrackKinematics trackKinematics;
132 
133  std::pair<double, double> vertexDist2D, vertexDist3D;
134  std::pair<double, double> tracksDist2D, tracksDist3D;
135 
136  unsigned int nVertices = 0;
137  unsigned int nVertexTracks = 0;
138  unsigned int nTracks = 0;
139  for(unsigned int i = 0; i < svInfo.nVertices(); i++) {
140  const Vertex &vertex = svInfo.secondaryVertex(i);
141  bool hasRefittedTracks = vertex.hasRefittedTracks();
143  unsigned int n = 0;
144  for(TrackRefVector::const_iterator track = tracks.begin();
145  track != tracks.end(); track++)
146  if (svInfo.trackWeight(i, *track) >= minTrackWeight)
147  n++;
148 
149  if (n < 1)
150  continue;
151  bool isTrackVertex = (n == 1);
152  ++*(isTrackVertex ? &nTracks : &nVertices);
153 
154  addMeas(*(isTrackVertex ? &tracksDist2D : &vertexDist2D),
155  svInfo.flightDistance(i, true));
156  addMeas(*(isTrackVertex ? &tracksDist3D : &vertexDist3D),
157  svInfo.flightDistance(i, false));
158 
159  TrackKinematics &kin = isTrackVertex ? trackKinematics
160  : vertexKinematics;
161  for(TrackRefVector::const_iterator track = tracks.begin();
162  track != tracks.end(); track++) {
163  float w = svInfo.trackWeight(i, *track);
164  if (w < minTrackWeight)
165  continue;
166  if (hasRefittedTracks) {
167  Track actualTrack =
168  vertex.refittedTrack(*track);
169  kin.add(actualTrack, w);
170  vars.insert(btau::trackEtaRel, etaRel(jetDir,
171  actualTrack.momentum()), true);
172  } else {
173  kin.add(**track, w);
174  vars.insert(btau::trackEtaRel, etaRel(jetDir,
175  (*track)->momentum()), true);
176  }
177  if (!isTrackVertex)
178  nVertexTracks++;
179  }
180  }
181 
182  Measurement1D dist2D, dist3D;
183  if (nVertices) {
184  vtxType = btag::Vertices::RecoVertex;
185 
186  if (nVertices == 1 && nTracks) {
187  vertexDist2D.first += tracksDist2D.first;
188  vertexDist2D.second += tracksDist2D.second;
189  vertexDist3D.first += tracksDist3D.first;
190  vertexDist3D.second += tracksDist3D.second;
191  vertexKinematics += trackKinematics;
192  }
193 
194  dist2D = Measurement1D(
195  vertexDist2D.first / vertexDist2D.second,
196  std::sqrt(1. / vertexDist2D.second));
197  dist3D = Measurement1D(
198  vertexDist3D.first / vertexDist3D.second,
199  std::sqrt(1. / vertexDist3D.second));
200 
201  vars.insert(btau::jetNSecondaryVertices, nVertices, true);
202  vars.insert(btau::vertexNTracks, nVertexTracks, true);
203  } else if (nTracks) {
205  vertexKinematics = trackKinematics;
206 
207  dist2D = Measurement1D(
208  tracksDist2D.first / tracksDist2D.second,
209  std::sqrt(1. / tracksDist2D.second));
210  dist3D = Measurement1D(
211  tracksDist3D.first / tracksDist3D.second,
212  std::sqrt(1. / tracksDist3D.second));
213  }
214 
215  if (nVertices || nTracks) {
216  vars.insert(btau::flightDistance2dVal, dist2D.value(), true);
217  vars.insert(btau::flightDistance2dSig, dist2D.significance(), true);
218  vars.insert(btau::flightDistance3dVal, dist3D.value(), true);
219  vars.insert(btau::flightDistance3dSig, dist3D.significance(), true);
220  vars.insert(btau::vertexJetDeltaR,
221  Geom::deltaR(svInfo.flightDirection(0), jetDir),true);
222  vars.insert(btau::jetNSingleTrackVertices, nTracks, true);
223  }
224 
225  std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
226  const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
227  ipInfo.impactParameterData();
229  ipInfo.selectedTracks();
230 
231  TrackRef trackPairV0Test[2];
232  for(unsigned int i = 0; i < indices.size(); i++) {
233  std::size_t idx = indices[i];
234  const TrackIPTagInfo::TrackIPData &data = ipData[idx];
235  const TrackRef &trackRef = tracks[idx];
236  const Track &track = *trackRef;
237 
238  // filter track
239 
240  if (!trackSelector(track, data, *jet, pv))
241  continue;
242 
243  // add track to kinematics for all tracks in jet
244 
245  allKinematics.add(track);
246 
247  // check against all other tracks for V0 track pairs
248 
249  trackPairV0Test[0] = tracks[idx];
250  bool ok = true;
251  for(unsigned int j = 0; j < indices.size(); j++) {
252  if (i == j)
253  continue;
254 
255  std::size_t pairIdx = indices[j];
256  const TrackIPTagInfo::TrackIPData &pairTrackData =
257  ipData[pairIdx];
258  const TrackRef &pairTrackRef = tracks[pairIdx];
259  const Track &pairTrack = *pairTrackRef;
260 
261  if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
262  continue;
263 
264  trackPairV0Test[1] = pairTrackRef;
265  if (!trackPairV0Filter(trackPairV0Test, 2)) {
266  ok = false;
267  break;
268  }
269  }
270  if (!ok)
271  continue;
272 
273  // add track variables
274 
275  math::XYZVector trackMom = track.momentum();
276  double trackMag = std::sqrt(trackMom.Mag2());
277 
278  vars.insert(btau::trackSip3dVal, data.ip3d.value(), true);
279  vars.insert(btau::trackSip3dSig, data.ip3d.significance(), true);
280  vars.insert(btau::trackSip2dVal, data.ip2d.value(), true);
281  vars.insert(btau::trackSip2dSig, data.ip2d.significance(), true);
282  vars.insert(btau::trackJetDistVal, data.distanceToJetAxis.value(), true);
283  vars.insert(btau::trackGhostTrackDistVal, data.distanceToGhostTrack.value(), true);
285  vars.insert(btau::trackGhostTrackWeight, data.ghostTrackWeight, true);
286  vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToGhostTrack - pv).mag() : -1.0, true);
287 
288  vars.insert(btau::trackMomentum, trackMag, true);
289  vars.insert(btau::trackEta, trackMom.Eta(), true);
290 
291  vars.insert(btau::trackChi2, track.normalizedChi2(), true);
292  vars.insert(btau::trackNPixelHits, track.hitPattern().pixelLayersWithMeasurement(), true);
293  vars.insert(btau::trackNTotalHits, track.hitPattern().trackerLayersWithMeasurement(), true);
294  vars.insert(btau::trackPtRel, VectorUtil::Perp(trackMom, jetDir), true);
295  vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
296  }
297 
298  vars.insert(btau::vertexCategory, vtxType, true);
299  vars.insert(btau::trackSumJetDeltaR,
300  VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
301  vars.insert(btau::trackSumJetEtRatio,
302  allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
303  vars.insert(btau::trackSip3dSigAboveCharm,
304  threshTrack(ipInfo, TrackIPTagInfo::IP3DSig, *jet, pv)
305  .ip3d.significance(),
306  true);
307  vars.insert(btau::trackSip2dSigAboveCharm,
308  threshTrack(ipInfo, TrackIPTagInfo::IP2DSig, *jet, pv)
309  .ip2d.significance(),
310  true);
311 
312  if (vtxType != btag::Vertices::NoVertex) {
313  math::XYZTLorentzVector allSum = allKinematics.vectorSum();
314  math::XYZTLorentzVector vertexSum = vertexKinematics.vectorSum();
315 
316  vars.insert(btau::vertexMass, vertexSum.M(), true);
317  if (allKinematics.numberOfTracks())
318  vars.insert(btau::vertexEnergyRatio,
319  vertexSum.E() / allSum.E(), true);
320  else
321  vars.insert(btau::vertexEnergyRatio, 1, true);
322  }
323 
324  vars.finalize();
325 
326  return vars;
327 }
virtual double et() const GCC11_FINAL
transverse energy
int i
Definition: DBlmapReader.cc:9
reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:148
reco::TrackIPTagInfo::SortCriteria getCriterium(const std::string &name)
Definition: TrackSorting.cc:11
Measurement1D flightDistance(unsigned int index, bool in2d=false) const
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:109
const edm::Ref< VertexCollection > & primaryVertex() const
void add(const reco::Track &track, double weight=1.0)
GhostTrackComputer(const edm::ParameterSet &params)
Track refittedTrack(const TrackBaseRef &track) const
Base class for all types of Jets.
Definition: Jet.h:20
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:120
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
int pixelLayersWithMeasurement() const
Definition: HitPattern.h:745
double error() const
Definition: Measurement1D.h:30
static edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset)
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
const GlobalVector & flightDirection(unsigned int index) const
const Vertex & secondaryVertex(unsigned int index) const
const reco::TrackIPTagInfo::TrackIPData & threshTrack(const reco::TrackIPTagInfo &trackIPTagInfo, const reco::TrackIPTagInfo::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
static void addMeas(std::pair< double, double > &sum, Measurement1D meas)
int trackerLayersWithMeasurement() const
Definition: HitPattern.h:740
T sqrt(T t)
Definition: SSEVec.h:48
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:142
int j
Definition: DBlmapReader.cc:9
float trackWeight(unsigned int svIndex, unsigned int trackindex) const
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:221
reco::V0Filter trackPairV0Filter
virtual edm::RefToBase< Jet > jet(void) const
returns a polymorphic reference to the tagged jet
Definition: JTATagInfo.h:20
TrackRefVector vertexTracks() const
std::vector< size_t > sortedIndexes(SortCriteria mode=IP3DSig) const
double significance() const
Definition: Measurement1D.h:32
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double value() const
Definition: Measurement1D.h:28
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double deltaR(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:84
const math::XYZTLorentzVector & vectorSum() const
const std::vector< TrackIPData > & impactParameterData() const
reco::TrackIPTagInfo::SortCriteria sortCriterium
reco::TrackSelector trackNoDeltaRSelector
reco::TrackSelector trackSelector
T w() const
dbl *** dir
Definition: mlp_gen.cc:35
int weight
Definition: histoStyle.py:50
const edm::RefVector< TrackCollection > & selectedTracks() const
static double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
void insert(const TaggingVariable &variable, bool delayed=false)
tuple trackSelector
Tracks selection.
Definition: valSkim_cff.py:4