src
TopQuarkAnalysis
TopJetCombination
plugins
TtSemiLepJetCombWMassMaxSumPt.cc
Go to the documentation of this file.
1
#include "
FWCore/Framework/interface/Event.h
"
2
#include "
FWCore/Framework/interface/global/EDProducer.h
"
3
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
4
5
#include "
AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h
"
6
#include "
DataFormats/PatCandidates/interface/Jet.h
"
7
8
class
TtSemiLepJetCombWMassMaxSumPt
:
public
edm::global::EDProducer
<> {
9
public
:
10
explicit
TtSemiLepJetCombWMassMaxSumPt
(
const
edm::ParameterSet
&);
11
12
private
:
13
void
produce
(
edm::StreamID
,
edm::Event
& evt,
const
edm::EventSetup
&
setup
)
const override
;
14
15
bool
isValid
(
const
int
&
idx
,
const
std::vector<pat::Jet>&
jets
)
const
{
16
return
(0 <=
idx
&&
idx
< (
int
)
jets
.size());
17
};
18
19
edm::EDGetTokenT<std::vector<pat::Jet>
>
jetsToken_
;
20
edm::EDGetTokenT<edm::View<reco::RecoCandidate>
>
lepsToken_
;
21
int
maxNJets_
;
22
double
wMass_
;
23
bool
useBTagging_
;
24
std::string
bTagAlgorithm_
;
25
double
minBDiscBJets_
;
26
double
maxBDiscLightJets_
;
27
};
28
29
TtSemiLepJetCombWMassMaxSumPt::TtSemiLepJetCombWMassMaxSumPt
(
const
edm::ParameterSet
&
cfg
)
30
: jetsToken_(consumes<
std
::
vector
<
pat
::
Jet
>>(
cfg
.getParameter<
edm
::
InputTag
>(
"jets"
))),
31
lepsToken_(consumes<
edm
::
View
<
reco
::RecoCandidate>>(
cfg
.getParameter<
edm
::
InputTag
>(
"leps"
))),
32
maxNJets_(
cfg
.getParameter<
int
>(
"maxNJets"
)),
33
wMass_(
cfg
.getParameter<double>(
"wMass"
)),
34
useBTagging_(
cfg
.getParameter<
bool
>(
"useBTagging"
)),
35
bTagAlgorithm_(
cfg
.getParameter<
std
::
string
>(
"bTagAlgorithm"
)),
36
minBDiscBJets_(
cfg
.getParameter<double>(
"minBDiscBJets"
)),
37
maxBDiscLightJets_(
cfg
.getParameter<double>(
"maxBDiscLightJets"
)) {
38
if
(
maxNJets_
< 4 &&
maxNJets_
!= -1)
39
throw
cms::Exception
(
"WrongConfig"
) <<
"Parameter maxNJets can not be set to "
<<
maxNJets_
<<
". \n"
40
<<
"It has to be larger than 4 or can be set to -1 to take all jets."
;
41
42
produces<std::vector<std::vector<int>>>();
43
produces<int>(
"NumberOfConsideredJets"
);
44
}
45
46
void
TtSemiLepJetCombWMassMaxSumPt::produce
(
edm::StreamID
,
edm::Event
& evt,
const
edm::EventSetup
&
setup
)
const
{
47
auto
pOut = std::make_unique<std::vector<std::vector<int>>>();
48
auto
pJetsConsidered = std::make_unique<int>(0);
49
50
std::vector<int>
match
;
51
for
(
unsigned
int
i
= 0;
i
< 4; ++
i
)
52
match
.push_back(-1);
53
54
// get jets
55
const
auto
&
jets
= evt.
get
(
jetsToken_
);
56
57
// get leptons
58
const
auto
&
leps
= evt.
get
(
lepsToken_
);
59
60
// skip events without lepton candidate or less than 4 jets or no MET
61
if
(
leps
.empty() ||
jets
.size() < 4) {
62
pOut->push_back(
match
);
63
evt.
put
(
std::move
(pOut));
64
*pJetsConsidered =
jets
.size();
65
evt.
put
(
std::move
(pJetsConsidered),
"NumberOfConsideredJets"
);
66
return
;
67
}
68
69
unsigned
maxNJets
=
maxNJets_
;
70
if
(
maxNJets_
== -1 || (
int
)
jets
.size() <
maxNJets_
)
71
maxNJets
=
jets
.size();
72
*pJetsConsidered =
maxNJets
;
73
evt.
put
(
std::move
(pJetsConsidered),
"NumberOfConsideredJets"
);
74
75
std::vector<bool> isBJet;
76
std::vector<bool> isLJet;
77
int
cntBJets = 0;
78
if
(
useBTagging_
) {
79
for
(
unsigned
int
idx
= 0;
idx
<
maxNJets
; ++
idx
) {
80
isBJet.push_back((
jets
[
idx
].bDiscriminator(
bTagAlgorithm_
) >
minBDiscBJets_
));
81
isLJet.push_back((
jets
[
idx
].bDiscriminator(
bTagAlgorithm_
) <
maxBDiscLightJets_
));
82
if
(
jets
[
idx
].bDiscriminator(
bTagAlgorithm_
) >
minBDiscBJets_
)
83
cntBJets++;
84
}
85
}
86
87
// -----------------------------------------------------
88
// associate those jets that get closest to the W mass
89
// with their invariant mass to the hadronic W boson
90
// -----------------------------------------------------
91
double
wDist = -1.;
92
std::vector<int> closestToWMassIndices;
93
closestToWMassIndices.push_back(-1);
94
closestToWMassIndices.push_back(-1);
95
for
(
unsigned
idx
= 0;
idx
<
maxNJets
; ++
idx
) {
96
if
(
useBTagging_
&& (!isLJet[
idx
] || (cntBJets <= 2 && isBJet[
idx
])))
97
continue
;
98
for
(
unsigned
jdx = (
idx
+ 1); jdx <
maxNJets
; ++jdx) {
99
if
(
useBTagging_
&&
100
(!isLJet[jdx] || (cntBJets <= 2 && isBJet[jdx]) || (cntBJets == 3 && isBJet[
idx
] && isBJet[jdx])))
101
continue
;
102
reco::Particle::LorentzVector
sum =
jets
[
idx
].p4() +
jets
[jdx].p4();
103
if
(wDist < 0. || wDist > fabs(sum.mass() -
wMass_
)) {
104
wDist = fabs(sum.mass() -
wMass_
);
105
closestToWMassIndices.clear();
106
closestToWMassIndices.push_back(
idx
);
107
closestToWMassIndices.push_back(jdx);
108
}
109
}
110
}
111
112
// -----------------------------------------------------
113
// associate those jets with maximum pt of the vectorial
114
// sum to the hadronic decay chain
115
// -----------------------------------------------------
116
double
maxPt
= -1.;
117
int
hadB = -1;
118
if
(
isValid
(closestToWMassIndices[0],
jets
) &&
isValid
(closestToWMassIndices[1],
jets
)) {
119
for
(
unsigned
idx
= 0;
idx
<
maxNJets
; ++
idx
) {
120
if
(
useBTagging_
&& !isBJet[
idx
])
121
continue
;
122
// make sure it's not used up already from the hadronic W
123
if
((
int
)
idx
!= closestToWMassIndices[0] && (
int
)
idx
!= closestToWMassIndices[1]) {
124
reco::Particle::LorentzVector
sum =
125
jets
[closestToWMassIndices[0]].p4() +
jets
[closestToWMassIndices[1]].p4() +
jets
[
idx
].p4();
126
if
(
maxPt
< 0. ||
maxPt
< sum.pt()) {
127
maxPt
= sum.pt();
128
hadB =
idx
;
129
}
130
}
131
}
132
}
133
134
// -----------------------------------------------------
135
// associate the remaining jet with maximum pt of the
136
// vectorial sum with the leading lepton with the
137
// leptonic b quark
138
// -----------------------------------------------------
139
maxPt
= -1.;
140
int
lepB = -1;
141
for
(
unsigned
idx
= 0;
idx
<
maxNJets
; ++
idx
) {
142
if
(
useBTagging_
&& !isBJet[
idx
])
143
continue
;
144
// make sure it's not used up already from the hadronic decay chain
145
if
((
int
)
idx
!= closestToWMassIndices[0] && (
int
)
idx
!= closestToWMassIndices[1] && (
int
)
idx
!= hadB) {
146
reco::Particle::LorentzVector
sum =
jets
[
idx
].p4() +
leps
[0].p4();
147
if
(
maxPt
< 0. ||
maxPt
< sum.pt()) {
148
maxPt
= sum.pt();
149
lepB =
idx
;
150
}
151
}
152
}
153
154
match
[
TtSemiLepEvtPartons::LightQ
] = closestToWMassIndices[0];
155
match
[
TtSemiLepEvtPartons::LightQBar
] = closestToWMassIndices[1];
156
match
[
TtSemiLepEvtPartons::HadB
] = hadB;
157
match
[
TtSemiLepEvtPartons::LepB
] = lepB;
158
159
pOut->push_back(
match
);
160
evt.
put
(
std::move
(pOut));
161
}
162
163
#include "
FWCore/Framework/interface/MakerMacros.h
"
164
DEFINE_FWK_MODULE
(
TtSemiLepJetCombWMassMaxSumPt
);
SUSYDQMAnalyzer_cfi.maxNJets
maxNJets
Definition:
SUSYDQMAnalyzer_cfi.py:13
EDProducer.h
edm::Event::put
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition:
Event.h:133
mps_fire.i
i
Definition:
mps_fire.py:428
TtSemiLepJetCombWMassMaxSumPt::wMass_
double wMass_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:22
TtSemiLepJetCombWMassMaxSumPt::TtSemiLepJetCombWMassMaxSumPt
TtSemiLepJetCombWMassMaxSumPt(const edm::ParameterSet &)
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:29
edm::Event::get
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition:
Event.h:346
TtSemiLepEvtPartons.h
Exception
Definition:
hltDiff.cc:245
heavyIonCSV_trainingSettings.idx
idx
Definition:
heavyIonCSV_trainingSettings.py:5
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
edm::StreamID
Definition:
StreamID.h:30
Event.h
HLT_2022v12_cff.InputTag
InputTag
Definition:
HLT_2022v12_cff.py:65221
TtSemiLepJetCombWMassMaxSumPt::minBDiscBJets_
double minBDiscBJets_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:25
MakerMacros.h
TtSemiLepJetCombWMassMaxSumPt::bTagAlgorithm_
std::string bTagAlgorithm_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:24
PDWG_EXODelayedJetMET_cff.jets
jets
Definition:
PDWG_EXODelayedJetMET_cff.py:14
std
Definition:
JetResolutionObject.h:76
edm::EDGetTokenT
Definition:
EDGetToken.h:33
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
pat
Definition:
HeavyIon.h:7
TtSemiLepEvtPartons::HadB
Definition:
TtSemiLepEvtPartons.h:25
ParameterSet.h
createfilelist.int
int
Definition:
createfilelist.py:10
Jet
Definition:
Jet.py:1
TtSemiLepJetCombWMassMaxSumPt::useBTagging_
bool useBTagging_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:23
TtSemiLepJetCombWMassMaxSumPt::jetsToken_
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:17
edm::EventSetup
Definition:
EventSetup.h:59
edm::global::EDProducer
Definition:
EDProducer.h:32
singleTopDQM_cfi.setup
setup
SETUP
Definition:
singleTopDQM_cfi.py:37
TtSemiLepEvtPartons::LepB
Definition:
TtSemiLepEvtPartons.h:25
looper.cfg
cfg
Definition:
looper.py:296
electrons_cff.bool
bool
Definition:
electrons_cff.py:381
L1TPhase2Offline_cfi.maxPt
maxPt
Definition:
L1TPhase2Offline_cfi.py:18
TtSemiLepHitFitProducer_Electrons_cfi.leps
leps
Definition:
TtSemiLepHitFitProducer_Electrons_cfi.py:5
TtSemiLepEvtPartons::LightQ
Definition:
TtSemiLepEvtPartons.h:25
TtSemiLepEvtPartons::LightQBar
Definition:
TtSemiLepEvtPartons.h:25
TtSemiLepJetCombWMassMaxSumPt::maxBDiscLightJets_
double maxBDiscLightJets_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:26
TtSemiLepJetCombWMassMaxSumPt
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:8
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:46
edm
HLT enums.
Definition:
AlignableModifier.h:19
trackerHitRTTI::vector
Definition:
trackerHitRTTI.h:21
TtSemiLepJetCombWMassMaxSumPt::maxNJets_
int maxNJets_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:21
Jet.h
edm::ParameterSet
Definition:
ParameterSet.h:47
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition:
Utils.h:10
TtSemiLepJetCombWMassMaxSumPt::lepsToken_
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:20
edm::Event
Definition:
Event.h:73
reco::Particle::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition:
Particle.h:21
sistrip::View
View
Definition:
ConstantsForView.h:26
TtSemiLepJetCombWMassMaxSumPt::produce
void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &setup) const override
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:46
eostools.move
def move(src, dest)
Definition:
eostools.py:511
TtSemiLepJetCombWMassMaxSumPt::isValid
bool isValid(const int &idx, const std::vector< pat::Jet > &jets) const
Definition:
TtSemiLepJetCombWMassMaxSumPt.cc:15
Generated for CMSSW Reference Manual by
1.8.14