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