src
TopQuarkAnalysis
TopTools
plugins
TtJetPartonMatch.h
Go to the documentation of this file.
1
#ifndef TtJetPartonMatch_h
2
#define TtJetPartonMatch_h
3
4
#include <memory>
5
#include <vector>
6
7
#include "
FWCore/Framework/interface/Event.h
"
8
#include "
FWCore/Framework/interface/global/EDProducer.h
"
9
#include "
FWCore/Framework/interface/Frameworkfwd.h
"
10
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
11
12
#include "
DataFormats/JetReco/interface/Jet.h
"
13
#include "
AnalysisDataFormats/TopObjects/interface/TtGenEvent.h
"
14
#include "
TopQuarkAnalysis/TopTools/interface/JetPartonMatching.h
"
15
16
// ----------------------------------------------------------------------
17
// template class for:
18
//
19
// * TtFullHadJetPartonMatch
20
// * TtSemiLepJetPartonMatch
21
// * TtFullLepJetPartonMatch
22
//
23
// the class provides plugins for jet-parton matching corresponding
24
// to the JetPartonMatching class; expected input is one of the
25
// classes in:
26
//
27
// AnalysisDataFormats/TopObjects/interface/TtFullHadEvtPartons.h
28
// AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h
29
// AnalysisDataFormats/TopObjects/interface/TtFullLepEvtPartons.h
30
//
31
// output is:
32
// * a vector of vectors containing the indices of the jets in the
33
// input collection matched to the partons in the order defined in
34
// the corresponding Tt*EvtPartons class
35
// * a set of vectors with quality parameters of the matching
36
//
37
// the matching can be performed on an arbitrary number of jet
38
// combinations; per default the combination which matches best
39
// according to the quality parameters will be stored; the length
40
// of the vectors will be 1 then
41
// ----------------------------------------------------------------------
42
43
template
<
typename
C>
44
class
TtJetPartonMatch
:
public
edm::global::EDProducer
<> {
45
public
:
47
explicit
TtJetPartonMatch
(
const
edm::ParameterSet
&);
49
void
produce
(
edm::StreamID
,
edm::Event
&,
const
edm::EventSetup
&)
const override
;
50
51
private
:
53
JetPartonMatching::algorithms
readAlgorithm
(
const
std::string
&
str
)
const
;
54
56
C
partons_
;
58
edm::EDGetTokenT<TtGenEvent>
genEvt_
;
60
edm::EDGetTokenT<edm::View<reco::Jet>
>
jets_
;
63
int
maxNJets_
;
66
int
maxNComb_
;
68
JetPartonMatching::algorithms
algorithm_
;
70
bool
useDeltaR_
;
73
bool
useMaxDist_
;
76
double
maxDist_
;
78
int
verbosity_
;
79
};
80
81
template
<
typename
C>
82
TtJetPartonMatch<C>::TtJetPartonMatch
(
const
edm::ParameterSet
&
cfg
)
83
: partons_(
cfg
.getParameter<
std
::
vector
<
std
::
string
>>(
"partonsToIgnore"
)),
84
genEvt_(consumes<
TtGenEvent
>(
edm
::
InputTag
(
"genEvt"
))),
85
jets_(consumes<
edm
::
View
<
reco
::
Jet
>>(
cfg
.getParameter<
edm
::
InputTag
>(
"jets"
))),
86
maxNJets_(
cfg
.getParameter<
int
>(
"maxNJets"
)),
87
maxNComb_(
cfg
.getParameter<
int
>(
"maxNComb"
)),
88
algorithm_(readAlgorithm(
cfg
.getParameter<
std
::
string
>(
"algorithm"
))),
89
useDeltaR_(
cfg
.getParameter<
bool
>(
"useDeltaR"
)),
90
useMaxDist_(
cfg
.getParameter<
bool
>(
"useMaxDist"
)),
91
maxDist_(
cfg
.getParameter<double>(
"maxDist"
)),
92
verbosity_(
cfg
.getParameter<
int
>(
"verbosity"
)) {
93
// produces a vector of jet/lepton indices in the order of
94
// * TtSemiLepEvtPartons
95
// * TtFullHadEvtPartons
96
// * TtFullLepEvtPartons
97
// and vectors of the corresponding quality parameters
98
produces<std::vector<std::vector<int>>>();
99
produces<std::vector<double>>(
"SumPt"
);
100
produces<std::vector<double>>(
"SumDR"
);
101
produces<int>(
"NumberOfConsideredJets"
);
102
}
103
104
template
<
typename
C>
105
void
TtJetPartonMatch<C>::produce
(
edm::StreamID
,
edm::Event
& evt,
const
edm::EventSetup
&
setup
)
const
{
106
// will write
107
// * parton match
108
// * sumPt
109
// * sumDR
110
// to the event
111
auto
match
= std::make_unique<std::vector<std::vector<int>>>();
112
auto
sumPt
= std::make_unique<std::vector<double>>();
113
auto
sumDR
= std::make_unique<std::vector<double>>();
114
std::unique_ptr<int> pJetsConsidered(
new
int
);
115
116
// get TtGenEvent and jet collection from the event
117
const
TtGenEvent
&
genEvt
= evt.
get
(genEvt_);
118
119
const
edm::View<reco::Jet>
& topJets = evt.
get
(jets_);
120
121
// fill vector of partons in the order of
122
// * TtFullLepEvtPartons
123
// * TtSemiLepEvtPartons
124
// * TtFullHadEvtPartons
125
std::vector<const reco::Candidate*>
partons
= partons_.vec(
genEvt
);
126
127
// prepare vector of jets
128
std::vector<const reco::Candidate*>
jets
;
129
for
(
unsigned
int
ij = 0; ij < topJets.
size
(); ++ij) {
130
// take all jets if maxNJets_ == -1; otherwise use
131
// maxNJets_ if maxNJets_ is big enough or use same
132
// number of jets as partons if maxNJets_ < number
133
// of partons
134
if
(maxNJets_ != -1) {
135
if
(maxNJets_ >= (
int
)
partons
.size()) {
136
if
((
int
)ij == maxNJets_)
137
break
;
138
}
else
{
139
if
(ij ==
partons
.size())
140
break
;
141
}
142
}
143
jets
.push_back(&topJets[ij]);
144
}
145
*pJetsConsidered =
jets
.size();
146
147
// do the matching with specified parameters
148
JetPartonMatching
jetPartonMatch(
partons
,
jets
, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
149
150
// print some info for each event
151
// if corresponding verbosity level set
152
if
(verbosity_ > 0)
153
jetPartonMatch.
print
();
154
155
for
(
unsigned
int
ic = 0; ic < jetPartonMatch.
getNumberOfAvailableCombinations
(); ++ic) {
156
if
((
int
)ic >= maxNComb_ && maxNComb_ >= 0)
157
break
;
158
std::vector<int>
matches
= jetPartonMatch.
getMatchesForPartons
(ic);
159
partons_.expand(
matches
);
// insert dummy indices for partons that were chosen to be ignored
160
match
->push_back(
matches
);
161
sumPt
->push_back(jetPartonMatch.
getSumDeltaPt
(ic));
162
sumDR
->push_back(jetPartonMatch.
getSumDeltaR
(ic));
163
}
164
evt.
put
(
std::move
(
match
));
165
evt.
put
(
std::move
(
sumPt
),
"SumPt"
);
166
evt.
put
(
std::move
(
sumDR
),
"SumDR"
);
167
evt.
put
(
std::move
(pJetsConsidered),
"NumberOfConsideredJets"
);
168
}
169
170
template
<
typename
C>
171
JetPartonMatching::algorithms
TtJetPartonMatch<C>::readAlgorithm
(
const
std::string
&
str
)
const
{
172
if
(
str
==
"totalMinDist"
)
173
return
JetPartonMatching::totalMinDist
;
174
else
if
(
str
==
"minSumDist"
)
175
return
JetPartonMatching::minSumDist
;
176
else
if
(
str
==
"ptOrderedMinDist"
)
177
return
JetPartonMatching::ptOrderedMinDist
;
178
else
if
(
str
==
"unambiguousOnly"
)
179
return
JetPartonMatching::unambiguousOnly
;
180
else
181
throw
cms::Exception
(
"Configuration"
) <<
"Chosen algorithm is not supported: "
<<
str
<<
"\n"
;
182
}
183
184
#endif
TtJetPartonMatch::readAlgorithm
JetPartonMatching::algorithms readAlgorithm(const std::string &str) const
convert string for algorithm into corresponding enumerator type
Definition:
TtJetPartonMatch.h:171
ProducerED_cfi.InputTag
InputTag
Definition:
ProducerED_cfi.py:5
TtJetPartonMatch::maxNJets_
int maxNJets_
Definition:
TtJetPartonMatch.h:63
edm::Event::put
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition:
Event.h:133
TtJetPartonMatch::partons_
C partons_
partons
Definition:
TtJetPartonMatch.h:56
edm::Event::get
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition:
Event.h:347
TtJetPartonMatch::genEvt_
edm::EDGetTokenT< TtGenEvent > genEvt_
TtGenEvent collection input.
Definition:
TtJetPartonMatch.h:58
Exception
Definition:
hltDiff.cc:245
edm::StreamID
Definition:
StreamID.h:30
Event.h
JetPartonMatching::minSumDist
Definition:
JetPartonMatching.h:16
TtGenEvent.h
PDWG_EXODelayedJetMET_cff.jets
jets
Definition:
PDWG_EXODelayedJetMET_cff.py:14
std
Definition:
JetResolutionObject.h:76
JetPartonMatching::algorithms
algorithms
Definition:
JetPartonMatching.h:16
edm::EDGetTokenT< TtGenEvent >
Frameworkfwd.h
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
JetPartonMatching::print
void print()
Definition:
JetPartonMatching.cc:377
edm::View::size
size_type size() const
JetPartonMatching::getNumberOfAvailableCombinations
unsigned int getNumberOfAvailableCombinations()
Definition:
JetPartonMatching.h:46
TtFullHadEvtBuilder_cfi.sumDR
sumDR
Definition:
TtFullHadEvtBuilder_cfi.py:39
ParameterSet.h
dqmAnalyzer_cff.partons
partons
Definition:
dqmAnalyzer_cff.py:25
JetPartonMatching::totalMinDist
Definition:
JetPartonMatching.h:16
edm::View
Definition:
CaloClusterFwd.h:14
TtJetPartonMatch::useMaxDist_
bool useMaxDist_
Definition:
TtJetPartonMatch.h:73
TtJetPartonMatch
Definition:
TtJetPartonMatch.h:44
TtGenEvent
Class derived from the TopGenEvent for ttbar events.
Definition:
TtGenEvent.h:18
createfilelist.int
int
Definition:
createfilelist.py:10
Jet
Definition:
Jet.py:1
JetPartonMatching::unambiguousOnly
Definition:
JetPartonMatching.h:16
JetPartonMatching::getMatchesForPartons
std::vector< int > getMatchesForPartons(const unsigned int comb=0)
Definition:
JetPartonMatching.cc:349
TtJetPartonMatch::useDeltaR_
bool useDeltaR_
switch to choose between deltaR/deltaTheta matching
Definition:
TtJetPartonMatch.h:70
edm::EventSetup
Definition:
EventSetup.h:56
Jet.h
edm::global::EDProducer
Definition:
EDProducer.h:32
singleTopDQM_cfi.setup
setup
SETUP
Definition:
singleTopDQM_cfi.py:37
correctionTermsCaloMet_cff.C
C
Definition:
correctionTermsCaloMet_cff.py:34
looper.cfg
cfg
Definition:
looper.py:296
JetPartonMatching
Definition:
JetPartonMatching.h:11
JetPartonMatching::ptOrderedMinDist
Definition:
JetPartonMatching.h:16
TtJetPartonMatch::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
write jet parton match objects into the event
Definition:
TtJetPartonMatch.h:105
TtJetPartonMatch::maxDist_
double maxDist_
Definition:
TtJetPartonMatch.h:76
nano_mu_local_reco_cff.bool
bool
Definition:
nano_mu_local_reco_cff.py:13
TtJetPartonMatch::algorithm_
JetPartonMatching::algorithms algorithm_
choice of algorithm
Definition:
TtJetPartonMatch.h:68
EDProducer.h
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:46
edm
HLT enums.
Definition:
AlignableModifier.h:19
trackerHitRTTI::vector
Definition:
trackerHitRTTI.h:21
TtJetPartonMatch::jets_
edm::EDGetTokenT< edm::View< reco::Jet > > jets_
jet collection input
Definition:
TtJetPartonMatch.h:60
TtJetPartonMatch::verbosity_
int verbosity_
verbosity level
Definition:
TtJetPartonMatch.h:78
edm::ParameterSet
Definition:
ParameterSet.h:47
JetPartonMatching::getSumDeltaR
double getSumDeltaR(const unsigned int comb=0)
Definition:
JetPartonMatching.h:58
oniaPATMuonsWithTrigger_cff.matches
matches
Definition:
oniaPATMuonsWithTrigger_cff.py:77
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition:
Utils.h:10
TtJetPartonMatch::TtJetPartonMatch
TtJetPartonMatch(const edm::ParameterSet &)
default conructor
Definition:
TtJetPartonMatch.h:82
TtJetPartonMatch::maxNComb_
int maxNComb_
Definition:
TtJetPartonMatch.h:66
edm::Event
Definition:
Event.h:73
str
#define str(s)
Definition:
TestProcessor.cc:56
JetPartonMatching.h
sistrip::View
View
Definition:
ConstantsForView.h:26
JetPartonMatching::getSumDeltaPt
double getSumDeltaPt(const unsigned int comb=0)
Definition:
JetPartonMatching.h:57
eostools.move
def move(src, dest)
Definition:
eostools.py:511
TtFullHadEvtBuilder_cfi.sumPt
sumPt
Definition:
TtFullHadEvtBuilder_cfi.py:38
TtGenEvtProducer_cfi.genEvt
genEvt
Definition:
TtGenEvtProducer_cfi.py:7
Generated for CMSSW Reference Manual by
1.8.14