RecoBTag
SecondaryVertex
interface
CombinedSVSoftLeptonComputer.h
Go to the documentation of this file.
1
#ifndef RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
2
#define RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
3
4
#include "
DataFormats/BTauReco/interface/CandSoftLeptonTagInfo.h
"
5
#include "
DataFormats/BTauReco/interface/TaggingVariable.h
"
6
#include "
DataFormats/BTauReco/interface/TemplatedSoftLeptonTagInfo.h
"
7
#include "
DataFormats/BTauReco/interface/VertexTypes.h
"
8
#include "
DataFormats/TrackReco/interface/TrackFwd.h
"
9
10
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
11
12
#include "
RecoBTag/SecondaryVertex/interface/CombinedSVComputer.h
"
13
14
class
CombinedSVSoftLeptonComputer
:
public
CombinedSVComputer
{
15
public
:
16
explicit
CombinedSVSoftLeptonComputer
(
const
edm::ParameterSet
&
params
);
17
~CombinedSVSoftLeptonComputer
()
override
=
default
;
18
double
flipSoftLeptonValue
(
double
value
)
const
;
19
20
template
<
class
IPTI,
class
SVTI>
21
reco::TaggingVariableList
operator()
(
const
IPTI &ipInfo,
22
const
SVTI &svInfo,
23
const
reco::CandSoftLeptonTagInfo
&muonInfo,
24
const
reco::CandSoftLeptonTagInfo
&elecInfo)
const
;
25
26
private
:
27
bool
SoftLeptonFlip
;
28
};
29
30
double
CombinedSVSoftLeptonComputer::flipSoftLeptonValue
(
double
value
)
const
{
return
SoftLeptonFlip
? -
value
:
value
; }
31
32
template
<
class
IPTI,
class
SVTI>
33
reco::TaggingVariableList
CombinedSVSoftLeptonComputer::operator()
(
const
IPTI &ipInfo,
34
const
SVTI &svInfo,
35
const
reco::CandSoftLeptonTagInfo
&muonInfo,
36
const
reco::CandSoftLeptonTagInfo
&elecInfo)
const
{
37
using namespace
reco
;
38
39
// call the inherited operator()
40
TaggingVariableList
vars
=
CombinedSVComputer::operator()
(ipInfo, svInfo);
41
42
//Jets with vtxCategory 99 cause problems
43
unsigned
int
vtxType =
44
(
vars
.checkTag(
reco::btau::vertexCategory
) ? (
unsigned
int
)(
vars
.get(
reco::btau::vertexCategory
)) : 99);
45
if
(vtxType == 99)
46
return
vars
;
47
48
// the following is specific to soft leptons
49
int
leptonCategory = 0;
// 0 = no lepton, 1 = muon, 2 = electron
50
51
for
(
unsigned
int
i
= 0;
i
< muonInfo.
leptons
();
52
++
i
)
// loop over all muons, not optimal -> find the best or use ranking from best to worst
53
{
54
leptonCategory = 1;
// muon category
55
const
SoftLeptonProperties
&propertiesMuon = muonInfo.
properties
(
i
);
56
vars
.insert(
btau::leptonPtRel
, propertiesMuon.
ptRel
,
true
);
57
vars
.insert(
btau::leptonSip3d
,
flipSoftLeptonValue
(propertiesMuon.
sip3d
),
true
);
58
vars
.insert(
btau::leptonDeltaR
, propertiesMuon.
deltaR
,
true
);
59
vars
.insert(
btau::leptonRatioRel
, propertiesMuon.
ratioRel
,
true
);
60
vars
.insert(
btau::leptonEtaRel
, propertiesMuon.
etaRel
,
true
);
61
vars
.insert(
btau::leptonRatio
, propertiesMuon.
ratio
,
true
);
62
}
63
64
if
(leptonCategory != 1)
// no soft muon found, try soft electron
65
{
66
for
(
unsigned
int
i
= 0;
i
< elecInfo.
leptons
();
67
++
i
)
// loop over all electrons, not optimal -> find the best or use ranking from best to worst
68
{
69
leptonCategory = 2;
// electron category
70
const
SoftLeptonProperties
&propertiesElec = elecInfo.
properties
(
i
);
71
vars
.insert(
btau::leptonPtRel
, propertiesElec.
ptRel
,
true
);
72
vars
.insert(
btau::leptonSip3d
,
flipSoftLeptonValue
(propertiesElec.
sip3d
),
true
);
73
vars
.insert(
btau::leptonDeltaR
, propertiesElec.
deltaR
,
true
);
74
vars
.insert(
btau::leptonRatioRel
, propertiesElec.
ratioRel
,
true
);
75
vars
.insert(
btau::leptonEtaRel
, propertiesElec.
etaRel
,
true
);
76
vars
.insert(
btau::leptonRatio
, propertiesElec.
ratio
,
true
);
77
}
78
}
79
80
// set the default value for vertexLeptonCategory to 2 (= NoVertexNoSoftLepton)
81
int
vertexLepCat = 2;
82
83
if
(leptonCategory == 0)
// no soft lepton
84
{
85
if
(vtxType == (
unsigned
int
)(
btag::Vertices::RecoVertex
))
86
vertexLepCat = 0;
87
else
if
(vtxType == (
unsigned
int
)(
btag::Vertices::PseudoVertex
))
88
vertexLepCat = 1;
89
else
90
vertexLepCat = 2;
91
}
else
if
(leptonCategory == 1)
// soft muon
92
{
93
if
(vtxType == (
unsigned
int
)(
btag::Vertices::RecoVertex
))
94
vertexLepCat = 3;
95
else
if
(vtxType == (
unsigned
int
)(
btag::Vertices::PseudoVertex
))
96
vertexLepCat = 4;
97
else
98
vertexLepCat = 5;
99
}
else
if
(leptonCategory == 2)
// soft electron
100
{
101
if
(vtxType == (
unsigned
int
)(
btag::Vertices::RecoVertex
))
102
vertexLepCat = 6;
103
else
if
(vtxType == (
unsigned
int
)(
btag::Vertices::PseudoVertex
))
104
vertexLepCat = 7;
105
else
106
vertexLepCat = 8;
107
}
108
vars
.insert(
btau::vertexLeptonCategory
, vertexLepCat,
true
);
109
110
vars
.finalize();
111
return
vars
;
112
}
113
114
#endif // RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
CombinedSVSoftLeptonComputer::flipSoftLeptonValue
double flipSoftLeptonValue(double value) const
Definition:
CombinedSVSoftLeptonComputer.h:30
CombinedSVComputer
Definition:
CombinedSVComputer.h:39
mps_fire.i
i
Definition:
mps_fire.py:428
reco::SoftLeptonProperties
Definition:
TemplatedSoftLeptonTagInfo.h:15
reco::SoftLeptonProperties::sip3d
float sip3d
Definition:
TemplatedSoftLeptonTagInfo.h:34
reco::TemplatedSoftLeptonTagInfo::leptons
unsigned int leptons(void) const
Definition:
TemplatedSoftLeptonTagInfo.h:118
reco::btau::vertexCategory
Definition:
TaggingVariable.h:69
CalibrationSummaryClient_cfi.params
params
Definition:
CalibrationSummaryClient_cfi.py:14
reco::btau::leptonEtaRel
Definition:
TaggingVariable.h:124
reco::TaggingVariableList
Definition:
TaggingVariable.h:194
CombinedSVComputer.h
reco::btag::Vertices::PseudoVertex
Definition:
VertexTypes.h:18
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
CombinedSVSoftLeptonComputer::~CombinedSVSoftLeptonComputer
~CombinedSVSoftLeptonComputer() override=default
VertexTypes.h
reco::btau::leptonDeltaR
Definition:
TaggingVariable.h:125
reco::btau::vertexLeptonCategory
Definition:
TaggingVariable.h:70
reco::SoftLeptonProperties::etaRel
float etaRel
Definition:
TemplatedSoftLeptonTagInfo.h:38
TrackFwd.h
vars
vars
Definition:
DeepTauId.cc:164
reco::btag::Vertices::RecoVertex
Definition:
VertexTypes.h:18
reco::SoftLeptonProperties::ratio
float ratio
Definition:
TemplatedSoftLeptonTagInfo.h:40
reco::SoftLeptonProperties::ptRel
float ptRel
Definition:
TemplatedSoftLeptonTagInfo.h:35
CombinedSVSoftLeptonComputer::operator()
reco::TaggingVariableList operator()(const IPTI &ipInfo, const SVTI &svInfo, const reco::CandSoftLeptonTagInfo &muonInfo, const reco::CandSoftLeptonTagInfo &elecInfo) const
Definition:
CombinedSVSoftLeptonComputer.h:33
reco::TemplatedSoftLeptonTagInfo
Definition:
TemplatedSoftLeptonTagInfo.h:108
edm::ParameterSet
Definition:
ParameterSet.h:47
CombinedSVSoftLeptonComputer::SoftLeptonFlip
bool SoftLeptonFlip
Definition:
CombinedSVSoftLeptonComputer.h:27
createfilelist.int
int
Definition:
createfilelist.py:10
value
Definition:
value.py:1
reco::SoftLeptonProperties::deltaR
float deltaR
Definition:
TemplatedSoftLeptonTagInfo.h:39
TemplatedSoftLeptonTagInfo.h
CombinedSVSoftLeptonComputer::CombinedSVSoftLeptonComputer
CombinedSVSoftLeptonComputer(const edm::ParameterSet ¶ms)
Definition:
CombinedSVSoftLeptonComputer.cc:6
relativeConstraints.value
value
Definition:
relativeConstraints.py:53
reco::TemplatedSoftLeptonTagInfo::properties
const SoftLeptonProperties & properties(size_t i) const
Definition:
TemplatedSoftLeptonTagInfo.h:122
reco::btau::leptonSip3d
Definition:
TaggingVariable.h:121
ParameterSet.h
reco::btau::leptonRatio
Definition:
TaggingVariable.h:126
CandSoftLeptonTagInfo.h
TaggingVariable.h
reco::SoftLeptonProperties::ratioRel
float ratioRel
Definition:
TemplatedSoftLeptonTagInfo.h:41
reco::btau::leptonPtRel
Definition:
TaggingVariable.h:122
CombinedSVComputer::operator()
virtual reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
Definition:
CombinedSVComputer.cc:133
reco::btau::leptonRatioRel
Definition:
TaggingVariable.h:127
CombinedSVSoftLeptonComputer
Definition:
CombinedSVSoftLeptonComputer.h:14
Generated for CMSSW Reference Manual by
1.8.16