RecoBTag
SecondaryVertex
interface
TemplatedSimpleSecondaryVertexComputer.h
Go to the documentation of this file.
1
#ifndef RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
2
#define RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
3
4
#include <cmath>
5
6
#include "Math/GenVector/VectorUtil.h"
7
8
#include "
DataFormats/TrackReco/interface/Track.h
"
9
#include "
DataFormats/VertexReco/interface/Vertex.h
"
10
#include "
DataFormats/BTauReco/interface/TemplatedSecondaryVertexTagInfo.h
"
11
12
#include "
RecoBTau/JetTagComputer/interface/JetTagComputer.h
"
13
14
#include "
RecoBTag/SecondaryVertex/interface/TrackKinematics.h
"
15
16
template
<
class
IPTI,
class
VTX>
17
class
TemplatedSimpleSecondaryVertexComputer
:
public
JetTagComputer
{
18
public
:
19
using
Tokens
=
void
;
20
21
typedef
reco::TemplatedSecondaryVertexTagInfo<IPTI, VTX>
TagInfo
;
22
23
TemplatedSimpleSecondaryVertexComputer
(
const
edm::ParameterSet
&
parameters
)
24
:
use2d
(!
parameters
.getParameter<
bool
>(
"use3d"
)),
25
useSig
(
parameters
.getParameter<
bool
>(
"useSignificance"
)),
26
unBoost
(
parameters
.getParameter<
bool
>(
"unBoost"
)),
27
minTracks
(
parameters
.getParameter<unsigned
int
>(
"minTracks"
)),
28
minVertices_
(1) {
29
uses
(
"svTagInfos"
);
30
minVertices_
=
31
parameters
.existsAs<
unsigned
int
>(
"minVertices"
) ?
parameters
.getParameter<
unsigned
int
>(
"minVertices"
) : 1;
32
}
33
34
float
discriminator
(
const
TagInfoHelper
&
tagInfos
)
const override
{
35
const
TagInfo
&
info
=
tagInfos
.get<
TagInfo
>();
36
if
(
info
.nVertices() <
minVertices_
)
37
return
-1;
38
unsigned
int
idx
= 0;
39
while
(
idx
<
info
.nVertices()) {
40
if
(
info
.nVertexTracks(
idx
) >=
minTracks
)
41
break
;
42
idx
++;
43
}
44
if
(
idx
>=
info
.nVertices())
45
return
-1.0;
46
47
double
gamma
;
48
if
(
unBoost
) {
49
reco::TrackKinematics
kinematics
(
info
.secondaryVertex(
idx
));
50
gamma
=
kinematics
.vectorSum().Gamma();
51
}
else
52
gamma
= 1.0;
53
54
double
value
;
55
if
(
useSig
)
56
value
=
info
.flightDistance(
idx
,
use2d
).significance();
57
else
58
value
=
info
.flightDistance(
idx
,
use2d
).value();
59
60
value
/=
gamma
;
61
62
if
(
useSig
)
63
value
= (
value
> 0) ? +
std::log
(1 +
value
) : -
std::log
(1 -
value
);
64
65
return
value
;
66
}
67
68
private
:
69
bool
use2d
;
70
bool
useSig
;
71
bool
unBoost
;
72
unsigned
int
minTracks
;
73
unsigned
int
minVertices_
;
74
};
75
76
#endif // RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
TemplatedSimpleSecondaryVertexComputer::use2d
bool use2d
Definition:
TemplatedSimpleSecondaryVertexComputer.h:69
reco::TrackKinematics
Definition:
TrackKinematics.h:16
HLT_2018_cff.tagInfos
tagInfos
Definition:
HLT_2018_cff.py:51519
electrons_cff.bool
bool
Definition:
electrons_cff.py:372
TemplatedSimpleSecondaryVertexComputer::unBoost
bool unBoost
Definition:
TemplatedSimpleSecondaryVertexComputer.h:71
reco::TemplatedSecondaryVertexTagInfo
Definition:
TemplatedSecondaryVertexTagInfo.h:47
TemplatedSecondaryVertexTagInfo.h
CustomPhysics_cfi.gamma
gamma
Definition:
CustomPhysics_cfi.py:17
charmTagsComputerCvsB_cfi.idx
idx
Definition:
charmTagsComputerCvsB_cfi.py:108
info
static const TGPicture * info(bool iBackgroundIsBlack)
Definition:
FWCollectionSummaryWidget.cc:152
JetTagComputer::TagInfoHelper
Definition:
JetTagComputer.h:16
parameters
parameters
Definition:
BeamSpot_PayloadInspector.cc:14
JetTagComputer.h
Track.h
TemplatedSimpleSecondaryVertexComputer::TagInfo
reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX > TagInfo
Definition:
TemplatedSimpleSecondaryVertexComputer.h:21
TemplatedSimpleSecondaryVertexComputer::Tokens
void Tokens
Definition:
TemplatedSimpleSecondaryVertexComputer.h:19
Vertex.h
TemplatedSimpleSecondaryVertexComputer::minVertices_
unsigned int minVertices_
Definition:
TemplatedSimpleSecondaryVertexComputer.h:73
edm::ParameterSet
Definition:
ParameterSet.h:36
TemplatedSimpleSecondaryVertexComputer::useSig
bool useSig
Definition:
TemplatedSimpleSecondaryVertexComputer.h:70
createfilelist.int
int
Definition:
createfilelist.py:10
value
Definition:
value.py:1
TrackKinematics.h
JetTagComputer
Definition:
JetTagComputer.h:14
TemplatedSimpleSecondaryVertexComputer::minTracks
unsigned int minTracks
Definition:
TemplatedSimpleSecondaryVertexComputer.h:72
JetTagComputer::uses
void uses(unsigned int id, const std::string &label)
Definition:
JetTagComputer.cc:17
TemplatedSimpleSecondaryVertexComputer::TemplatedSimpleSecondaryVertexComputer
TemplatedSimpleSecondaryVertexComputer(const edm::ParameterSet ¶meters)
Definition:
TemplatedSimpleSecondaryVertexComputer.h:23
relativeConstraints.value
value
Definition:
relativeConstraints.py:53
TemplatedSimpleSecondaryVertexComputer
Definition:
TemplatedSimpleSecondaryVertexComputer.h:17
dqm-mbProfile.log
log
Definition:
dqm-mbProfile.py:17
funct::void
TEMPL(T2) struct Divides void
Definition:
Factorize.h:29
Pyquen_Zmumu_2760GeV_dimuonAcc_cfi.kinematics
kinematics
Definition:
Pyquen_Zmumu_2760GeV_dimuonAcc_cfi.py:19
TemplatedSimpleSecondaryVertexComputer::discriminator
float discriminator(const TagInfoHelper &tagInfos) const override
Definition:
TemplatedSimpleSecondaryVertexComputer.h:34
Generated for CMSSW Reference Manual by
1.8.16