TopQuarkAnalysis
TopJetCombination
plugins
TtSemiLepJetCombWMassDeltaTopMass.cc
Go to the documentation of this file.
1
#include "
TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepJetCombWMassDeltaTopMass.h
"
2
3
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
4
5
#include "
TopQuarkAnalysis/TopTools/interface/MEzCalculator.h
"
6
7
TtSemiLepJetCombWMassDeltaTopMass::TtSemiLepJetCombWMassDeltaTopMass
(
const
edm::ParameterSet
&
cfg
)
8
: jetsToken_(consumes<
std
::
vector
<
pat
::
Jet
> >(
cfg
.getParameter<
edm
::
InputTag
>(
"jets"
))),
9
lepsToken_(consumes<
edm
::
View
<
reco
::RecoCandidate> >(
cfg
.getParameter<
edm
::
InputTag
>(
"leps"
))),
10
metsToken_(consumes<
std
::
vector
<
pat
::
MET
> >(
cfg
.getParameter<
edm
::
InputTag
>(
"mets"
))),
11
maxNJets_(
cfg
.getParameter<
int
>(
"maxNJets"
)),
12
wMass_(
cfg
.getParameter<double>(
"wMass"
)),
13
useBTagging_(
cfg
.getParameter<
bool
>(
"useBTagging"
)),
14
bTagAlgorithm_(
cfg
.getParameter<
std
::
string
>(
"bTagAlgorithm"
)),
15
minBDiscBJets_(
cfg
.getParameter<double>(
"minBDiscBJets"
)),
16
maxBDiscLightJets_(
cfg
.getParameter<double>(
"maxBDiscLightJets"
)),
17
neutrinoSolutionType_(
cfg
.getParameter<
int
>(
"neutrinoSolutionType"
)) {
18
if
(
maxNJets_
< 4 &&
maxNJets_
!= -1)
19
throw
cms::Exception
(
"WrongConfig"
) <<
"Parameter maxNJets can not be set to "
<<
maxNJets_
<<
". \n"
20
<<
"It has to be larger than 4 or can be set to -1 to take all jets."
;
21
22
produces<std::vector<std::vector<int> > >();
23
produces<int>(
"NumberOfConsideredJets"
);
24
}
25
26
TtSemiLepJetCombWMassDeltaTopMass::~TtSemiLepJetCombWMassDeltaTopMass
() {}
27
28
void
TtSemiLepJetCombWMassDeltaTopMass::produce
(
edm::Event
& evt,
const
edm::EventSetup
&
setup
) {
29
std::unique_ptr<std::vector<std::vector<int> > > pOut(
new
std::vector
<std::vector<int> >);
30
std::unique_ptr<int> pJetsConsidered(
new
int
);
31
32
std::vector<int>
match
;
33
for
(
unsigned
int
i
= 0;
i
< 4; ++
i
)
34
match
.push_back(-1);
35
36
// get jets
37
edm::Handle<std::vector<pat::Jet>
>
jets
;
38
evt.
getByToken
(
jetsToken_
,
jets
);
39
40
// get leptons
41
edm::Handle<edm::View<reco::RecoCandidate>
>
leps
;
42
evt.
getByToken
(
lepsToken_
,
leps
);
43
44
// get MET
45
edm::Handle<std::vector<pat::MET>
>
mets
;
46
evt.
getByToken
(
metsToken_
,
mets
);
47
48
// skip events without lepton candidate or less than 4 jets or no MET
49
if
(
leps
->empty() ||
jets
->size() < 4 ||
mets
->empty()) {
50
pOut->push_back(
match
);
51
evt.
put
(
std::move
(pOut));
52
*pJetsConsidered =
jets
->size();
53
evt.
put
(
std::move
(pJetsConsidered),
"NumberOfConsideredJets"
);
54
return
;
55
}
56
57
unsigned
maxNJets
=
maxNJets_
;
58
if
(
maxNJets_
== -1 || (
int
)
jets
->size() <
maxNJets_
)
59
maxNJets
=
jets
->size();
60
*pJetsConsidered =
maxNJets
;
61
evt.
put
(
std::move
(pJetsConsidered),
"NumberOfConsideredJets"
);
62
63
std::vector<bool> isBJet;
64
std::vector<bool> isLJet;
65
int
cntBJets = 0;
66
if
(
useBTagging_
) {
67
for
(
unsigned
int
idx
= 0;
idx
<
maxNJets
; ++
idx
) {
68
isBJet.push_back(((*
jets
)[
idx
].bDiscriminator(
bTagAlgorithm_
) >
minBDiscBJets_
));
69
isLJet.push_back(((*
jets
)[
idx
].bDiscriminator(
bTagAlgorithm_
) <
maxBDiscLightJets_
));
70
if
((*
jets
)[
idx
].bDiscriminator(
bTagAlgorithm_
) >
minBDiscBJets_
)
71
cntBJets++;
72
}
73
}
74
75
// -----------------------------------------------------
76
// associate those jets that get closest to the W mass
77
// with their invariant mass to the hadronic W boson
78
// -----------------------------------------------------
79
double
wDist = -1.;
80
std::vector<int> closestToWMassIndices;
81
closestToWMassIndices.push_back(-1);
82
closestToWMassIndices.push_back(-1);
83
for
(
unsigned
idx
= 0;
idx
<
maxNJets
; ++
idx
) {
84
if
(
useBTagging_
&& (!isLJet[
idx
] || (cntBJets <= 2 && isBJet[
idx
])))
85
continue
;
86
for
(
unsigned
jdx = (
idx
+ 1); jdx <
maxNJets
; ++jdx) {
87
if
(
useBTagging_
&&
88
(!isLJet[jdx] || (cntBJets <= 2 && isBJet[jdx]) || (cntBJets == 3 && isBJet[
idx
] && isBJet[jdx])))
89
continue
;
90
reco::Particle::LorentzVector
sum = (*jets)[
idx
].p4() + (*jets)[jdx].p4();
91
if
(wDist < 0. || wDist > fabs(sum.mass() -
wMass_
)) {
92
wDist = fabs(sum.mass() -
wMass_
);
93
closestToWMassIndices.clear();
94
closestToWMassIndices.push_back(
idx
);
95
closestToWMassIndices.push_back(jdx);
96
}
97
}
98
}
99
100
// -----------------------------------------------------
101
// build a leptonic W boson from the lepton and the MET
102
// -----------------------------------------------------
103
double
neutrino_pz = 0.;
104
if
(
neutrinoSolutionType_
!= -1) {
105
MEzCalculator
mez;
106
mez.
SetMET
(*
mets
->begin());
107
if
(dynamic_cast<const reco::Muon*>(&(
leps
->front())))
108
mez.
SetLepton
((*
leps
)[0],
true
);
109
else
if
(dynamic_cast<const reco::GsfElectron*>(&(
leps
->front())))
110
mez.
SetLepton
((*
leps
)[0],
false
);
111
else
112
throw
cms::Exception
(
"UnimplementedFeature"
)
113
<<
"Type of lepton given together with MET for solving neutrino kinematics is neither muon nor electron.\n"
;
114
neutrino_pz = mez.
Calculate
(
neutrinoSolutionType_
);
115
}
116
const
math::XYZTLorentzVector
neutrino(
mets
->begin()->px(),
117
mets
->begin()->py(),
118
neutrino_pz,
119
sqrt
(
mets
->begin()->px() *
mets
->begin()->px() +
120
mets
->begin()->py() *
mets
->begin()->py() + neutrino_pz * neutrino_pz));
121
const
reco::Particle::LorentzVector
lepW = neutrino +
leps
->front().p4();
122
123
// -----------------------------------------------------
124
// associate those jets to the hadronic and the leptonic
125
// b quark that minimize the difference between
126
// hadronic and leptonic top mass
127
// -----------------------------------------------------
128
double
deltaTop = -1.;
129
int
hadB = -1;
130
int
lepB = -1;
131
if
(
isValid
(closestToWMassIndices[0],
jets
) &&
isValid
(closestToWMassIndices[1],
jets
)) {
132
const
reco::Particle::LorentzVector
hadW =
133
(*jets)[closestToWMassIndices[0]].p4() + (*jets)[closestToWMassIndices[1]].p4();
134
// find hadronic b candidate
135
for
(
unsigned
idx
= 0;
idx
<
maxNJets
; ++
idx
) {
136
if
(
useBTagging_
&& !isBJet[
idx
])
137
continue
;
138
// make sure it's not used up already from the hadronic W
139
if
((
int
)
idx
!= closestToWMassIndices[0] && (
int
)
idx
!= closestToWMassIndices[1]) {
140
reco::Particle::LorentzVector
hadTop = hadW + (*jets)[
idx
].p4();
141
// find leptonic b candidate
142
for
(
unsigned
jdx = 0; jdx <
maxNJets
; ++jdx) {
143
if
(
useBTagging_
&& !isBJet[jdx])
144
continue
;
145
// make sure it's not used up already from the hadronic branch
146
if
((
int
)jdx != closestToWMassIndices[0] && (
int
)jdx != closestToWMassIndices[1] && jdx !=
idx
) {
147
reco::Particle::LorentzVector
lepTop = lepW + (*jets)[jdx].p4();
148
if
(deltaTop < 0. || deltaTop > fabs(hadTop.mass() - lepTop.mass())) {
149
deltaTop = fabs(hadTop.mass() - lepTop.mass());
150
hadB =
idx
;
151
lepB = jdx;
152
}
153
}
154
}
155
}
156
}
157
}
158
159
match
[
TtSemiLepEvtPartons::LightQ
] = closestToWMassIndices[0];
160
match
[
TtSemiLepEvtPartons::LightQBar
] = closestToWMassIndices[1];
161
match
[
TtSemiLepEvtPartons::HadB
] = hadB;
162
match
[
TtSemiLepEvtPartons::LepB
] = lepB;
163
164
pOut->push_back(
match
);
165
evt.
put
(
std::move
(pOut));
166
}
electrons_cff.bool
bool
Definition:
electrons_cff.py:393
mps_fire.i
i
Definition:
mps_fire.py:428
TtSemiLepEvtPartons::LightQBar
Definition:
TtSemiLepEvtPartons.h:25
TtSemiLepJetCombWMassDeltaTopMass::useBTagging_
bool useBTagging_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:32
sistrip::View
View
Definition:
ConstantsForView.h:26
TtSemiLepEvtPartons::LepB
Definition:
TtSemiLepEvtPartons.h:25
TtSemiLepJetCombWMassDeltaTopMass::produce
void produce(edm::Event &evt, const edm::EventSetup &setup) override
Definition:
TtSemiLepJetCombWMassDeltaTopMass.cc:28
edm
HLT enums.
Definition:
AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition:
HLT_FULL_cff.py:89287
TtSemiLepJetCombWMassDeltaTopMass::maxBDiscLightJets_
double maxBDiscLightJets_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:35
singleTopDQM_cfi.jets
jets
Definition:
singleTopDQM_cfi.py:42
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
singleTopDQM_cfi.mets
mets
Definition:
singleTopDQM_cfi.py:43
TtSemiLepJetCombWMassDeltaTopMass::neutrinoSolutionType_
int neutrinoSolutionType_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:36
MEzCalculator
Definition:
MEzCalculator.h:18
edm::Handle
Definition:
AssociativeIterator.h:50
singleTopDQM_cfi.setup
setup
Definition:
singleTopDQM_cfi.py:37
reco::Particle::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition:
Particle.h:21
heavyIonCSV_trainingSettings.idx
idx
Definition:
heavyIonCSV_trainingSettings.py:5
TtSemiLepHitFitProducer_Electrons_cfi.leps
leps
Definition:
TtSemiLepHitFitProducer_Electrons_cfi.py:5
TtSemiLepJetCombWMassDeltaTopMass::maxNJets_
int maxNJets_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:30
TtSemiLepJetCombWMassDeltaTopMass::bTagAlgorithm_
std::string bTagAlgorithm_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:33
Jet
Definition:
Jet.py:1
TtSemiLepJetCombWMassDeltaTopMass::metsToken_
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:29
mathSSE::sqrt
T sqrt(T t)
Definition:
SSEVec.h:19
MEzCalculator.h
TtSemiLepJetCombWMassDeltaTopMass::isValid
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:23
MEzCalculator::SetLepton
void SetLepton(const pat::Particle &lepton, bool isMuon=true)
Set lepton.
Definition:
MEzCalculator.h:31
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
Event.h:531
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
TtSemiLepJetCombWMassDeltaTopMass::minBDiscBJets_
double minBDiscBJets_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:34
edm::ParameterSet
Definition:
ParameterSet.h:47
MEzCalculator::Calculate
double Calculate(int type=1)
member functions
Definition:
MEzCalculator.cc:14
TtSemiLepEvtPartons::HadB
Definition:
TtSemiLepEvtPartons.h:25
TtSemiLepJetCombWMassDeltaTopMass::wMass_
double wMass_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:31
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition:
Utils.h:10
TtSemiLepEvtPartons::LightQ
Definition:
TtSemiLepEvtPartons.h:25
createfilelist.int
int
Definition:
createfilelist.py:10
TtSemiLepJetCombWMassDeltaTopMass.h
SUSYDQMAnalyzer_cfi.maxNJets
maxNJets
Definition:
SUSYDQMAnalyzer_cfi.py:13
edm::Event::put
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition:
Event.h:133
trackerHitRTTI::vector
Definition:
trackerHitRTTI.h:21
TtSemiLepJetCombWMassDeltaTopMass::TtSemiLepJetCombWMassDeltaTopMass
TtSemiLepJetCombWMassDeltaTopMass(const edm::ParameterSet &)
Definition:
TtSemiLepJetCombWMassDeltaTopMass.cc:7
edm::EventSetup
Definition:
EventSetup.h:57
pat
Definition:
HeavyIon.h:7
TtSemiLepJetCombWMassDeltaTopMass::lepsToken_
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:28
TtSemiLepJetCombWMassDeltaTopMass::~TtSemiLepJetCombWMassDeltaTopMass
~TtSemiLepJetCombWMassDeltaTopMass() override
Definition:
TtSemiLepJetCombWMassDeltaTopMass.cc:26
looper.cfg
cfg
Definition:
looper.py:297
MEzCalculator::SetMET
void SetMET(const pat::MET &MET)
Set MET.
Definition:
MEzCalculator.h:25
eostools.move
def move(src, dest)
Definition:
eostools.py:511
std
Definition:
JetResolutionObject.h:76
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition:
LorentzVector.h:29
Exception
Definition:
hltDiff.cc:246
MET
ParameterSet.h
edm::Event
Definition:
Event.h:73
TtSemiLepJetCombWMassDeltaTopMass::jetsToken_
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
Definition:
TtSemiLepJetCombWMassDeltaTopMass.h:25
Generated for CMSSW Reference Manual by
1.8.16