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