RecoBTag
FeatureTools
src
SecondaryVertexConverter.cc
Go to the documentation of this file.
1
#include "
RecoBTag/FeatureTools/interface/deep_helpers.h
"
2
#include "
DataFormats/BTauReco/interface/SecondaryVertexFeatures.h
"
3
4
#include "
DataFormats/JetReco/interface/Jet.h
"
5
#include "
DataFormats/Candidate/interface/VertexCompositePtrCandidate.h
"
6
#include "
DataFormats/VertexReco/interface/Vertex.h
"
7
8
#include "
RecoBTag/FeatureTools/interface/SecondaryVertexConverter.h
"
9
10
namespace
btagbtvdeep
{
11
12
void
svToFeatures
(
const
reco::VertexCompositePtrCandidate
&
sv
,
13
const
reco::Vertex
&
pv
,
14
const
reco::Jet
&
jet
,
15
SecondaryVertexFeatures
&
sv_features
,
16
const
bool
flip
) {
17
math::XYZVector
jet_dir =
jet
.momentum().Unit();
18
sv_features
.pt =
sv
.pt();
19
sv_features
.deltaR =
catch_infs_and_bound
(std::fabs(
reco::deltaR
(
sv
, jet_dir)) - 0.5, 0, -2, 0);
20
sv_features
.mass =
sv
.mass();
21
sv_features
.ntracks =
sv
.numberOfDaughters();
22
sv_features
.chi2 =
sv
.vertexChi2();
23
sv_features
.normchi2 =
catch_infs_and_bound
(
sv_features
.chi2 /
sv
.vertexNdof(), 1000, -1000, 1000);
24
const
auto
& dxy_meas =
vertexDxy
(
sv
,
pv
);
25
sv_features
.dxy = dxy_meas.value();
26
sv_features
.dxysig =
catch_infs_and_bound
(dxy_meas.value() / dxy_meas.error(), 0, -1, 800);
27
const
auto
& d3d_meas =
vertexD3d
(
sv
,
pv
);
28
sv_features
.d3d = d3d_meas.value();
29
sv_features
.d3dsig =
catch_infs_and_bound
(d3d_meas.value() / d3d_meas.error(), 0, -1, 800);
30
sv_features
.costhetasvpv = (
flip
? -1.f : 1.f) *
vertexDdotP
(
sv
,
pv
);
31
sv_features
.enratio =
sv
.energy() /
jet
.energy();
32
}
33
34
}
// namespace btagbtvdeep
btagbtvdeep
Definition:
BoostedDoubleSVTagInfoFeatures.h:4
SecondaryVertexConverter.h
reco::Jet
Base class for all types of Jets.
Definition:
Jet.h:20
SecondaryVertexFeatures.h
btagbtvdeep::vertexDxy
Measurement1D vertexDxy(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv)
Definition:
deep_helpers.cc:49
Jet.h
reco::VertexCompositePtrCandidate
Definition:
VertexCompositePtrCandidate.h:16
pfNegativeDeepFlavourTagInfos_cfi.flip
flip
Definition:
pfNegativeDeepFlavourTagInfos_cfi.py:8
btagbtvdeep::SecondaryVertexFeatures
Definition:
SecondaryVertexFeatures.h:6
pfDeepBoostedJetPreprocessParams_cfi.sv
sv
Definition:
pfDeepBoostedJetPreprocessParams_cfi.py:226
VertexCompositePtrCandidate.h
btagbtvdeep::catch_infs_and_bound
const float catch_infs_and_bound(const float in, const float replace_value, const float lowerbound, const float upperbound, const float offset=0., const bool use_offsets=true)
Definition:
deep_helpers.cc:32
Vertex.h
pfParticleNetPreprocessParams_cfi.sv_features
sv_features
Definition:
pfParticleNetPreprocessParams_cfi.py:174
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition:
Vector3D.h:31
MetAnalyzer.pv
def pv(vc)
Definition:
MetAnalyzer.py:7
btagbtvdeep::svToFeatures
void svToFeatures(const reco::VertexCompositePtrCandidate &sv, const reco::Vertex &pv, const reco::Jet &jet, SecondaryVertexFeatures &sv_features, const bool flip=false)
Definition:
SecondaryVertexConverter.cc:12
metsig::jet
Definition:
SignAlgoResolutions.h:47
btagbtvdeep::vertexD3d
Measurement1D vertexD3d(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv)
Definition:
deep_helpers.cc:58
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition:
deltaR.h:30
deep_helpers.h
reco::Vertex
Definition:
Vertex.h:35
btagbtvdeep::vertexDdotP
float vertexDdotP(const reco::VertexCompositePtrCandidate &sv, const reco::Vertex &pv)
Definition:
deep_helpers.cc:67
Generated for CMSSW Reference Manual by
1.8.16