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 
38 #define range_for(i, x) \
39  for(int i = (x).begin; i != (x).end; i += (x).increment)
40 
41 
43  public:
44  explicit CombinedSVComputer(const edm::ParameterSet &params);
45  virtual ~CombinedSVComputer() = default;
47  operator () (const reco::TrackIPTagInfo &ipInfo,
48  const reco::SecondaryVertexTagInfo &svInfo) const;
50  operator () (const reco::CandIPTagInfo &ipInfo,
51  const reco::CandSecondaryVertexTagInfo &svInfo) const;
52 
53  struct IterationRange {
55  };
56  double flipValue(double value, bool vertex) const;
57  IterationRange flipIterate(int size, bool vertex) const;
58 
60 
62  threshTrack(const reco::CandIPTagInfo &trackIPTagInfo,
63  const reco::btag::SortCriteria sort,
64  const reco::Jet &jet,
65  const GlobalPoint &pv) const;
67  threshTrack(const reco::TrackIPTagInfo &trackIPTagInfo,
68  const reco::btag::SortCriteria sort,
69  const reco::Jet &jet,
70  const GlobalPoint &pv) const;
71  template <class SVTI,class IPTI>
73  const IPTI & ipInfo,const SVTI & svInfo,
74  double & vtx_track_ptSum, double & vtx_track_ESum) const;
75 
76  private:
77  bool trackFlip;
78  bool vertexFlip;
79  double charmCut;
84  unsigned int pseudoMultiplicityMin;
85  unsigned int trackMultiplicityMin;
91  std::vector<reco::btau::TaggingVariableName> taggingVariables;
92 };
93 
94 template <class SVTI,class IPTI>
96  const IPTI & ipInfo, const SVTI & svInfo,
97  double & vtx_track_ptSum, double & vtx_track_ESum) const
98 {
99  using namespace ROOT::Math;
100  using namespace reco;
101 
102  typedef typename IPTI::input_container Container;
103  typedef typename Container::value_type TrackRef;
104 
105  edm::RefToBase<Jet> jet = ipInfo.jet();
106  math::XYZVector jetDir = jet->momentum().Unit();
107  bool havePv = ipInfo.primaryVertex().isNonnull();
108  GlobalPoint pv;
109  if (havePv)
110  pv = GlobalPoint(ipInfo.primaryVertex()->x(),
111  ipInfo.primaryVertex()->y(),
112  ipInfo.primaryVertex()->z());
113 
115 
116 
117  vars.insert(btau::jetPt, jet->pt(), true);
118  vars.insert(btau::jetEta, jet->eta(), true);
119  vars.insert(btau::jetAbsEta, fabs(jet->eta()), true);
120 
121  if (ipInfo.selectedTracks().size() < trackMultiplicityMin)
122  return;
123 
124  vars.insert(btau::jetNTracks, ipInfo.selectedTracks().size(), true);
125 
126  TrackKinematics allKinematics;
127  TrackKinematics trackJetKinematics;
128 
129  double jet_track_ESum= 0.;
130 
131  int vtx = -1;
132 
133  IterationRange range = flipIterate(svInfo.nVertices(), true);
134  range_for(i , range) {
135  if (vtx < 0) vtx = i;
136  }
137 
138  if (vtx >= 0) {
139  vtxType = btag::Vertices::RecoVertex;
140  vars.insert(btau::flightDistance1dVal,flipValue(svInfo.flightDistance(vtx, 1).value(),true),true);
141  vars.insert(btau::flightDistance1dSig,flipValue(svInfo.flightDistance(vtx, 1).significance(),true),true);
142  vars.insert(btau::flightDistance2dVal,flipValue(svInfo.flightDistance(vtx, 2).value(),true),true);
143  vars.insert(btau::flightDistance2dSig,flipValue(svInfo.flightDistance(vtx, 2).significance(),true),true);
144  vars.insert(btau::flightDistance3dVal,flipValue(svInfo.flightDistance(vtx, 3).value(),true),true);
145  vars.insert(btau::flightDistance3dSig,flipValue(svInfo.flightDistance(vtx, 3).significance(),true),true);
146  vars.insert(btau::vertexJetDeltaR,Geom::deltaR(svInfo.flightDirection(vtx), vertexFlip ? -jetDir : jetDir),true);
147  vars.insert(btau::jetNSecondaryVertices, svInfo.nVertices(), true);
148  }
149 
150  std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
151  const std::vector<reco::btag::TrackIPData> &ipData = ipInfo.impactParameterData();
152 
153  const Container &tracks = ipInfo.selectedTracks();
154  std::vector<TrackRef> pseudoVertexTracks;
155 
156  std::vector<TrackRef> trackPairV0Test(2);
157  range = flipIterate(indices.size(), false);
158  range_for(i, range) {
159  std::size_t idx = indices[i];
160  const reco::btag::TrackIPData &data = ipData[idx];
161  const TrackRef &track = tracks[idx];
162 
163  jet_track_ESum += std::sqrt(track->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
164 
165  // add track to kinematics for all tracks in jet
166  //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...
167 
168  // filter track -> this track selection can be tighter than the vertex track selection (used to fill the track related variables...)
169  if (!trackSelector(track, data, *jet, pv))
170  continue;
171 
172  // add track to kinematics for all tracks in jet
173  allKinematics.add(track);
174 
175  // if no vertex was reconstructed, attempt pseudo vertex
176  if (vtxType == btag::Vertices::NoVertex && trackPseudoSelector(track, data, *jet, pv)) {
177  pseudoVertexTracks.push_back(track);
178  vertexKinematics.add(track);
179  }
180 
181  // check against all other tracks for V0 track pairs
182  trackPairV0Test[0] = track;
183  bool ok = true;
184  range_for(j, range) {
185  if (i == j)
186  continue;
187 
188  std::size_t pairIdx = indices[j];
189  const reco::btag::TrackIPData &pairTrackData = ipData[pairIdx];
190  const TrackRef &pairTrack = tracks[pairIdx];
191 
192  if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
193  continue;
194 
195  trackPairV0Test[1] = pairTrack;
196  if (!trackPairV0Filter(trackPairV0Test)) {
197  ok = false;
198  break;
199  }
200  }
201  if (!ok)
202  continue;
203 
204  trackJetKinematics.add(track);
205 
206  // add track variables
207  math::XYZVector trackMom = track->momentum();
208  double trackMag = std::sqrt(trackMom.Mag2());
209 
210  vars.insert(btau::trackSip3dVal, flipValue(data.ip3d.value(), false), true);
211  vars.insert(btau::trackSip3dSig, flipValue(data.ip3d.significance(), false), true);
212  vars.insert(btau::trackSip2dVal, flipValue(data.ip2d.value(), false), true);
213  vars.insert(btau::trackSip2dSig, flipValue(data.ip2d.significance(), false), true);
215 // vars.insert(btau::trackJetDistSig, data.distanceToJetAxis.significance(), true);
216 // vars.insert(btau::trackFirstTrackDist, data.distanceToFirstTrack, true);
217 // vars.insert(btau::trackGhostTrackVal, data.distanceToGhostTrack.value(), true);
218 // vars.insert(btau::trackGhostTrackSig, data.distanceToGhostTrack.significance(), true);
219  vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToJetAxis - pv).mag() : -1.0, true);
220 
221  vars.insert(btau::trackMomentum, trackMag, true);
222  vars.insert(btau::trackEta, trackMom.Eta(), true);
223 
224  // check for erroneous Perp^2 values before evaluating Perp (fix for DeepCSV NaN outputs)
225  double perp_trackMom_jetDir = VectorUtil::Perp2(trackMom, jetDir);
226  perp_trackMom_jetDir = (perp_trackMom_jetDir > 0. ? std::sqrt(perp_trackMom_jetDir ) : 0. );
227 
228  vars.insert(btau::trackPtRel, perp_trackMom_jetDir, true);
229  vars.insert(btau::trackPPar, jetDir.Dot(trackMom), true);
230  vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
231  vars.insert(btau::trackPtRatio, perp_trackMom_jetDir / trackMag, true);
232  vars.insert(btau::trackPParRatio, jetDir.Dot(trackMom) / trackMag, true);
233  }
234 
235  if (vtxType == btag::Vertices::NoVertex && vertexKinematics.numberOfTracks() >= pseudoMultiplicityMin && pseudoVertexV0Filter(pseudoVertexTracks))
236  {
238  for(typename std::vector<TrackRef>::const_iterator trkIt = pseudoVertexTracks.begin(); trkIt != pseudoVertexTracks.end(); ++trkIt)
239  {
240  vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir,(*trkIt)->momentum()), true);
241  vtx_track_ptSum += std::sqrt((*trkIt)->momentum().Perp2());
242  vtx_track_ESum += std::sqrt((*trkIt)->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
243  }
244  }
245 
246  vars.insert(btau::vertexCategory, vtxType, true);
247 
248  vars.insert(btau::trackJetPt, trackJetKinematics.vectorSum().Pt(), true);
249  vars.insert(btau::trackSumJetDeltaR,VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
250  vars.insert(btau::trackSumJetEtRatio,allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
251 
252  vars.insert(btau::trackSip3dSigAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.significance(),false),true);
253  vars.insert(btau::trackSip3dValAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.value(),false),true);
254  vars.insert(btau::trackSip2dSigAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.significance(),false),true);
255  vars.insert(btau::trackSip2dValAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.value(),false),true);
256 
257  if (vtxType != btag::Vertices::NoVertex) {
259  ? allKinematics.weightedVectorSum()
260  : allKinematics.vectorSum();
262  ? vertexKinematics.weightedVectorSum()
263  : vertexKinematics.vectorSum();
264 
265  if (vtxType != btag::Vertices::RecoVertex) {
266  vars.insert(btau::vertexNTracks,vertexKinematics.numberOfTracks(), true);
267  vars.insert(btau::vertexJetDeltaR,VectorUtil::DeltaR(vertexSum, jetDir), true);
268  }
269 
270  double vertexMass = vertexSum.M();
271  if (vtxType == btag::Vertices::RecoVertex &&
273  const GlobalVector& dir = svInfo.flightDirection(vtx);
274  double vertexPt2 = math::XYZVector(dir.x(), dir.y(), dir.z()).Cross(vertexSum).Mag2() / dir.mag2();
275  vertexMass = std::sqrt(vertexMass * vertexMass + vertexPt2) + std::sqrt(vertexPt2);
276  }
277  vars.insert(btau::vertexMass, vertexMass, true);
278 
279  double varPi = (vertexMass/5.2794) * (vtx_track_ESum /jet_track_ESum); // 5.2794 should be the average B meson mass of PDG! CHECK!!!
280  vars.insert(btau::massVertexEnergyFraction, varPi / (varPi + 0.04), true);
281  double varB = (std::sqrt(5.2794) * vtx_track_ptSum) / ( vertexMass * std::sqrt(jet->pt()));
282  vars.insert(btau::vertexBoostOverSqrtJetPt,varB*varB/(varB*varB + 10.), true);
283 
284  if (allKinematics.numberOfTracks()) {
285  vars.insert(btau::vertexEnergyRatio, vertexSum.E() / allSum.E(), true);
286  }
287  else {
288  vars.insert(btau::vertexEnergyRatio, 1, true);
289  }
290  }
291 
292  reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>( &* jet ) ;
293  pat::Jet const * patJet = dynamic_cast<pat::Jet const *>( &* jet ) ;
294  if ( pfJet != nullptr )
295  {
305  vars.insert(btau::muonMultiplicity,pfJet->muonMultiplicity(), true);
309 
310  }
311  else if( patJet != nullptr && patJet->isPFJet() )
312  {
317  vars.insert(btau::muonEnergyFraction,patJet->muonEnergyFraction(), true);
320  vars.insert(btau::photonMultiplicity,patJet->photonMultiplicity(), true);
322  vars.insert(btau::muonMultiplicity,patJet->muonMultiplicity(), true);
326 
327  }
328  else
329  {
332  vars.insert(btau::photonEnergyFraction,0., true);
333  vars.insert(btau::electronEnergyFraction,0., true);
334  vars.insert(btau::muonEnergyFraction,0., true);
337  vars.insert(btau::photonMultiplicity,0, true);
338  vars.insert(btau::electronMultiplicity,0, true);
339  vars.insert(btau::muonMultiplicity,0, true);
340  vars.insert(btau::hadronMultiplicity,0, true);
342  vars.insert(btau::totalMultiplicity,0, true);
343  }
344 }
345 
346 
347 #endif // RecoBTag_SecondaryVertex_CombinedSVComputer_h
size
Write out results.
reco::TrackSelector trackSelector
const double piPlus
Definition: ParticleMasses.h:9
T mag2() const
Definition: PV3DBase.h:66
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:131
const math::XYZTLorentzVector & weightedVectorSum() const
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:374
float photonEnergyFraction() const
photonEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:384
int electronMultiplicity() const
electronMultiplicity
Definition: Jet.h:409
virtual ~CombinedSVComputer()=default
unsigned int pseudoMultiplicityMin
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:372
void add(const reco::Track &track, double weight=1.0)
int photonMultiplicity() const
photonMultiplicity
Definition: Jet.h:407
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:100
Base class for all types of Jets.
Definition: Jet.h:20
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
Measurement1D ip2d
Definition: IPTagInfo.h:31
reco::TrackSelector trackNoDeltaRSelector
reco::V0Filter pseudoVertexV0Filter
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
int muonMultiplicity() const
muonMultiplicity
Definition: Jet.h:683
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:340
virtual reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: Jet.h:403
CombinedSVComputer(const edm::ParameterSet &params)
float photonEnergyFraction() const
photonEnergyFraction
Definition: PFJet.h:108
Jets made from PFObjects.
Definition: PFJet.h:21
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
reco::btag::SortCriteria sortCriterium
reco::TrackSelector trackPseudoSelector
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
def pv(vc)
Definition: MetAnalyzer.py:7
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:129
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:104
Definition: value.py:1
Measurement1D ip3d
Definition: IPTagInfo.h:32
IterationRange flipIterate(int size, bool vertex) const
bool isPFJet() const
check to see if the jet is a reco::PFJet
Definition: Jet.h:255
reco::V0Filter trackPairV0Filter
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
std::vector< reco::btau::TaggingVariableName > taggingVariables
double significance() const
Definition: Measurement1D.h:32
unsigned int numberOfTracks() const
float electronEnergyFraction() const
electronEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:388
#define range_for(i, x)
GlobalPoint closestToJetAxis
Definition: IPTagInfo.h:29
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
double flipValue(double value, bool vertex) const
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
float electronEnergyFraction() const
electronEnergyFraction
Definition: PFJet.h:112
Analysis-level calorimeter jet class.
Definition: Jet.h:80
edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset) const
double value() const
Definition: Measurement1D.h:28
const reco::btag::TrackIPData & threshTrack(const reco::CandIPTagInfo &trackIPTagInfo, const reco::btag::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: PFJet.h:127
Measurement1D distanceToJetAxis
Definition: IPTagInfo.h:33
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:135
void fillCommonVariables(reco::TaggingVariableList &vars, reco::TrackKinematics &vertexKinematics, const IPTI &ipInfo, const SVTI &svInfo, double &vtx_track_ptSum, double &vtx_track_ESum) const
fixed size matrix
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: Jet.h:405
const math::XYZTLorentzVector & vectorSum() const
float muonEnergyFraction() const
muonEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:392
dbl *** dir
Definition: mlp_gen.cc:35
float muonEnergyFraction() const
muonEnergyFraction
Definition: PFJet.h:116
T x() const
Definition: PV3DBase.h:62
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:133
void insert(const TaggingVariable &variable, bool delayed=false)
unsigned int trackMultiplicityMin