CMS 3D CMS Logo

CombinedSVComputer.h
Go to the documentation of this file.
1 #ifndef RecoBTag_SecondaryVertex_CombinedSVComputer_h
2 #define RecoBTag_SecondaryVertex_CombinedSVComputer_h
3 
4 #include <iostream>
5 #include <cstddef>
6 #include <string>
7 #include <cmath>
8 #include <vector>
9 
10 #include <Math/VectorUtil.h>
11 
15 
31 
36 
37 #define range_for(i, x) for (int i = (x).begin; i != (x).end; i += (x).increment)
38 
40 public:
42  virtual ~CombinedSVComputer() = default;
44  const reco::SecondaryVertexTagInfo &svInfo) const;
46  const reco::CandSecondaryVertexTagInfo &svInfo) const;
47 
48  struct IterationRange {
50  };
51  double flipValue(double value, bool vertex) const;
52  IterationRange flipIterate(int size, bool vertex) const;
53 
55 
56  const reco::btag::TrackIPData &threshTrack(const reco::CandIPTagInfo &trackIPTagInfo,
58  const reco::Jet &jet,
59  const GlobalPoint &pv) const;
60  const reco::btag::TrackIPData &threshTrack(const reco::TrackIPTagInfo &trackIPTagInfo,
62  const reco::Jet &jet,
63  const GlobalPoint &pv) const;
64  template <class SVTI, class IPTI>
66  reco::TrackKinematics &vertexKinematics,
67  const IPTI &ipInfo,
68  const SVTI &svInfo,
69  double &vtx_track_ptSum,
70  double &vtx_track_ESum) const;
71 
72 private:
73  bool trackFlip;
74  bool vertexFlip;
75  double charmCut;
80  unsigned int pseudoMultiplicityMin;
81  unsigned int trackMultiplicityMin;
87  std::vector<reco::btau::TaggingVariableName> taggingVariables;
88 };
89 
90 template <class SVTI, class IPTI>
92  reco::TrackKinematics &vertexKinematics,
93  const IPTI &ipInfo,
94  const SVTI &svInfo,
95  double &vtx_track_ptSum,
96  double &vtx_track_ESum) const {
97  using namespace ROOT::Math;
98  using namespace reco;
99 
100  typedef typename IPTI::input_container Container;
101  typedef typename Container::value_type TrackRef;
102 
103  edm::RefToBase<Jet> jet = ipInfo.jet();
104  math::XYZVector jetDir = jet->momentum().Unit();
105  bool havePv = ipInfo.primaryVertex().isNonnull();
106  GlobalPoint pv;
107  if (havePv)
108  pv = GlobalPoint(ipInfo.primaryVertex()->x(), ipInfo.primaryVertex()->y(), ipInfo.primaryVertex()->z());
109 
111 
112  vars.insert(btau::jetPt, jet->pt(), true);
113  vars.insert(btau::jetEta, jet->eta(), true);
114  vars.insert(btau::jetAbsEta, fabs(jet->eta()), true);
115 
116  if (ipInfo.selectedTracks().size() < trackMultiplicityMin)
117  return;
118 
119  vars.insert(btau::jetNTracks, ipInfo.selectedTracks().size(), true);
120 
121  TrackKinematics allKinematics;
122  TrackKinematics trackJetKinematics;
123 
124  double jet_track_ESum = 0.;
125 
126  int vtx = -1;
127 
128  IterationRange range = flipIterate(svInfo.nVertices(), true);
129  range_for(i, range) {
130  if (vtx < 0)
131  vtx = i;
132  }
133 
134  if (vtx >= 0) {
135  vtxType = btag::Vertices::RecoVertex;
136  vars.insert(btau::flightDistance1dVal, flipValue(svInfo.flightDistance(vtx, 1).value(), true), true);
137  vars.insert(btau::flightDistance1dSig, flipValue(svInfo.flightDistance(vtx, 1).significance(), true), true);
138  vars.insert(btau::flightDistance2dVal, flipValue(svInfo.flightDistance(vtx, 2).value(), true), true);
139  vars.insert(btau::flightDistance2dSig, flipValue(svInfo.flightDistance(vtx, 2).significance(), true), true);
140  vars.insert(btau::flightDistance3dVal, flipValue(svInfo.flightDistance(vtx, 3).value(), true), true);
141  vars.insert(btau::flightDistance3dSig, flipValue(svInfo.flightDistance(vtx, 3).significance(), true), true);
142  vars.insert(btau::vertexJetDeltaR, Geom::deltaR(svInfo.flightDirection(vtx), vertexFlip ? -jetDir : jetDir), true);
143  vars.insert(btau::jetNSecondaryVertices, svInfo.nVertices(), true);
144  }
145 
146  std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
147  const std::vector<reco::btag::TrackIPData> &ipData = ipInfo.impactParameterData();
148 
149  const Container &tracks = ipInfo.selectedTracks();
150  std::vector<TrackRef> pseudoVertexTracks;
151 
152  std::vector<TrackRef> trackPairV0Test(2);
153  range = flipIterate(indices.size(), false);
154  range_for(i, range) {
155  std::size_t idx = indices[i];
156  const reco::btag::TrackIPData &data = ipData[idx];
157  const TrackRef &track = tracks[idx];
158 
159  jet_track_ESum += std::sqrt(track->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
160 
161  // add track to kinematics for all tracks in jet
162  //allKinematics.add(track); // would make more sense for some variables, e.g. vertexEnergyRatio nicely between 0 and 1, but not necessarily the best option for the discriminating power...
163 
164  // filter track -> this track selection can be tighter than the vertex track selection (used to fill the track related variables...)
165  if (!trackSelector(track, data, *jet, pv))
166  continue;
167 
168  // add track to kinematics for all tracks in jet
169  allKinematics.add(track);
170 
171  // if no vertex was reconstructed, attempt pseudo vertex
172  if (vtxType == btag::Vertices::NoVertex && trackPseudoSelector(track, data, *jet, pv)) {
173  pseudoVertexTracks.push_back(track);
174  vertexKinematics.add(track);
175  }
176 
177  // check against all other tracks for V0 track pairs
178  trackPairV0Test[0] = track;
179  bool ok = true;
180  range_for(j, range) {
181  if (i == j)
182  continue;
183 
184  std::size_t pairIdx = indices[j];
185  const reco::btag::TrackIPData &pairTrackData = ipData[pairIdx];
186  const TrackRef &pairTrack = tracks[pairIdx];
187 
188  if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
189  continue;
190 
191  trackPairV0Test[1] = pairTrack;
192  if (!trackPairV0Filter(trackPairV0Test)) {
193  ok = false;
194  break;
195  }
196  }
197  if (!ok)
198  continue;
199 
200  trackJetKinematics.add(track);
201 
202  // add track variables
203  math::XYZVector trackMom = track->momentum();
204  double trackMag = std::sqrt(trackMom.Mag2());
205 
206  vars.insert(btau::trackSip3dVal, flipValue(data.ip3d.value(), false), true);
207  vars.insert(btau::trackSip3dSig, flipValue(data.ip3d.significance(), false), true);
208  vars.insert(btau::trackSip2dVal, flipValue(data.ip2d.value(), false), true);
209  vars.insert(btau::trackSip2dSig, flipValue(data.ip2d.significance(), false), true);
210  vars.insert(btau::trackJetDistVal, data.distanceToJetAxis.value(), true);
211  // vars.insert(btau::trackJetDistSig, data.distanceToJetAxis.significance(), true);
212  // vars.insert(btau::trackFirstTrackDist, data.distanceToFirstTrack, true);
213  // vars.insert(btau::trackGhostTrackVal, data.distanceToGhostTrack.value(), true);
214  // vars.insert(btau::trackGhostTrackSig, data.distanceToGhostTrack.significance(), true);
215  vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToJetAxis - pv).mag() : -1.0, true);
216 
217  vars.insert(btau::trackMomentum, trackMag, true);
218  vars.insert(btau::trackEta, trackMom.Eta(), true);
219 
220  // check for erroneous Perp^2 values before evaluating Perp (fix for DeepCSV NaN outputs)
221  double perp_trackMom_jetDir = VectorUtil::Perp2(trackMom, jetDir);
222  perp_trackMom_jetDir = (perp_trackMom_jetDir > 0. ? std::sqrt(perp_trackMom_jetDir) : 0.);
223 
224  vars.insert(btau::trackPtRel, perp_trackMom_jetDir, true);
225  vars.insert(btau::trackPPar, jetDir.Dot(trackMom), true);
226  vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
227  vars.insert(btau::trackPtRatio, perp_trackMom_jetDir / trackMag, true);
228  vars.insert(btau::trackPParRatio, jetDir.Dot(trackMom) / trackMag, true);
229  }
230 
231  if (vtxType == btag::Vertices::NoVertex && vertexKinematics.numberOfTracks() >= pseudoMultiplicityMin &&
232  pseudoVertexV0Filter(pseudoVertexTracks)) {
234  for (typename std::vector<TrackRef>::const_iterator trkIt = pseudoVertexTracks.begin();
235  trkIt != pseudoVertexTracks.end();
236  ++trkIt) {
237  vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir, (*trkIt)->momentum()), true);
238  vtx_track_ptSum += std::sqrt((*trkIt)->momentum().Perp2());
239  vtx_track_ESum += std::sqrt((*trkIt)->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
240  }
241  }
242 
243  vars.insert(btau::vertexCategory, vtxType, true);
244 
245  vars.insert(btau::trackJetPt, trackJetKinematics.vectorSum().Pt(), true);
246  vars.insert(btau::trackSumJetDeltaR, VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
247  vars.insert(btau::trackSumJetEtRatio, allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
248 
250  flipValue(threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.significance(), false),
251  true);
253  flipValue(threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.value(), false),
254  true);
256  flipValue(threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.significance(), false),
257  true);
259  flipValue(threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.value(), false),
260  true);
261 
262  if (vtxType != btag::Vertices::NoVertex) {
263  math::XYZTLorentzVector allSum = useTrackWeights ? allKinematics.weightedVectorSum() : allKinematics.vectorSum();
264  math::XYZTLorentzVector vertexSum =
265  useTrackWeights ? vertexKinematics.weightedVectorSum() : vertexKinematics.vectorSum();
266 
267  if (vtxType != btag::Vertices::RecoVertex) {
268  vars.insert(btau::vertexNTracks, vertexKinematics.numberOfTracks(), true);
269  vars.insert(btau::vertexJetDeltaR, VectorUtil::DeltaR(vertexSum, jetDir), true);
270  }
271 
272  double vertexMass = vertexSum.M();
274  const GlobalVector &dir = svInfo.flightDirection(vtx);
275  double vertexPt2 = math::XYZVector(dir.x(), dir.y(), dir.z()).Cross(vertexSum).Mag2() / dir.mag2();
276  vertexMass = std::sqrt(vertexMass * vertexMass + vertexPt2) + std::sqrt(vertexPt2);
277  }
278  vars.insert(btau::vertexMass, vertexMass, true);
279 
280  double varPi = (vertexMass / 5.2794) *
281  (vtx_track_ESum / jet_track_ESum); // 5.2794 should be the average B meson mass of PDG! CHECK!!!
282  vars.insert(btau::massVertexEnergyFraction, varPi / (varPi + 0.04), true);
283  double varB = (std::sqrt(5.2794) * vtx_track_ptSum) / (vertexMass * std::sqrt(jet->pt()));
284  vars.insert(btau::vertexBoostOverSqrtJetPt, varB * varB / (varB * varB + 10.), true);
285 
286  if (allKinematics.numberOfTracks()) {
287  vars.insert(btau::vertexEnergyRatio, vertexSum.E() / allSum.E(), true);
288  } else {
289  vars.insert(btau::vertexEnergyRatio, 1, true);
290  }
291  }
292 
293  reco::PFJet const *pfJet = dynamic_cast<reco::PFJet const *>(&*jet);
294  pat::Jet const *patJet = dynamic_cast<pat::Jet const *>(&*jet);
295  if (pfJet != nullptr) {
298  vars.insert(btau::photonEnergyFraction, pfJet->photonEnergyFraction(), true);
300  vars.insert(btau::muonEnergyFraction, pfJet->muonEnergyFraction(), true);
303  vars.insert(btau::photonMultiplicity, pfJet->photonMultiplicity(), true);
304  vars.insert(btau::electronMultiplicity, pfJet->electronMultiplicity(), true);
305  vars.insert(btau::muonMultiplicity, pfJet->muonMultiplicity(), true);
306  vars.insert(
310  true);
313  pfJet->electronMultiplicity() + pfJet->muonMultiplicity(),
314  true);
315 
316  } else if (patJet != nullptr && patJet->isPFJet()) {
319  vars.insert(btau::photonEnergyFraction, patJet->photonEnergyFraction(), true);
321  vars.insert(btau::muonEnergyFraction, patJet->muonEnergyFraction(), true);
324  vars.insert(btau::photonMultiplicity, patJet->photonMultiplicity(), true);
325  vars.insert(btau::electronMultiplicity, patJet->electronMultiplicity(), true);
326  vars.insert(btau::muonMultiplicity, patJet->muonMultiplicity(), true);
327  vars.insert(
329  vars.insert(
332  true);
335  patJet->photonMultiplicity() + patJet->electronMultiplicity() + patJet->muonMultiplicity(),
336  true);
337 
338  } else {
339  vars.insert(btau::chargedHadronEnergyFraction, 0., true);
340  vars.insert(btau::neutralHadronEnergyFraction, 0., true);
341  vars.insert(btau::photonEnergyFraction, 0., true);
342  vars.insert(btau::electronEnergyFraction, 0., true);
343  vars.insert(btau::muonEnergyFraction, 0., true);
344  vars.insert(btau::chargedHadronMultiplicity, 0, true);
345  vars.insert(btau::neutralHadronMultiplicity, 0, true);
346  vars.insert(btau::photonMultiplicity, 0, true);
347  vars.insert(btau::electronMultiplicity, 0, true);
348  vars.insert(btau::muonMultiplicity, 0, true);
349  vars.insert(btau::hadronMultiplicity, 0, true);
350  vars.insert(btau::hadronPhotonMultiplicity, 0, true);
351  vars.insert(btau::totalMultiplicity, 0, true);
352  }
353 }
354 
355 #endif // RecoBTag_SecondaryVertex_CombinedSVComputer_h
size
Write out results.
reco::TrackSelector trackSelector
virtual reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
const double piPlus
Definition: ParticleMasses.h:9
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: Jet.h:445
virtual ~CombinedSVComputer()=default
unsigned int pseudoMultiplicityMin
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:124
float electronEnergyFraction() const
electronEnergyFraction
Definition: PFJet.h:109
void add(const reco::Track &track, double weight=1.0)
Base class for all types of Jets.
Definition: Jet.h:20
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
reco::TrackSelector trackNoDeltaRSelector
reco::V0Filter pseudoVertexV0Filter
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
double flipValue(double value, bool vertex) const
CombinedSVComputer(const edm::ParameterSet &params)
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:402
float photonEnergyFraction() const
photonEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:418
unsigned int numberOfTracks() const
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:128
Jets made from PFObjects.
Definition: PFJet.h:20
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
reco::btag::SortCriteria sortCriterium
float photonEnergyFraction() const
photonEnergyFraction
Definition: PFJet.h:105
IterationRange flipIterate(int size, bool vertex) const
reco::TrackSelector trackPseudoSelector
T sqrt(T t)
Definition: SSEVec.h:19
edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset) const
const math::XYZTLorentzVector & vectorSum() const
def pv(vc)
Definition: MetAnalyzer.py:7
int photonMultiplicity() const
photonMultiplicity
Definition: Jet.h:447
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:130
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:97
Definition: value.py:1
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:126
bool isPFJet() const
check to see if the jet is a reco::PFJet
Definition: Jet.h:275
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: Jet.h:443
reco::V0Filter trackPairV0Filter
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:132
std::vector< reco::btau::TaggingVariableName > taggingVariables
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:398
int electronMultiplicity() const
electronMultiplicity
Definition: Jet.h:449
float muonEnergyFraction() const
muonEnergyFraction
Definition: PFJet.h:113
#define range_for(i, x)
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:101
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
int muonMultiplicity() const
muonMultiplicity
Definition: Jet.h:730
Analysis-level calorimeter jet class.
Definition: Jet.h:77
float muonEnergyFraction() const
muonEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:430
fixed size matrix
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
const reco::btag::TrackIPData & threshTrack(const reco::CandIPTagInfo &trackIPTagInfo, const reco::btag::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
const math::XYZTLorentzVector & weightedVectorSum() const
vars
Definition: DeepTauId.cc:30
void fillCommonVariables(reco::TaggingVariableList &vars, reco::TrackKinematics &vertexKinematics, const IPTI &ipInfo, const SVTI &svInfo, double &vtx_track_ptSum, double &vtx_track_ESum) const
float electronEnergyFraction() const
electronEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:424
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
unsigned int trackMultiplicityMin