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/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::EDProducer
{
45
public
:
47
explicit
TtJetPartonMatch
(
const
edm::ParameterSet
&);
49
~TtJetPartonMatch
()
override
;
51
void
produce
(
edm::Event
&,
const
edm::EventSetup
&)
override
;
52
53
private
:
55
JetPartonMatching::algorithms
readAlgorithm
(
const
std::string
&
str
);
56
58
C
partons_
;
60
edm::EDGetTokenT<TtGenEvent>
genEvt_
;
62
edm::EDGetTokenT<edm::View<reco::Jet>
>
jets_
;
65
int
maxNJets_
;
68
int
maxNComb_
;
70
JetPartonMatching::algorithms
algorithm_
;
72
bool
useDeltaR_
;
75
bool
useMaxDist_
;
78
double
maxDist_
;
80
int
verbosity_
;
81
};
82
83
template
<
typename
C>
84
TtJetPartonMatch<C>::TtJetPartonMatch
(
const
edm::ParameterSet
&
cfg
)
85
: partons_(
cfg
.getParameter<
std
::
vector
<
std
::
string
> >(
"partonsToIgnore"
)),
86
genEvt_(consumes<
TtGenEvent
>(
edm
::
InputTag
(
"genEvt"
))),
87
jets_(consumes<
edm
::
View
<
reco
::
Jet
> >(
cfg
.getParameter<
edm
::
InputTag
>(
"jets"
))),
88
maxNJets_(
cfg
.getParameter<
int
>(
"maxNJets"
)),
89
maxNComb_(
cfg
.getParameter<
int
>(
"maxNComb"
)),
90
algorithm_(readAlgorithm(
cfg
.getParameter<
std
::
string
>(
"algorithm"
))),
91
useDeltaR_(
cfg
.getParameter<
bool
>(
"useDeltaR"
)),
92
useMaxDist_(
cfg
.getParameter<
bool
>(
"useMaxDist"
)),
93
maxDist_(
cfg
.getParameter<double>(
"maxDist"
)),
94
verbosity_(
cfg
.getParameter<
int
>(
"verbosity"
)) {
95
// produces a vector of jet/lepton indices in the order of
96
// * TtSemiLepEvtPartons
97
// * TtFullHadEvtPartons
98
// * TtFullLepEvtPartons
99
// and vectors of the corresponding quality parameters
100
produces<std::vector<std::vector<int> > >();
101
produces<std::vector<double> >(
"SumPt"
);
102
produces<std::vector<double> >(
"SumDR"
);
103
produces<int>(
"NumberOfConsideredJets"
);
104
}
105
106
template
<
typename
C>
107
TtJetPartonMatch<C>::~TtJetPartonMatch
() {}
108
109
template
<
typename
C>
110
void
TtJetPartonMatch<C>::produce
(
edm::Event
& evt,
const
edm::EventSetup
&
setup
) {
111
// will write
112
// * parton match
113
// * sumPt
114
// * sumDR
115
// to the event
116
std::unique_ptr<std::vector<std::vector<int> > >
match
(
new
std::vector
<std::vector<int> >);
117
std::unique_ptr<std::vector<double> >
sumPt
(
new
std::vector<double>);
118
std::unique_ptr<std::vector<double> >
sumDR
(
new
std::vector<double>);
119
std::unique_ptr<int> pJetsConsidered(
new
int
);
120
121
// get TtGenEvent and jet collection from the event
122
edm::Handle<TtGenEvent>
genEvt
;
123
evt.
getByToken
(genEvt_,
genEvt
);
124
125
edm::Handle<edm::View<reco::Jet>
> topJets;
126
evt.
getByToken
(jets_, topJets);
127
128
// fill vector of partons in the order of
129
// * TtFullLepEvtPartons
130
// * TtSemiLepEvtPartons
131
// * TtFullHadEvtPartons
132
std::vector<const reco::Candidate*>
partons
= partons_.vec(*
genEvt
);
133
134
// prepare vector of jets
135
std::vector<const reco::Candidate*>
jets
;
136
for
(
unsigned
int
ij = 0; ij < topJets->size(); ++ij) {
137
// take all jets if maxNJets_ == -1; otherwise use
138
// maxNJets_ if maxNJets_ is big enough or use same
139
// number of jets as partons if maxNJets_ < number
140
// of partons
141
if
(maxNJets_ != -1) {
142
if
(maxNJets_ >= (
int
)
partons
.size()) {
143
if
((
int
)ij == maxNJets_)
144
break
;
145
}
else
{
146
if
(ij ==
partons
.size())
147
break
;
148
}
149
}
150
jets
.push_back((
const
reco::Candidate
*)&(*topJets)[ij]);
151
}
152
*pJetsConsidered =
jets
.size();
153
154
// do the matching with specified parameters
155
JetPartonMatching
jetPartonMatch(
partons
,
jets
, algorithm_, useMaxDist_, useDeltaR_, maxDist_);
156
157
// print some info for each event
158
// if corresponding verbosity level set
159
if
(verbosity_ > 0)
160
jetPartonMatch.
print
();
161
162
for
(
unsigned
int
ic = 0; ic < jetPartonMatch.
getNumberOfAvailableCombinations
(); ++ic) {
163
if
((
int
)ic >= maxNComb_ && maxNComb_ >= 0)
164
break
;
165
std::vector<int>
matches
= jetPartonMatch.
getMatchesForPartons
(ic);
166
partons_.expand(
matches
);
// insert dummy indices for partons that were chosen to be ignored
167
match
->push_back(
matches
);
168
sumPt
->push_back(jetPartonMatch.
getSumDeltaPt
(ic));
169
sumDR
->push_back(jetPartonMatch.
getSumDeltaR
(ic));
170
}
171
evt.
put
(
std::move
(
match
));
172
evt.
put
(
std::move
(
sumPt
),
"SumPt"
);
173
evt.
put
(
std::move
(
sumDR
),
"SumDR"
);
174
evt.
put
(
std::move
(pJetsConsidered),
"NumberOfConsideredJets"
);
175
}
176
177
template
<
typename
C>
178
JetPartonMatching::algorithms
TtJetPartonMatch<C>::readAlgorithm
(
const
std::string
&
str
) {
179
if
(
str
==
"totalMinDist"
)
180
return
JetPartonMatching::totalMinDist
;
181
else
if
(
str
==
"minSumDist"
)
182
return
JetPartonMatching::minSumDist
;
183
else
if
(
str
==
"ptOrderedMinDist"
)
184
return
JetPartonMatching::ptOrderedMinDist
;
185
else
if
(
str
==
"unambiguousOnly"
)
186
return
JetPartonMatching::unambiguousOnly
;
187
else
188
throw
cms::Exception
(
"Configuration"
) <<
"Chosen algorithm is not supported: "
<<
str
<<
"\n"
;
189
}
190
191
#endif
TtJetPartonMatch::jets_
edm::EDGetTokenT< edm::View< reco::Jet > > jets_
jet collection input
Definition:
TtJetPartonMatch.h:62
JetPartonMatching::getSumDeltaR
double getSumDeltaR(const unsigned int comb=0)
Definition:
JetPartonMatching.h:58
electrons_cff.bool
bool
Definition:
electrons_cff.py:366
EDProducer.h
sistrip::View
View
Definition:
ConstantsForView.h:26
edm::EDGetTokenT< TtGenEvent >
edm
HLT enums.
Definition:
AlignableModifier.h:19
JetPartonMatching::unambiguousOnly
Definition:
JetPartonMatching.h:16
oniaPATMuonsWithTrigger_cff.matches
matches
Definition:
oniaPATMuonsWithTrigger_cff.py:77
TtJetPartonMatch::readAlgorithm
JetPartonMatching::algorithms readAlgorithm(const std::string &str)
convert string for algorithm into corresponding enumerator type
Definition:
TtJetPartonMatch.h:178
HLT_FULL_cff.InputTag
InputTag
Definition:
HLT_FULL_cff.py:89285
TtJetPartonMatch::~TtJetPartonMatch
~TtJetPartonMatch() override
default destructor
Definition:
TtJetPartonMatch.h:107
JetPartonMatching::print
void print()
Definition:
JetPartonMatching.cc:377
JetPartonMatching::getNumberOfAvailableCombinations
unsigned int getNumberOfAvailableCombinations()
Definition:
JetPartonMatching.h:46
TtJetPartonMatch::verbosity_
int verbosity_
verbosity level
Definition:
TtJetPartonMatch.h:80
TtJetPartonMatch
Definition:
TtJetPartonMatch.h:44
TtGenEvent
Class derived from the TopGenEvent for ttbar events.
Definition:
TtGenEvent.h:18
Jet.h
singleTopDQM_cfi.jets
jets
Definition:
singleTopDQM_cfi.py:42
TtJetPartonMatch::useMaxDist_
bool useMaxDist_
Definition:
TtJetPartonMatch.h:75
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
edm::Handle
Definition:
AssociativeIterator.h:50
singleTopDQM_cfi.setup
setup
Definition:
singleTopDQM_cfi.py:37
TtFullHadEvtBuilder_cfi.sumPt
sumPt
Definition:
TtFullHadEvtBuilder_cfi.py:38
JetPartonMatching::getMatchesForPartons
std::vector< int > getMatchesForPartons(const unsigned int comb=0)
Definition:
JetPartonMatching.cc:349
TtJetPartonMatch::produce
void produce(edm::Event &, const edm::EventSetup &) override
write jet parton match objects into the event
Definition:
TtJetPartonMatch.h:110
Jet
Definition:
Jet.py:1
TtJetPartonMatch::maxNJets_
int maxNJets_
Definition:
TtJetPartonMatch.h:65
str
#define str(s)
Definition:
TestProcessor.cc:52
TtJetPartonMatch::maxNComb_
int maxNComb_
Definition:
TtJetPartonMatch.h:68
dqmAnalyzer_cff.partons
partons
Definition:
dqmAnalyzer_cff.py:27
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
Event.h:535
TtJetPartonMatch::useDeltaR_
bool useDeltaR_
switch to choose between deltaR/deltaTheta matching
Definition:
TtJetPartonMatch.h:72
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
JetPartonMatching::algorithms
algorithms
Definition:
JetPartonMatching.h:16
TtFullHadEvtBuilder_cfi.sumDR
sumDR
Definition:
TtFullHadEvtBuilder_cfi.py:39
edm::ParameterSet
Definition:
ParameterSet.h:47
Event.h
TtGenEvent.h
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition:
Utils.h:10
createfilelist.int
int
Definition:
createfilelist.py:10
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
edm::EventSetup
Definition:
EventSetup.h:58
JetPartonMatching::getSumDeltaPt
double getSumDeltaPt(const unsigned int comb=0)
Definition:
JetPartonMatching.h:57
looper.cfg
cfg
Definition:
looper.py:297
reco::Candidate
Definition:
Candidate.h:27
TtJetPartonMatch::maxDist_
double maxDist_
Definition:
TtJetPartonMatch.h:78
eostools.move
def move(src, dest)
Definition:
eostools.py:511
std
Definition:
JetResolutionObject.h:76
JetPartonMatching::totalMinDist
Definition:
JetPartonMatching.h:16
JetPartonMatching.h
gen::C
C
Definition:
PomwigHadronizer.cc:78
Frameworkfwd.h
JetPartonMatching::minSumDist
Definition:
JetPartonMatching.h:16
Exception
Definition:
hltDiff.cc:245
TtJetPartonMatch::TtJetPartonMatch
TtJetPartonMatch(const edm::ParameterSet &)
default conructor
Definition:
TtJetPartonMatch.h:84
edm::EDProducer
Definition:
EDProducer.h:35
TtJetPartonMatch::partons_
C partons_
partons
Definition:
TtJetPartonMatch.h:58
JetPartonMatching::ptOrderedMinDist
Definition:
JetPartonMatching.h:16
TtJetPartonMatch::algorithm_
JetPartonMatching::algorithms algorithm_
choice of algorithm
Definition:
TtJetPartonMatch.h:70
TtJetPartonMatch::genEvt_
edm::EDGetTokenT< TtGenEvent > genEvt_
TtGenEvent collection input.
Definition:
TtJetPartonMatch.h:60
ParameterSet.h
edm::Event
Definition:
Event.h:73
JetPartonMatching
Definition:
JetPartonMatching.h:11
TtGenEvtProducer_cfi.genEvt
genEvt
Definition:
TtGenEvtProducer_cfi.py:7
Generated for CMSSW Reference Manual by
1.8.16